ROMS
Loading...
Searching...
No Matches
obc_adjust_mod Module Reference

Functions/Subroutines

subroutine, public obc_adjust (ng, tile, linp)
 
subroutine obc_adjust_tile (ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, imins, imaxs, jmins, jmaxs, linp)
 
subroutine, public load_obc (ng, tile, lout)
 
subroutine load_obc_tile (ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, imins, imaxs, jmins, jmaxs, lout)
 

Function/Subroutine Documentation

◆ load_obc()

subroutine, public obc_adjust_mod::load_obc ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lout )

Definition at line 502 of file obc_adjust.F.

503!
504!=======================================================================
505! !
506! This routine loads open boundaries into nonlinear storage arrays. !
507! In 4DVAR open boundary adjustment, the boundary values are stored !
508! in arrays with extra dimensions to aid minimization at other times !
509! in addition to initialization time. !
510! !
511! On Input: !
512! !
513! ng Nested grid number. !
514! tile Domain partition. !
515! Lout Time index to process in storage arrays. !
516! !
517!=======================================================================
518!
519 USE mod_param
520!
521! Imported variable declarations.
522!
523 integer, intent(in) :: ng, tile, Lout
524!
525! Local variable declarations.
526!
527 character (len=*), parameter :: MyFile = &
528 & __FILE__//", load_obc"
529!
530# include "tile.h"
531!
532# ifdef PROFILE
533 CALL wclock_on (ng, inlm, 8, __line__, myfile)
534# endif
535 CALL load_obc_tile (ng, tile, &
536 & lbi, ubi, lbj, ubj, lbij, ubij, &
537 & imins, imaxs, jmins, jmaxs, &
538 & lout)
539# ifdef PROFILE
540 CALL wclock_off (ng, inlm, 8, __line__, myfile)
541# endif
542!
543 RETURN
integer, parameter inlm
Definition mod_param.F:662
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_param::inlm, load_obc_tile(), wclock_off(), and wclock_on().

Referenced by load_obc_tile(), and main3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ load_obc_tile()

subroutine obc_adjust_mod::load_obc_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) lout )
private

Definition at line 547 of file obc_adjust.F.

551!***********************************************************************
552!
553 USE mod_param
554 USE mod_boundary
555 USE mod_ncparam
556 USE mod_scalars
557!
558! Imported variable declarations.
559!
560 integer, intent(in) :: ng, tile
561 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
562 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
563 integer, intent(in) :: Lout
564!
565! Local variable declarations.
566!
567 integer :: i, ic, j
568# ifdef SOLVE3D
569 integer :: it, k
570# endif
571
572# include "set_bounds.h"
573!
574!-----------------------------------------------------------------------
575! Load nonlinear open boundary fields into storage arrays.
576!-----------------------------------------------------------------------
577!
578 load_obc : IF (mod(iic(ng)-1,nobc(ng)).eq.0) THEN
579 obccount(ng)=obccount(ng)+1
580 ic=obccount(ng)
581!
582! Free-surface open boundaries.
583!
584 IF (lbc(iwest,isfsur,ng)%acquire.and. &
585 & lobc(iwest,isfsur,ng).and. &
586 & domain(ng)%Western_Edge(tile)) THEN
587 DO j=jstr,jend
588 boundary(ng)%zeta_obc(j,iwest,ic,lout)= &
589 & boundary(ng)%zeta_west(j)
590 END DO
591 END IF
592
593 IF (lbc(ieast,isfsur,ng)%acquire.and. &
594 & lobc(ieast,isfsur,ng).and. &
595 & domain(ng)%Eastern_Edge(tile)) THEN
596 DO j=jstr,jend
597 boundary(ng)%zeta_obc(j,ieast,ic,lout)= &
598 & boundary(ng)%zeta_east(j)
599 END DO
600 END IF
601
602 IF (lbc(isouth,isfsur,ng)%acquire.and. &
603 & lobc(isouth,isfsur,ng).and. &
604 & domain(ng)%Southern_Edge(tile)) THEN
605 DO i=istr,iend
606 boundary(ng)%zeta_obc(i,isouth,ic,lout)= &
607 & boundary(ng)%zeta_south(i)
608 END DO
609 END IF
610
611 IF (lbc(inorth,isfsur,ng)%acquire.and. &
612 & lobc(inorth,isfsur,ng).and. &
613 & domain(ng)%Northern_Edge(tile)) THEN
614 DO i=istr,iend
615 boundary(ng)%zeta_obc(i,inorth,ic,lout)= &
616 & boundary(ng)%zeta_north(i)
617 END DO
618 END IF
619!
620! 2D U-momentum open boundaries.
621!
622 IF (lbc(iwest,isubar,ng)%acquire.and. &
623 & lobc(iwest,isubar,ng).and. &
624 & domain(ng)%Western_Edge(tile)) THEN
625 DO j=jstr,jend
626 boundary(ng)%ubar_obc(j,iwest,ic,lout)= &
627 & boundary(ng)%ubar_west(j)
628 END DO
629 END IF
630
631 IF (lbc(ieast,isubar,ng)%acquire.and. &
632 & lobc(ieast,isubar,ng).and. &
633 & domain(ng)%Eastern_Edge(tile)) THEN
634 DO j=jstr,jend
635 boundary(ng)%ubar_obc(j,ieast,ic,lout)= &
636 & boundary(ng)%ubar_east(j)
637 END DO
638 END IF
639
640 IF (lbc(isouth,isubar,ng)%acquire.and. &
641 & lobc(isouth,isubar,ng).and. &
642 & domain(ng)%Southern_Edge(tile)) THEN
643 DO i=istru,iend
644 boundary(ng)%ubar_obc(i,isouth,ic,lout)= &
645 & boundary(ng)%ubar_south(i)
646 END DO
647 END IF
648
649 IF (lbc(inorth,isubar,ng)%acquire.and. &
650 & lobc(inorth,isubar,ng).and. &
651 & domain(ng)%Northern_Edge(tile)) THEN
652 DO i=istru,iend
653 boundary(ng)%ubar_obc(i,inorth,ic,lout)= &
654 & boundary(ng)%ubar_north(i)
655 END DO
656 END IF
657!
658! 2D V-momentum open boundaries.
659!
660 IF (lbc(iwest,isvbar,ng)%acquire.and. &
661 & lobc(iwest,isvbar,ng).and. &
662 & domain(ng)%Western_Edge(tile)) THEN
663 DO j=jstrv,jend
664 boundary(ng)%vbar_obc(j,iwest,ic,lout)= &
665 & boundary(ng)%vbar_west(j)
666 END DO
667 END IF
668
669 IF (lbc(ieast,isvbar,ng)%acquire.and. &
670 & lobc(ieast,isvbar,ng).and. &
671 & domain(ng)%Eastern_Edge(tile)) THEN
672 DO j=jstrv,jend
673 boundary(ng)%vbar_obc(j,ieast,ic,lout)= &
674 & boundary(ng)%vbar_east(j)
675 END DO
676 END IF
677
678 IF (lbc(isouth,isvbar,ng)%acquire.and. &
679 & lobc(isouth,isvbar,ng).and. &
680 & domain(ng)%Southern_Edge(tile)) THEN
681 DO i=istr,iend
682 boundary(ng)%vbar_obc(i,isouth,ic,lout)= &
683 & boundary(ng)%vbar_south(i)
684 END DO
685 END IF
686
687 IF (lbc(inorth,isvbar,ng)%acquire.and. &
688 & lobc(inorth,isvbar,ng).and. &
689 & domain(ng)%Northern_Edge(tile)) THEN
690 DO i=istr,iend
691 boundary(ng)%vbar_obc(i,inorth,ic,lout)= &
692 & boundary(ng)%vbar_north(i)
693 END DO
694 END IF
695
696# ifdef SOLVE3D
697!
698! 3D U-momentum open boundaries.
699!
700 IF (lbc(iwest,isuvel,ng)%acquire.and. &
701 & lobc(iwest,isuvel,ng).and. &
702 & domain(ng)%Western_Edge(tile)) THEN
703 DO k=1,n(ng)
704 DO j=jstr,jend
705 boundary(ng)%u_obc(j,k,iwest,ic,lout)= &
706 & boundary(ng)%u_west(j,k)
707 END DO
708 END DO
709 END IF
710
711 IF (lbc(ieast,isuvel,ng)%acquire.and. &
712 & lobc(ieast,isuvel,ng).and. &
713 & domain(ng)%Eastern_Edge(tile)) THEN
714 DO k=1,n(ng)
715 DO j=jstr,jend
716 boundary(ng)%u_obc(j,k,ieast,ic,lout)= &
717 & boundary(ng)%u_east(j,k)
718 END DO
719 END DO
720 END IF
721
722 IF (lbc(isouth,isuvel,ng)%acquire.and. &
723 & lobc(isouth,isuvel,ng).and. &
724 & domain(ng)%Southern_Edge(tile)) THEN
725 DO k=1,n(ng)
726 DO i=istru,iend
727 boundary(ng)%u_obc(i,k,isouth,ic,lout)= &
728 & boundary(ng)%u_south(i,k)
729 END DO
730 END DO
731 END IF
732
733 IF (lbc(inorth,isuvel,ng)%acquire.and. &
734 & lobc(inorth,isuvel,ng).and. &
735 & domain(ng)%Northern_Edge(tile)) THEN
736 DO k=1,n(ng)
737 DO i=istru,iend
738 boundary(ng)%u_obc(i,k,inorth,ic,lout)= &
739 & boundary(ng)%u_north(i,k)
740 END DO
741 END DO
742 END IF
743!
744! 3D V-momentum open boundaries.
745!
746 IF (lbc(iwest,isvvel,ng)%acquire.and. &
747 & lobc(iwest,isvvel,ng).and. &
748 & domain(ng)%Western_Edge(tile)) THEN
749 DO k=1,n(ng)
750 DO j=jstrv,jend
751 boundary(ng)%v_obc(j,k,iwest,ic,lout)= &
752 & boundary(ng)%v_west(j,k)
753 END DO
754 END DO
755 END IF
756
757 IF (lbc(ieast,isvvel,ng)%acquire.and. &
758 & lobc(ieast,isvvel,ng).and. &
759 & domain(ng)%Eastern_Edge(tile)) THEN
760 DO k=1,n(ng)
761 DO j=jstrv,jend
762 boundary(ng)%v_obc(j,k,ieast,ic,lout)= &
763 & boundary(ng)%v_east(j,k)
764 END DO
765 END DO
766 END IF
767
768 IF (lbc(isouth,isvvel,ng)%acquire.and. &
769 & lobc(isouth,isvvel,ng).and. &
770 & domain(ng)%Southern_Edge(tile)) THEN
771 DO k=1,n(ng)
772 DO i=istr,iend
773 boundary(ng)%v_obc(i,k,isouth,ic,lout)= &
774 & boundary(ng)%v_south(i,k)
775 END DO
776 END DO
777 END IF
778
779 IF (lbc(inorth,isvvel,ng)%acquire.and. &
780 & lobc(inorth,isvvel,ng).and. &
781 & domain(ng)%Northern_Edge(tile)) THEN
782 DO k=1,n(ng)
783 DO i=istr,iend
784 boundary(ng)%v_obc(i,k,inorth,ic,lout)= &
785 & boundary(ng)%v_north(i,k)
786 END DO
787 END DO
788 END IF
789!
790! Tracers open boundaries.
791!
792 DO it=1,nt(ng)
793 IF (lbc(iwest,istvar(it),ng)%acquire.and. &
794 & lobc(iwest,istvar(it),ng).and. &
795 & domain(ng)%Western_Edge(tile)) THEN
796 DO k=1,n(ng)
797 DO j=jstr,jend
798 boundary(ng)%t_obc(j,k,iwest,ic,lout,it)= &
799 & boundary(ng)%t_west(j,k,it)
800 END DO
801 END DO
802 END IF
803
804 IF (lbc(ieast,istvar(it),ng)%acquire.and. &
805 & lobc(ieast,istvar(it),ng).and. &
806 & domain(ng)%Eastern_Edge(tile)) THEN
807 DO k=1,n(ng)
808 DO j=jstr,jend
809 boundary(ng)%t_obc(j,k,ieast,ic,lout,it)= &
810 & boundary(ng)%t_east(j,k,it)
811 END DO
812 END DO
813 END IF
814
815 IF (lbc(isouth,istvar(it),ng)%acquire.and. &
816 & lobc(isouth,istvar(it),ng).and. &
817 & domain(ng)%Southern_Edge(tile)) THEN
818 DO k=1,n(ng)
819 DO i=istr,iend
820 boundary(ng)%t_obc(i,k,isouth,ic,lout,it)= &
821 & boundary(ng)%t_south(i,k,it)
822 END DO
823 END DO
824 END IF
825
826 IF (lbc(inorth,istvar(it),ng)%acquire.and. &
827 & lobc(inorth,istvar(it),ng).and. &
828 & domain(ng)%Northern_Edge(tile)) THEN
829 DO k=1,n(ng)
830 DO i=istr,iend
831 boundary(ng)%t_obc(i,k,inorth,ic,lout,it)= &
832 & boundary(ng)%t_north(i,k,it)
833 END DO
834 END DO
835 END IF
836 END DO
837# endif
838 END IF load_obc
839!
840 RETURN
type(t_boundary), dimension(:), allocatable boundary
integer isvvel
integer isvbar
integer, dimension(:), allocatable istvar
integer isuvel
integer isfsur
integer isubar
integer, dimension(:), allocatable n
Definition mod_param.F:479
type(t_lbc), dimension(:,:,:), allocatable lbc
Definition mod_param.F:375
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer, dimension(:), allocatable obccount
integer, dimension(:), allocatable iic
integer, dimension(:), allocatable nobc
logical, dimension(:,:,:), allocatable lobc
integer, parameter iwest
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth

References mod_boundary::boundary, mod_param::domain, mod_scalars::ieast, mod_scalars::iic, mod_scalars::inorth, mod_ncparam::isfsur, mod_scalars::isouth, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::iwest, mod_param::lbc, load_obc(), mod_scalars::lobc, mod_param::n, mod_scalars::nobc, mod_param::nt, and mod_scalars::obccount.

Referenced by load_obc().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ obc_adjust()

subroutine, public obc_adjust_mod::obc_adjust ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) linp )

Definition at line 35 of file obc_adjust.F.

36!***********************************************************************
37!
38 USE mod_param
39!
40! Imported variable declarations.
41!
42 integer, intent(in) :: ng, tile, Linp
43!
44! Local variable declarations.
45!
46 character (len=*), parameter :: MyFile = &
47 & __FILE__
48!
49# include "tile.h"
50!
51# ifdef PROFILE
52 CALL wclock_on (ng, inlm, 7, __line__, myfile)
53# endif
54 CALL obc_adjust_tile (ng, tile, &
55 & lbi, ubi, lbj, ubj, lbij, ubij, &
56 & imins, imaxs, jmins, jmaxs, &
57 & linp)
58# ifdef PROFILE
59 CALL wclock_off (ng, inlm, 7, __line__, myfile)
60# endif
61!
62 RETURN

References mod_param::inlm, obc_adjust_tile(), wclock_off(), and wclock_on().

Referenced by main3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ obc_adjust_tile()

subroutine obc_adjust_mod::obc_adjust_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) linp )
private

Definition at line 66 of file obc_adjust.F.

70!***********************************************************************
71!
72 USE mod_param
73 USE mod_boundary
74 USE mod_ncparam
75 USE mod_scalars
76!
77! Imported variable declarations.
78!
79 integer, intent(in) :: ng, tile
80 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
81 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
82 integer, intent(in) :: Linp
83!
84! Local variable declarations.
85!
86 integer :: i, it1, it2, j
87# ifdef SOLVE3D
88 integer :: k, it
89# endif
90 real(r8) :: fac, fac1, fac2
91
92# include "set_bounds.h"
93!
94!-----------------------------------------------------------------------
95! Adjust nonlinear open boundaries with 4DVar increments.
96!
97! HGA: need to set fac2=0, when the boundary condition in the NLM
98! are radiations until we know what to do.
99!-----------------------------------------------------------------------
100!
101! Set time records and interpolation factor, if any.
102!
103 IF (nbrec(ng).eq.1) THEN
104 it1=1
105 it2=1
106 fac1=1.0_r8
107 fac2=0.0_r8
108 ELSE
109# ifdef GENERIC_DSTART
110 it1=max(0,(iic(ng)-ntstart(ng))/nobc(ng))+1
111# else
112 it1=max(0,(iic(ng)-1)/nobc(ng))+1
113# endif
114 it2=min(it1+1,nbrec(ng))
115 fac1=obc_time(it2,ng)-(time(ng)+dt(ng))
116 fac2=(time(ng)+dt(ng))-obc_time(it1,ng)
117 fac=1.0_r8/(fac1+fac2)
118 fac1=fac*fac1
119 fac2=fac*fac2
120 END IF
121!
122! Free-surface open boundaries.
123!
124 IF (lbc(iwest,isfsur,ng)%acquire.and. &
125 & lobc(iwest,isfsur,ng).and. &
126 & domain(ng)%Western_Edge(tile)) THEN
127 DO j=jstr,jend
128 boundary(ng)%zeta_west(j)=boundary(ng)%zeta_west(j)+ &
129 & fac1* &
130 & boundary(ng)%tl_zeta_obc(j, &
131 & iwest,it1,linp)+ &
132 & fac2* &
133 & boundary(ng)%tl_zeta_obc(j, &
134 & iwest,it2,linp)
135 END DO
136 END IF
137
138 IF (lbc(ieast,isfsur,ng)%acquire.and. &
139 & lobc(ieast,isfsur,ng).and. &
140 & domain(ng)%Eastern_Edge(tile)) THEN
141 DO j=jstr,jend
142 boundary(ng)%zeta_east(j)=boundary(ng)%zeta_east(j)+ &
143 & fac1* &
144 & boundary(ng)%tl_zeta_obc(j, &
145 & ieast,it1,linp)+ &
146 & fac2* &
147 & boundary(ng)%tl_zeta_obc(j, &
148 & ieast,it2,linp)
149 END DO
150 END IF
151
152 IF (lbc(isouth,isfsur,ng)%acquire.and. &
153 & lobc(isouth,isfsur,ng).and. &
154 & domain(ng)%Southern_Edge(tile)) THEN
155 DO i=istr,iend
156 boundary(ng)%zeta_south(i)=boundary(ng)%zeta_south(i)+ &
157 & fac1* &
158 & boundary(ng)%tl_zeta_obc(i, &
159 & isouth,it1,linp)+ &
160 & fac2* &
161 & boundary(ng)%tl_zeta_obc(i, &
162 & isouth,it2,linp)
163 END DO
164 END IF
165
166 IF (lbc(inorth,isfsur,ng)%acquire.and. &
167 & lobc(inorth,isfsur,ng).and. &
168 & domain(ng)%Northern_Edge(tile)) THEN
169 DO i=istr,iend
170 boundary(ng)%zeta_north(i)=boundary(ng)%zeta_north(i)+ &
171 & fac1* &
172 & boundary(ng)%tl_zeta_obc(i, &
173 & inorth,it1,linp)+ &
174 & fac2* &
175 & boundary(ng)%tl_zeta_obc(i, &
176 & inorth,it2,linp)
177 END DO
178 END IF
179!
180! 2D U-momentum open boundaries.
181!
182 IF (lbc(iwest,isubar,ng)%acquire.and. &
183 & lobc(iwest,isubar,ng).and. &
184 & domain(ng)%Western_Edge(tile)) THEN
185 DO j=jstr,jend
186 boundary(ng)%ubar_west(j)=boundary(ng)%ubar_west(j)+ &
187 & fac1* &
188 & boundary(ng)%tl_ubar_obc(j, &
189 & iwest,it1,linp)+ &
190 & fac2* &
191 & boundary(ng)%tl_ubar_obc(j, &
192 & iwest,it2,linp)
193 END DO
194 END IF
195
196 IF (lbc(ieast,isubar,ng)%acquire.and. &
197 & lobc(ieast,isubar,ng).and. &
198 & domain(ng)%Eastern_Edge(tile)) THEN
199 DO j=jstr,jend
200 boundary(ng)%ubar_east(j)=boundary(ng)%ubar_east(j)+ &
201 & fac1* &
202 & boundary(ng)%tl_ubar_obc(j, &
203 & ieast,it1,linp)+ &
204 & fac2* &
205 & boundary(ng)%tl_ubar_obc(j, &
206 & ieast,it2,linp)
207 END DO
208 END IF
209
210 IF (lbc(isouth,isubar,ng)%acquire.and. &
211 & lobc(isouth,isubar,ng).and. &
212 & domain(ng)%Southern_Edge(tile)) THEN
213 DO i=istru,iend
214 boundary(ng)%ubar_south(i)=boundary(ng)%ubar_south(i)+ &
215 & fac1* &
216 & boundary(ng)%tl_ubar_obc(i, &
217 & isouth,it1,linp)+ &
218 & fac2* &
219 & boundary(ng)%tl_ubar_obc(i, &
220 & isouth,it2,linp)
221 END DO
222 END IF
223
224 IF (lbc(inorth,isubar,ng)%acquire.and. &
225 & lobc(inorth,isubar,ng).and. &
226 & domain(ng)%Northern_Edge(tile)) THEN
227 DO i=istru,iend
228 boundary(ng)%ubar_north(i)=boundary(ng)%ubar_north(i)+ &
229 & fac1* &
230 & boundary(ng)%tl_ubar_obc(i, &
231 & inorth,it1,linp)+ &
232 & fac2* &
233 & boundary(ng)%tl_ubar_obc(i, &
234 & inorth,it2,linp)
235 END DO
236 END IF
237!
238! 2D V-momentum open boundaries.
239!
240 IF (lbc(iwest,isvbar,ng)%acquire.and. &
241 & lobc(iwest,isvbar,ng).and. &
242 & domain(ng)%Western_Edge(tile)) THEN
243 DO j=jstrv,jend
244 boundary(ng)%vbar_west(j)=boundary(ng)%vbar_west(j)+ &
245 & fac1* &
246 & boundary(ng)%tl_vbar_obc(j, &
247 & iwest,it1,linp)+ &
248 & fac2* &
249 & boundary(ng)%tl_vbar_obc(j, &
250 & iwest,it2,linp)
251 END DO
252 END IF
253
254 IF (lbc(ieast,isvbar,ng)%acquire.and. &
255 & lobc(ieast,isvbar,ng).and. &
256 & domain(ng)%Eastern_Edge(tile)) THEN
257 DO j=jstrv,jend
258 boundary(ng)%vbar_east(j)=boundary(ng)%vbar_east(j)+ &
259 & fac1* &
260 & boundary(ng)%tl_vbar_obc(j, &
261 & ieast,it1,linp)+ &
262 & fac2* &
263 & boundary(ng)%tl_vbar_obc(j, &
264 & ieast,it2,linp)
265 END DO
266 END IF
267
268 IF (lbc(isouth,isvbar,ng)%acquire.and. &
269 & lobc(isouth,isvbar,ng).and. &
270 & domain(ng)%Southern_Edge(tile)) THEN
271 DO i=istr,iend
272 boundary(ng)%vbar_south(i)=boundary(ng)%vbar_south(i)+ &
273 & fac1* &
274 & boundary(ng)%tl_vbar_obc(i, &
275 & isouth,it1,linp)+ &
276 & fac2* &
277 & boundary(ng)%tl_vbar_obc(i, &
278 & isouth,it2,linp)
279 END DO
280 END IF
281
282 IF (lbc(inorth,isvbar,ng)%acquire.and. &
283 & lobc(inorth,isvbar,ng).and. &
284 & domain(ng)%Northern_Edge(tile)) THEN
285 DO i=istr,iend
286 boundary(ng)%vbar_north(i)=boundary(ng)%vbar_north(i)+ &
287 & fac1* &
288 & boundary(ng)%tl_vbar_obc(i, &
289 & inorth,it1,linp)+ &
290 & fac2* &
291 & boundary(ng)%tl_vbar_obc(i, &
292 & inorth,it2,linp)
293 END DO
294 END IF
295
296# ifdef SOLVE3D
297!
298! 3D U-momentum open boundaries.
299!
300 IF (lbc(iwest,isuvel,ng)%acquire.and. &
301 & lobc(iwest,isuvel,ng).and. &
302 & domain(ng)%Western_Edge(tile)) THEN
303 DO k=1,n(ng)
304 DO j=jstr,jend
305 boundary(ng)%u_west(j,k)=boundary(ng)%u_west(j,k)+ &
306 & fac1* &
307 & boundary(ng)%tl_u_obc(j,k, &
308 & iwest,it1,linp)+ &
309 & fac2* &
310 & boundary(ng)%tl_u_obc(j,k, &
311 & iwest,it2,linp)
312 END DO
313 END DO
314 END IF
315
316 IF (lbc(ieast,isuvel,ng)%acquire.and. &
317 & lobc(ieast,isuvel,ng).and. &
318 & domain(ng)%Eastern_Edge(tile)) THEN
319 DO k=1,n(ng)
320 DO j=jstr,jend
321 boundary(ng)%u_east(j,k)=boundary(ng)%u_east(j,k)+ &
322 & fac1* &
323 & boundary(ng)%tl_u_obc(j,k, &
324 & ieast,it1,linp)+ &
325 & fac2* &
326 & boundary(ng)%tl_u_obc(j,k, &
327 & ieast,it2,linp)
328 END DO
329 END DO
330 END IF
331
332 IF (lbc(isouth,isuvel,ng)%acquire.and. &
333 & lobc(isouth,isuvel,ng).and. &
334 & domain(ng)%Southern_Edge(tile)) THEN
335 DO k=1,n(ng)
336 DO i=istru,iend
337 boundary(ng)%u_south(i,k)=boundary(ng)%u_south(i,k)+ &
338 & fac1* &
339 & boundary(ng)%tl_u_obc(i,k, &
340 & isouth,it1,linp)+ &
341 & fac2* &
342 & boundary(ng)%tl_u_obc(i,k, &
343 & isouth,it2,linp)
344 END DO
345 END DO
346 END IF
347
348 IF (lbc(inorth,isuvel,ng)%acquire.and. &
349 & lobc(inorth,isuvel,ng).and. &
350 & domain(ng)%Northern_Edge(tile)) THEN
351 DO k=1,n(ng)
352 DO i=istru,iend
353 boundary(ng)%u_north(i,k)=boundary(ng)%u_north(i,k)+ &
354 & fac1* &
355 & boundary(ng)%tl_u_obc(i,k, &
356 & inorth,it1,linp)+ &
357 & fac2* &
358 & boundary(ng)%tl_u_obc(i,k, &
359 & inorth,it2,linp)
360 END DO
361 END DO
362 END IF
363!
364! 3D V-momentum open boundaries.
365!
366 IF (lbc(iwest,isvvel,ng)%acquire.and. &
367 & lobc(iwest,isvvel,ng).and. &
368 & domain(ng)%Western_Edge(tile)) THEN
369 DO k=1,n(ng)
370 DO j=jstrv,jend
371 boundary(ng)%v_west(j,k)=boundary(ng)%v_west(j,k)+ &
372 & fac1* &
373 & boundary(ng)%tl_v_obc(j,k, &
374 & iwest,it1,linp)+ &
375 & fac2* &
376 & boundary(ng)%tl_v_obc(j,k, &
377 & iwest,it2,linp)
378 END DO
379 END DO
380 END IF
381
382 IF (lbc(ieast,isvvel,ng)%acquire.and. &
383 & lobc(ieast,isvvel,ng).and. &
384 & domain(ng)%Eastern_Edge(tile)) THEN
385 DO k=1,n(ng)
386 DO j=jstrv,jend
387 boundary(ng)%v_east(j,k)=boundary(ng)%v_east(j,k)+ &
388 & fac1* &
389 & boundary(ng)%tl_v_obc(j,k, &
390 & ieast,it1,linp)+ &
391 & fac2* &
392 & boundary(ng)%tl_v_obc(j,k, &
393 & ieast,it2,linp)
394 END DO
395 END DO
396 END IF
397
398 IF (lbc(isouth,isvvel,ng)%acquire.and. &
399 & lobc(isouth,isvvel,ng).and. &
400 & domain(ng)%Southern_Edge(tile)) THEN
401 DO k=1,n(ng)
402 DO i=istr,iend
403 boundary(ng)%v_south(i,k)=boundary(ng)%v_south(i,k)+ &
404 & fac1* &
405 & boundary(ng)%tl_v_obc(i,k, &
406 & isouth,it1,linp)+ &
407 & fac2* &
408 & boundary(ng)%tl_v_obc(i,k, &
409 & isouth,it2,linp)
410 END DO
411 END DO
412 END IF
413
414 IF (lbc(inorth,isvvel,ng)%acquire.and. &
415 & lobc(inorth,isvvel,ng).and. &
416 & domain(ng)%Northern_Edge(tile)) THEN
417 DO k=1,n(ng)
418 DO i=istr,iend
419 boundary(ng)%v_north(i,k)=boundary(ng)%v_north(i,k)+ &
420 & fac1* &
421 & boundary(ng)%tl_v_obc(i,k, &
422 & inorth,it1,linp)+ &
423 & fac2* &
424 & boundary(ng)%tl_v_obc(i,k, &
425 & inorth,it2,linp)
426 END DO
427 END DO
428 END IF
429!
430! Tracers open boundaries.
431!
432 DO it=1,nt(ng)
433 IF (lbc(iwest,istvar(it),ng)%acquire.and. &
434 & lobc(iwest,istvar(it),ng).and. &
435 & domain(ng)%Western_Edge(tile)) THEN
436 DO k=1,n(ng)
437 DO j=jstr,jend
438 boundary(ng)%t_west(j,k,it)=boundary(ng)%t_west(j,k,it)+ &
439 & fac1* &
440 & boundary(ng)%tl_t_obc(j,k, &
441 & iwest,it1,linp,it)+ &
442 & fac2* &
443 & boundary(ng)%tl_t_obc(j,k, &
444 & iwest,it2,linp,it)
445 END DO
446 END DO
447 END IF
448
449 IF (lbc(ieast,istvar(it),ng)%acquire.and. &
450 & lobc(ieast,istvar(it),ng).and. &
451 & domain(ng)%Eastern_Edge(tile)) THEN
452 DO k=1,n(ng)
453 DO j=jstr,jend
454 boundary(ng)%t_east(j,k,it)=boundary(ng)%t_east(j,k,it)+ &
455 & fac1* &
456 & boundary(ng)%tl_t_obc(j,k, &
457 & ieast,it1,linp,it)+ &
458 & fac2* &
459 & boundary(ng)%tl_t_obc(j,k, &
460 & ieast,it2,linp,it)
461 END DO
462 END DO
463 END IF
464
465 IF (lbc(isouth,istvar(it),ng)%acquire.and. &
466 & lobc(isouth,istvar(it),ng).and. &
467 & domain(ng)%Southern_Edge(tile)) THEN
468 DO k=1,n(ng)
469 DO i=istr,iend
470 boundary(ng)%t_south(i,k,it)=boundary(ng)%t_south(i,k,it)+&
471 & fac1* &
472 & boundary(ng)%tl_t_obc(i,k, &
473 & isouth,it1,linp,it)+ &
474 & fac2* &
475 & boundary(ng)%tl_t_obc(i,k, &
476 & isouth,it2,linp,it)
477 END DO
478 END DO
479 END IF
480
481 IF (lbc(inorth,istvar(it),ng)%acquire.and. &
482 & lobc(inorth,istvar(it),ng).and. &
483 & domain(ng)%Northern_Edge(tile)) THEN
484 DO k=1,n(ng)
485 DO i=istr,iend
486 boundary(ng)%t_north(i,k,it)=boundary(ng)%t_north(i,k,it)+&
487 & fac1* &
488 & boundary(ng)%tl_t_obc(i,k, &
489 & inorth,it1,linp,it)+ &
490 & fac2* &
491 & boundary(ng)%tl_t_obc(i,k, &
492 & inorth,it2,linp,it)
493 END DO
494 END DO
495 END IF
496 END DO
497# endif
498!
499 RETURN
real(dp), dimension(:), allocatable dt
real(dp), dimension(:,:), allocatable obc_time
real(dp), dimension(:), allocatable time
integer, dimension(:), allocatable ntstart
integer, dimension(:), allocatable nbrec

References mod_boundary::boundary, mod_param::domain, mod_scalars::dt, mod_scalars::ieast, mod_scalars::iic, mod_scalars::inorth, mod_ncparam::isfsur, mod_scalars::isouth, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::iwest, mod_param::lbc, mod_scalars::lobc, mod_param::n, mod_scalars::nbrec, mod_scalars::nobc, mod_param::nt, mod_scalars::ntstart, mod_scalars::obc_time, and mod_scalars::time.

Referenced by obc_adjust().

Here is the caller graph for this function: