ROMS
Loading...
Searching...
No Matches
tl_obc_adjust.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#ifdef ADJUST_BOUNDARY
5!
6!git $Id$
7!================================================== Hernan G. Arango ===
8! Copyright (c) 2002-2025 The ROMS Group !
9! Licensed under a MIT/X style license !
10! See License_ROMS.md !
11!=======================================================================
12! !
13! This routine time-interpolates 4DVar tangent linear model open !
14! boundary increments. The increments can be constant (Nbrec=1) !
15! or time interpolated between snapshots (Nbrec>1). !
16! !
17! On Input: !
18! !
19! ng Nested grid number. !
20! tile Domain partition. !
21! Linp 4DVar increment time index to process. !
22! !
23!=======================================================================
24!
25 implicit none
26!
27 PRIVATE
28 PUBLIC :: tl_obc_adjust
29# ifdef SOLVE3D
30 PUBLIC :: tl_obc2d_adjust
31# endif
32!
33 CONTAINS
34!
35!***********************************************************************
36 SUBROUTINE tl_obc_adjust (ng, tile, Linp)
37!***********************************************************************
38!
39 USE mod_param
40!
41! Imported variable declarations.
42!
43 integer, intent(in) :: ng, tile, linp
44!
45! Local variable declarations.
46!
47 character (len=*), parameter :: myfile = &
48 & __FILE__
49!
50# include "tile.h"
51!
52# ifdef PROFILE
53 CALL wclock_on (ng, itlm, 7, __line__, myfile)
54# endif
55 CALL tl_obc_adjust_tile (ng, tile, &
56 & lbi, ubi, lbj, ubj, lbij, ubij, &
57 & imins, imaxs, jmins, jmaxs, &
58 & linp)
59# ifdef PROFILE
60 CALL wclock_off (ng, itlm, 7, __line__, myfile)
61# endif
62!
63 RETURN
64 END SUBROUTINE tl_obc_adjust
65!
66!***********************************************************************
67 SUBROUTINE tl_obc_adjust_tile (ng, tile, &
68 & LBi, UBi, LBj, UBj, LBij, UBij, &
69 & IminS, ImaxS, JminS, JmaxS, &
70 & Linp)
71!***********************************************************************
72!
73 USE mod_param
74 USE mod_boundary
75 USE mod_ncparam
76 USE mod_scalars
77!
78! Imported variable declarations.
79!
80 integer, intent(in) :: ng, tile
81 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
82 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
83 integer, intent(in) :: Linp
84!
85! Local variable declarations.
86!
87 integer :: i, it1, it2, j
88# ifdef SOLVE3D
89 integer :: it, k
90# endif
91 real(r8) :: fac, fac1, fac2
92
93# include "set_bounds.h"
94!
95!-----------------------------------------------------------------------
96! Adjust tangent linear open boundary fields with 4DVar increments.
97!-----------------------------------------------------------------------
98!
99! Set time records and interpolation factor, if any.
100!
101 IF (nbrec(ng).eq.1) THEN
102 it1=1
103 it2=1
104 fac1=1.0_r8
105 fac2=0.0_r8
106 ELSE
107# ifdef GENERIC_DSTART
108 it1=max(0,(iic(ng)-ntstart(ng))/nobc(ng))+1
109# else
110 it1=max(0,(iic(ng)-1)/nobc(ng))+1
111# endif
112 it2=min(it1+1,nbrec(ng))
113 fac1=obc_time(it2,ng)-(time(ng)+dt(ng))
114 fac2=(time(ng)+dt(ng))-obc_time(it1,ng)
115 fac=1.0_r8/(fac1+fac2)
116 fac1=fac*fac1
117 fac2=fac*fac2
118 END IF
119
120# ifdef SOLVE3D
121!
122! 3D U-momentum open boundaries.
123!
124 IF (tl_lbc(iwest,isuvel,ng)%acquire.and. &
125 & lobc(iwest,isuvel,ng).and. &
126 & domain(ng)%Western_Edge(tile)) THEN
127 DO k=1,n(ng)
128 DO j=jstr,jend
129 boundary(ng)%tl_u_west(j,k)=fac1* &
130 & boundary(ng)%tl_u_obc(j,k, &
131 & iwest,it1,linp)+ &
132 & fac2* &
133 & boundary(ng)%tl_u_obc(j,k, &
134 & iwest,it2,linp)
135 END DO
136 END DO
137 END IF
138
139 IF (tl_lbc(ieast,isuvel,ng)%acquire.and. &
140 & lobc(ieast,isuvel,ng).and. &
141 & domain(ng)%Eastern_Edge(tile)) THEN
142 DO k=1,n(ng)
143 DO j=jstr,jend
144 boundary(ng)%tl_u_east(j,k)=fac1* &
145 & boundary(ng)%tl_u_obc(j,k, &
146 & ieast,it1,linp)+ &
147 & fac2* &
148 & boundary(ng)%tl_u_obc(j,k, &
149 & ieast,it2,linp)
150 END DO
151 END DO
152 END IF
153
154 IF (tl_lbc(isouth,isuvel,ng)%acquire.and. &
155 & lobc(isouth,isuvel,ng).and. &
156 & domain(ng)%Southern_Edge(tile)) THEN
157 DO k=1,n(ng)
158 DO i=istru,iend
159 boundary(ng)%tl_u_south(i,k)=fac1* &
160 & boundary(ng)%tl_u_obc(i,k, &
161 & isouth,it1,linp)+ &
162 & fac2* &
163 & boundary(ng)%tl_u_obc(i,k, &
164 & isouth,it2,linp)
165 END DO
166 END DO
167 END IF
168
169 IF (tl_lbc(inorth,isuvel,ng)%acquire.and. &
170 & lobc(inorth,isuvel,ng).and. &
171 & domain(ng)%Northern_Edge(tile)) THEN
172 DO k=1,n(ng)
173 DO i=istru,iend
174 boundary(ng)%tl_u_north(i,k)=fac1* &
175 & boundary(ng)%tl_u_obc(i,k, &
176 & inorth,it1,linp)+ &
177 & fac2* &
178 & boundary(ng)%tl_u_obc(i,k, &
179 & inorth,it2,linp)
180 END DO
181 END DO
182 END IF
183!
184! 3D V-momentum open boundaries.
185!
186 IF (tl_lbc(iwest,isvvel,ng)%acquire.and. &
187 & lobc(iwest,isvvel,ng).and. &
188 & domain(ng)%Western_Edge(tile)) THEN
189 DO k=1,n(ng)
190 DO j=jstrv,jend
191 boundary(ng)%tl_v_west(j,k)=fac1* &
192 & boundary(ng)%tl_v_obc(j,k, &
193 & iwest,it1,linp)+ &
194 & fac2* &
195 & boundary(ng)%tl_v_obc(j,k, &
196 & iwest,it2,linp)
197 END DO
198 END DO
199 END IF
200
201 IF (tl_lbc(ieast,isvvel,ng)%acquire.and. &
202 & lobc(ieast,isvvel,ng).and. &
203 & domain(ng)%Eastern_Edge(tile)) THEN
204 DO k=1,n(ng)
205 DO j=jstrv,jend
206 boundary(ng)%tl_v_east(j,k)=fac1* &
207 & boundary(ng)%tl_v_obc(j,k, &
208 & ieast,it1,linp)+ &
209 & fac2* &
210 & boundary(ng)%tl_v_obc(j,k, &
211 & ieast,it2,linp)
212 END DO
213 END DO
214 END IF
215
216 IF (tl_lbc(isouth,isvvel,ng)%acquire.and. &
217 & lobc(isouth,isvvel,ng).and. &
218 & domain(ng)%Southern_Edge(tile)) THEN
219 DO k=1,n(ng)
220 DO i=istr,iend
221 boundary(ng)%tl_v_south(i,k)=fac1* &
222 & boundary(ng)%tl_v_obc(i,k, &
223 & isouth,it1,linp)+ &
224 & fac2* &
225 & boundary(ng)%tl_v_obc(i,k, &
226 & isouth,it2,linp)
227 END DO
228 END DO
229 END IF
230
231 IF (tl_lbc(inorth,isvvel,ng)%acquire.and. &
232 & lobc(inorth,isvvel,ng).and. &
233 & domain(ng)%Northern_Edge(tile)) THEN
234 DO k=1,n(ng)
235 DO i=istr,iend
236 boundary(ng)%tl_v_north(i,k)=fac1* &
237 & boundary(ng)%tl_v_obc(i,k, &
238 & inorth,it1,linp)+ &
239 & fac2* &
240 & boundary(ng)%tl_v_obc(i,k, &
241 & inorth,it2,linp)
242 END DO
243 END DO
244 END IF
245!
246! Tracers open boundaries.
247!
248 DO it=1,nt(ng)
249 IF (tl_lbc(iwest,istvar(it),ng)%acquire.and. &
250 & lobc(iwest,istvar(it),ng).and. &
251 & domain(ng)%Western_Edge(tile)) THEN
252 DO k=1,n(ng)
253 DO j=jstr,jend
254 boundary(ng)%tl_t_west(j,k,it)=fac1* &
255 & boundary(ng)%tl_t_obc(j, &
256 & k,iwest,it1,linp,it)+ &
257 & fac2* &
258 & boundary(ng)%tl_t_obc(j, &
259 & k,iwest,it2,linp,it)
260 END DO
261 END DO
262 END IF
263
264 IF (tl_lbc(ieast,istvar(it),ng)%acquire.and. &
265 & lobc(ieast,istvar(it),ng).and. &
266 & domain(ng)%Eastern_Edge(tile)) THEN
267 DO k=1,n(ng)
268 DO j=jstr,jend
269 boundary(ng)%tl_t_east(j,k,it)=fac1* &
270 & boundary(ng)%tl_t_obc(j, &
271 & k,ieast,it1,linp,it)+ &
272 & fac2* &
273 & boundary(ng)%tl_t_obc(j, &
274 & k,ieast,it2,linp,it)
275 END DO
276 END DO
277 END IF
278
279 IF (tl_lbc(isouth,istvar(it),ng)%acquire.and. &
280 & lobc(isouth,istvar(it),ng).and. &
281 & domain(ng)%Southern_Edge(tile)) THEN
282 DO k=1,n(ng)
283 DO i=istr,iend
284 boundary(ng)%tl_t_south(i,k,it)=fac1* &
285 & boundary(ng)%tl_t_obc(i, &
286 & k,isouth,it1,linp,it)+ &
287 & fac2* &
288 & boundary(ng)%tl_t_obc(i, &
289 & k,isouth,it2,linp,it)
290 END DO
291 END DO
292 END IF
293
294 IF (tl_lbc(inorth,istvar(it),ng)%acquire.and. &
295 & lobc(inorth,istvar(it),ng).and. &
296 & domain(ng)%Northern_Edge(tile)) THEN
297 DO k=1,n(ng)
298 DO i=istr,iend
299 boundary(ng)%tl_t_north(i,k,it)=fac1* &
300 & boundary(ng)%tl_t_obc(i, &
301 & k,inorth,it1,linp,it)+ &
302 & fac2* &
303 & boundary(ng)%tl_t_obc(i, &
304 & k,inorth,it2,linp,it)
305 END DO
306 END DO
307 END IF
308 END DO
309# endif
310!
311! Free-surface open boundaries.
312!
313 IF (tl_lbc(iwest,isfsur,ng)%acquire.and. &
314 & lobc(iwest,isfsur,ng).and. &
315 & domain(ng)%Western_Edge(tile)) THEN
316 DO j=jstr,jend
317 boundary(ng)%tl_zeta_west(j)=fac1* &
318 & boundary(ng)%tl_zeta_obc(j, &
319 & iwest,it1,linp)+ &
320 & fac2* &
321 & boundary(ng)%tl_zeta_obc(j, &
322 & iwest,it2,linp)
323 END DO
324 END IF
325
326 IF (tl_lbc(ieast,isfsur,ng)%acquire.and. &
327 & lobc(ieast,isfsur,ng).and. &
328 & domain(ng)%Eastern_Edge(tile)) THEN
329 DO j=jstr,jend
330 boundary(ng)%tl_zeta_east(j)=fac1* &
331 & boundary(ng)%tl_zeta_obc(j, &
332 & ieast,it1,linp)+ &
333 & fac2* &
334 & boundary(ng)%tl_zeta_obc(j, &
335 & ieast,it2,linp)
336 END DO
337 END IF
338
339 IF (tl_lbc(isouth,isfsur,ng)%acquire.and. &
340 & lobc(isouth,isfsur,ng).and. &
341 & domain(ng)%Southern_Edge(tile)) THEN
342 DO i=istr,iend
343 boundary(ng)%tl_zeta_south(i)=fac1* &
344 & boundary(ng)%tl_zeta_obc(i, &
345 & isouth,it1,linp)+ &
346 & fac2* &
347 & boundary(ng)%tl_zeta_obc(i, &
348 & isouth,it2,linp)
349 END DO
350 END IF
351
352 IF (tl_lbc(inorth,isfsur,ng)%acquire.and. &
353 & lobc(inorth,isfsur,ng).and. &
354 & domain(ng)%Northern_Edge(tile)) THEN
355 DO i=istr,iend
356 boundary(ng)%tl_zeta_north(i)=fac1* &
357 & boundary(ng)%tl_zeta_obc(i, &
358 & inorth,it1,linp)+ &
359 & fac2* &
360 & boundary(ng)%tl_zeta_obc(i, &
361 & inorth,it2,linp)
362 END DO
363 END IF
364
365# ifndef SOLVE3D
366!
367! 2D U-momentum open boundaries.
368!
369 IF (tl_lbc(iwest,isubar,ng)%acquire.and. &
370 & lobc(iwest,isubar,ng).and. &
371 & domain(ng)%Western_Edge(tile)) THEN
372 DO j=jstr,jend
373 boundary(ng)%tl_ubar_west(j)=fac1* &
374 & boundary(ng)%tl_ubar_obc(j, &
375 & iwest,it1,linp)+ &
376 & fac2* &
377 & boundary(ng)%tl_ubar_obc(j, &
378 & iwest,it2,linp)
379 END DO
380 END IF
381
382 IF (tl_lbc(ieast,isubar,ng)%acquire.and. &
383 & lobc(ieast,isubar,ng).and. &
384 & domain(ng)%Eastern_Edge(tile)) THEN
385 DO j=jstr,jend
386 boundary(ng)%tl_ubar_east(j)=fac1* &
387 & boundary(ng)%tl_ubar_obc(j, &
388 & ieast,it1,linp)+ &
389 & fac2* &
390 & boundary(ng)%tl_ubar_obc(j, &
391 & ieast,it2,linp)
392 END DO
393 END IF
394
395 IF (tl_lbc(isouth,isubar,ng)%acquire.and. &
396 & lobc(isouth,isubar,ng).and. &
397 & domain(ng)%Southern_Edge(tile)) THEN
398 DO i=istru,iend
399 boundary(ng)%tl_ubar_south(i)=fac1* &
400 & boundary(ng)%tl_ubar_obc(i, &
401 & isouth,it1,linp)+ &
402 & fac2* &
403 & boundary(ng)%tl_ubar_obc(i, &
404 & isouth,it2,linp)
405 END DO
406 END IF
407
408 IF (tl_lbc(inorth,isubar,ng)%acquire.and. &
409 & lobc(inorth,isubar,ng).and. &
410 & domain(ng)%Northern_Edge(tile)) THEN
411 DO i=istru,iend
412 boundary(ng)%tl_ubar_north(i)=fac1* &
413 & boundary(ng)%tl_ubar_obc(i, &
414 & inorth,it1,linp)+ &
415 & fac2* &
416 & boundary(ng)%tl_ubar_obc(i, &
417 & inorth,it2,linp)
418 END DO
419 END IF
420!
421! 2D V-momentum open boundaries.
422!
423 IF (tl_lbc(iwest,isvbar,ng)%acquire.and. &
424 & lobc(iwest,isvbar,ng).and. &
425 & domain(ng)%Western_Edge(tile)) THEN
426 DO j=jstrv,jend
427 boundary(ng)%tl_vbar_west(j)=fac1* &
428 & boundary(ng)%tl_vbar_obc(j, &
429 & iwest,it1,linp)+ &
430 & fac2* &
431 & boundary(ng)%tl_vbar_obc(j, &
432 & iwest,it2,linp)
433 END DO
434 END IF
435
436 IF (tl_lbc(ieast,isvbar,ng)%acquire.and. &
437 & lobc(ieast,isvbar,ng).and. &
438 & domain(ng)%Eastern_Edge(tile)) THEN
439 DO j=jstrv,jend
440 boundary(ng)%tl_vbar_east(j)=fac1* &
441 & boundary(ng)%tl_vbar_obc(j, &
442 & ieast,it1,linp)+ &
443 & fac2* &
444 & boundary(ng)%tl_vbar_obc(j, &
445 & ieast,it2,linp)
446 END DO
447 END IF
448
449 IF (tl_lbc(isouth,isvbar,ng)%acquire.and. &
450 & lobc(isouth,isvbar,ng).and. &
451 & domain(ng)%Southern_Edge(tile)) THEN
452 DO i=istr,iend
453 boundary(ng)%tl_vbar_south(i)=fac1* &
454 & boundary(ng)%tl_vbar_obc(i, &
455 & isouth,it1,linp)+ &
456 & fac2* &
457 & boundary(ng)%tl_vbar_obc(i, &
458 & isouth,it2,linp)
459 END DO
460 END IF
461
462 IF (tl_lbc(inorth,isvbar,ng)%acquire.and. &
463 & lobc(inorth,isvbar,ng).and. &
464 & domain(ng)%Northern_Edge(tile)) THEN
465 DO i=istr,iend
466 boundary(ng)%tl_vbar_north(i)=fac1* &
467 & boundary(ng)%tl_vbar_obc(i, &
468 & inorth,it1,linp)+ &
469 & fac2* &
470 & boundary(ng)%tl_vbar_obc(i, &
471 & inorth,it2,linp)
472 END DO
473 END IF
474# endif
475!
476 RETURN
477 END SUBROUTINE tl_obc_adjust_tile
478
479# ifdef SOLVE3D
480!
481!***********************************************************************
482 SUBROUTINE tl_obc2d_adjust (ng, tile, Linp)
483!***********************************************************************
484!
485 USE mod_param
486 USE mod_grid
487!
488! Imported variable declarations.
489!
490 integer, intent(in) :: ng, tile, linp
491!
492! Local variable declarations.
493!
494 character (len=*), parameter :: myfile = &
495 & __FILE__//", tl_obc2d_adjust"
496!
497# include "tile.h"
498!
499# ifdef PROFILE
500 CALL wclock_on (ng, itlm, 7, __line__, myfile)
501# endif
502 CALL tl_obc2d_adjust_tile (ng, tile, &
503 & lbi, ubi, lbj, ubj, lbij, ubij, &
504 & imins, imaxs, jmins, jmaxs, &
505 & linp, &
506# ifdef MASKING
507 & grid(ng) % umask, &
508 & grid(ng) % vmask, &
509# endif
510 & grid(ng) % Hz, &
511 & grid(ng) % Hz_bry, &
512 & grid(ng) % tl_Hz, &
513 & grid(ng) % tl_Hz_bry)
514# ifdef PROFILE
515 CALL wclock_off (ng, itlm, 7, __line__, myfile)
516# endif
517!
518 RETURN
519 END SUBROUTINE tl_obc2d_adjust
520!
521!***********************************************************************
522 SUBROUTINE tl_obc2d_adjust_tile (ng, tile, &
523 & LBi, UBi, LBj, UBj, LBij, UBij, &
524 & IminS, ImaxS, JminS, JmaxS, &
525 & Linp, &
526# ifdef MASKING
527 & umask, vmask, &
528# endif
529 & Hz, Hz_bry, &
530 & tl_Hz, tl_Hz_bry)
531!***********************************************************************
532!
533 USE mod_param
534 USE mod_boundary
535 USE mod_ncparam
536 USE mod_scalars
537!
538! Imported variable declarations.
539!
540 integer, intent(in) :: ng, tile
541 integer, intent(in) :: LBi, UBi, LBj, UBj, LBij, UBij
542 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
543 integer, intent(in) :: Linp
544!
545# ifdef ASSUMED_SHAPE
546# ifdef MASKING
547 real(r8), intent(in) :: umask(LBi:,LBj:)
548 real(r8), intent(in) :: vmask(LBi:,LBj:)
549# endif
550 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
551 real(r8), intent(in) :: Hz_bry(LBij:,:,:)
552 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
553 real(r8), intent(in) :: tl_Hz_bry(LBij:,:,:)
554
555# else
556
557# ifdef MASKING
558 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
559 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
560# endif
561 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
562 real(r8), intent(in) :: Hz_bry(LBij:UBij,N(ng),4)
563 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
564 real(r8), intent(in) :: tl_Hz_bry(LBij:UBij,N(ng),4)
565# endif
566!
567! Local variable declarations.
568!
569 integer :: i, it1, it2, j, k
570
571 real(r8) :: fac, fac1, fac2
572 real(r8) :: cff1, cff2, tl_cff1, tl_cff2
573
574 real(r8), dimension(0:N(ng)) :: CF
575 real(r8), dimension(0:N(ng)) :: DC
576
577 real(r8), dimension(0:N(ng)) :: tl_CF
578 real(r8), dimension(0:N(ng)) :: tl_DC
579
580# include "set_bounds.h"
581!
582!-----------------------------------------------------------------------
583! Adjust tangent linear open boundary fields with 4DVar increments.
584!-----------------------------------------------------------------------
585!
586! Set time records and interpolation factor, if any.
587!
588 IF (nbrec(ng).eq.1) THEN
589 it1=1
590 it2=1
591 fac1=1.0_r8
592 fac2=0.0_r8
593 ELSE
594# ifdef GENERIC_DSTART
595 it1=max(0,(iic(ng)-ntstart(ng))/nobc(ng))+1
596# else
597 it1=max(0,(iic(ng)-1)/nobc(ng))+1
598# endif
599 it2=min(it1+1,nbrec(ng))
600 fac1=obc_time(it2,ng)-(time(ng)+dt(ng))
601 fac2=(time(ng)+dt(ng))-obc_time(it1,ng)
602 fac=1.0_r8/(fac1+fac2)
603 fac1=fac*fac1
604 fac2=fac*fac2
605 END IF
606!
607! 2D U-momentum open boundaries: integrate 3D U-momentum at the open
608! boundaries.
609!
610 IF (tl_lbc(iwest,isubar,ng)%acquire.and. &
611 & lobc(iwest,isubar,ng).and. &
612 & domain(ng)%Western_Edge(tile)) THEN
613 i=bounds(ng)%edge(iwest,r2dvar)
614 DO j=jstr,jend
615 dc(0)=0.0_r8
616 tl_dc(0)=0.0_r8
617 cf(0)=0.0_r8
618 tl_cf(0)=0.0_r8
619 DO k=1,n(ng)
620 dc(k)=0.5_r8*(hz_bry(j,k,iwest)+ &
621 & hz(i+1,j,k))
622 tl_dc(k)=0.5_r8*(tl_hz_bry(j,k,iwest)+ &
623 & tl_hz(i+1,j,k))
624 dc(0)=dc(0)+dc(k)
625 tl_dc(0)=tl_dc(0)+tl_dc(k)
626 cf(0)=cf(0)+dc(k)*boundary(ng)%u_west(j,k)
627 tl_cf(0)=tl_cf(0)+ &
628 & tl_dc(k)*boundary(ng)%u_west(j,k)+ &
629 & dc(k)*boundary(ng)%tl_u_west(j,k)
630 END DO
631 cff1=1.0_r8/dc(0)
632 tl_cff1=-cff1*cff1*tl_dc(0)
633 tl_cff2=tl_cf(0)*cff1+cf(0)*tl_cff1
634# ifdef MASKING
635 tl_cff2=tl_cff2*umask(i,j)
636# endif
637 boundary(ng)%tl_ubar_west(j)=tl_cff2
638 END DO
639 END IF
640
641 IF (tl_lbc(ieast,isubar,ng)%acquire.and. &
642 & lobc(ieast,isubar,ng).and. &
643 & domain(ng)%Eastern_Edge(tile)) THEN
644 i=bounds(ng)%edge(ieast,r2dvar)
645 DO j=jstr,jend
646 dc(0)=0.0_r8
647 tl_dc(0)=0.0_r8
648 cf(0)=0.0_r8
649 tl_cf(0)=0.0_r8
650 DO k=1,n(ng)
651 dc(k)=0.5_r8*(hz(i-1,j,k)+ &
652 & hz_bry(j,k,ieast))
653 tl_dc(k)=0.5_r8*(tl_hz(i-1,j,k)+ &
654 & tl_hz_bry(j,k,ieast))
655 dc(0)=dc(0)+dc(k)
656 tl_dc(0)=tl_dc(0)+tl_dc(k)
657 cf(0)=cf(0)+dc(k)*boundary(ng)%u_east(j,k)
658 tl_cf(0)=tl_cf(0)+ &
659 & tl_dc(k)*boundary(ng)%u_east(j,k)+ &
660 & dc(k)*boundary(ng)%tl_u_east(j,k)
661 END DO
662 cff1=1.0_r8/dc(0)
663 tl_cff1=-cff1*cff1*tl_dc(0)
664 tl_cff2=tl_cf(0)*cff1+cf(0)*tl_cff1
665# ifdef MASKING
666 tl_cff2=tl_cff2*umask(i,j)
667# endif
668 boundary(ng)%tl_ubar_east(j)=tl_cff2
669 END DO
670 END IF
671
672 IF (tl_lbc(isouth,isubar,ng)%acquire.and. &
673 & lobc(isouth,isubar,ng).and. &
674 & domain(ng)%Southern_Edge(tile)) THEN
675 j=bounds(ng)%edge(isouth,r2dvar)
676 DO i=istr,iend
677 dc(0)=0.0_r8
678 tl_dc(0)=0.0_r8
679 cf(0)=0.0_r8
680 tl_cf(0)=0.0_r8
681 DO k=1,n(ng)
682 dc(k)=0.5_r8*(hz_bry(i-1,k,isouth)+ &
683 & hz_bry(i ,k,isouth))
684 tl_dc(k)=0.5_r8*(tl_hz_bry(i-1,k,isouth)+ &
685 & tl_hz_bry(i ,k,isouth))
686 dc(0)=dc(0)+dc(k)
687 tl_dc(0)=tl_dc(0)+tl_dc(k)
688 cf(0)=cf(0)+dc(k)*boundary(ng)%u_south(i,k)
689 tl_cf(0)=tl_cf(0)+ &
690 & tl_dc(k)*boundary(ng)%u_south(i,k)+ &
691 & dc(k)*boundary(ng)%tl_u_south(i,k)
692 END DO
693 cff1=1.0_r8/dc(0)
694 tl_cff1=-cff1*cff1*tl_dc(0)
695 tl_cff2=tl_cf(0)*cff1+cf(0)*tl_cff1
696# ifdef MASKING
697 tl_cff2=tl_cff2*umask(i,j)
698# endif
699 boundary(ng)%tl_ubar_south(i)=tl_cff2
700 END DO
701 END IF
702
703 IF (tl_lbc(inorth,isubar,ng)%acquire.and. &
704 & lobc(inorth,isubar,ng).and. &
705 & domain(ng)%Northern_Edge(tile)) THEN
706 j=bounds(ng)%edge(inorth,r2dvar)
707 DO i=istr,iend
708 dc(0)=0.0_r8
709 tl_dc(0)=0.0_r8
710 cf(0)=0.0_r8
711 tl_cf(0)=0.0_r8
712 DO k=1,n(ng)
713 dc(k)=0.5_r8*(hz_bry(i-1,k,inorth)+ &
714 & hz_bry(i ,k,inorth))
715 tl_dc(k)=0.5_r8*(tl_hz_bry(i-1,k,inorth)+ &
716 & tl_hz_bry(i ,k,inorth))
717 dc(0)=dc(0)+dc(k)
718 tl_dc(0)=tl_dc(0)+tl_dc(k)
719 cf(0)=cf(0)+dc(k)*boundary(ng)%u_north(i,k)
720 tl_cf(0)=tl_cf(0)+ &
721 & tl_dc(k)*boundary(ng)%u_north(i,k)+ &
722 & dc(k)*boundary(ng)%tl_u_north(i,k)
723 END DO
724 cff1=1.0_r8/dc(0)
725 tl_cff1=-cff1*cff1*tl_dc(0)
726 tl_cff2=tl_cf(0)*cff1+cf(0)*tl_cff1
727# ifdef MASKING
728 tl_cff2=tl_cff2*umask(i,j)
729# endif
730 boundary(ng)%tl_ubar_north(i)=tl_cff2
731 END DO
732 END IF
733!
734! 2D V-momentum open boundaries: integrate 3D V-momentum at the open
735! boundaries.
736!
737 IF (tl_lbc(iwest,isvbar,ng)%acquire.and. &
738 & lobc(iwest,isvbar,ng).and. &
739 & domain(ng)%Western_Edge(tile)) THEN
740 i=bounds(ng)%edge(iwest,r2dvar)
741 DO j=jstrv,jend
742 dc(0)=0.0_r8
743 tl_dc(0)=0.0_r8
744 cf(0)=0.0_r8
745 tl_cf(0)=0.0_r8
746 DO k=1,n(ng)
747 dc(k)=0.5_r8*(hz_bry(j-1,k,iwest)+ &
748 & hz_bry(j ,k,iwest))
749 tl_dc(k)=0.5_r8*(tl_hz_bry(j-1,k,iwest)+ &
750 & tl_hz_bry(j ,k,iwest))
751 dc(0)=dc(0)+dc(k)
752 tl_dc(0)=tl_dc(0)+tl_dc(k)
753 cf(0)=cf(0)+dc(k)*boundary(ng)%v_west(j,k)
754 tl_cf(0)=tl_cf(0)+ &
755 & tl_dc(k)*boundary(ng)%v_west(j,k)+ &
756 & dc(k)*boundary(ng)%tl_v_west(j,k)
757 END DO
758 cff1=1.0_r8/dc(0)
759 tl_cff1=-cff1*cff1*tl_dc(0)
760 tl_cff2=tl_cf(0)*cff1+cf(0)*tl_cff1
761# ifdef MASKING
762 tl_cff2=tl_cff2*vmask(i,j)
763# endif
764 boundary(ng)%tl_vbar_west(j)=tl_cff2
765 END DO
766 END IF
767
768 IF (tl_lbc(ieast,isvbar,ng)%acquire.and. &
769 & lobc(ieast,isvbar,ng).and. &
770 & domain(ng)%Eastern_Edge(tile)) THEN
771 i=bounds(ng)%edge(ieast,r2dvar)
772 DO j=jstrv,jend
773 dc(0)=0.0_r8
774 tl_dc(0)=0.0_r8
775 cf(0)=0.0_r8
776 tl_cf(0)=0.0_r8
777 DO k=1,n(ng)
778 dc(k)=0.5_r8*(hz_bry(j-1,k,ieast)+ &
779 & hz_bry(j ,k,ieast))
780 tl_dc(k)=0.5_r8*(tl_hz_bry(j-1,k,ieast)+ &
781 & tl_hz_bry(j ,k,ieast))
782 dc(0)=dc(0)+dc(k)
783 tl_dc(0)=tl_dc(0)+tl_dc(k)
784 cf(0)=cf(0)+dc(k)*boundary(ng)%v_east(j,k)
785 tl_cf(0)=tl_cf(0)+ &
786 & tl_dc(k)*boundary(ng)%v_east(j,k)+ &
787 & dc(k)*boundary(ng)%tl_v_east(j,k)
788 END DO
789 cff1=1.0_r8/dc(0)
790 tl_cff1=-cff1*cff1*tl_dc(0)
791 tl_cff2=tl_cf(0)*cff1+cf(0)*tl_cff1
792# ifdef MASKING
793 tl_cff2=tl_cff2*vmask(i,j)
794# endif
795 boundary(ng)%tl_vbar_east(j)=tl_cff2
796 END DO
797 END IF
798
799 IF (tl_lbc(isouth,isvbar,ng)%acquire.and. &
800 & lobc(isouth,isvbar,ng).and. &
801 & domain(ng)%Southern_Edge(tile)) THEN
802 j=bounds(ng)%edge(isouth,r2dvar)
803 DO i=istr,iend
804 dc(0)=0.0_r8
805 tl_dc(0)=0.0_r8
806 cf(0)=0.0_r8
807 tl_cf(0)=0.0_r8
808 DO k=1,n(ng)
809 dc(k)=0.5_r8*(hz_bry(i,k,isouth)+ &
810 & hz(i+1,j,k))
811 tl_dc(k)=0.5_r8*(tl_hz_bry(i,k,isouth)+ &
812 & tl_hz(i+1,j,k))
813 dc(0)=dc(0)+dc(k)
814 tl_dc(0)=tl_dc(0)+tl_dc(k)
815 cf(0)=cf(0)+dc(k)*boundary(ng)%v_south(i,k)
816 tl_cf(0)=tl_cf(0)+ &
817 & tl_dc(k)*boundary(ng)%v_south(i,k)+ &
818 & dc(k)*boundary(ng)%tl_v_south(i,k)
819 END DO
820 cff1=1.0_r8/dc(0)
821 tl_cff1=-cff1*cff1*tl_dc(0)
822 tl_cff2=tl_cf(0)*cff1+cf(0)*tl_cff1
823# ifdef MASKING
824 tl_cff2=tl_cff2*vmask(i,j)
825# endif
826 boundary(ng)%tl_vbar_south(i)=tl_cff2
827 END DO
828 END IF
829
830 IF (tl_lbc(inorth,isvbar,ng)%acquire.and. &
831 & lobc(inorth,isvbar,ng).and. &
832 & domain(ng)%Northern_Edge(tile)) THEN
833 j=bounds(ng)%edge(inorth,r2dvar)
834 DO i=istr,iend
835 dc(0)=0.0_r8
836 tl_dc(0)=0.0_r8
837 cf(0)=0.0_r8
838 tl_cf(0)=0.0_r8
839 DO k=1,n(ng)
840 dc(k)=0.5_r8*(hz(i,j-1,k)+ &
841 & hz_bry(i,k,inorth))
842 tl_dc(k)=0.5_r8*(tl_hz(i,j-1,k)+ &
843 & tl_hz_bry(i,k,inorth))
844 dc(0)=dc(0)+dc(k)
845 tl_dc(0)=tl_dc(0)+tl_dc(k)
846 cf(0)=cf(0)+dc(k)*boundary(ng)%v_north(i,k)
847 tl_cf(0)=tl_cf(0)+ &
848 & tl_dc(k)*boundary(ng)%v_north(i,k)+ &
849 & dc(k)*boundary(ng)%tl_v_north(i,k)
850 END DO
851 cff1=1.0_r8/dc(0)
852 tl_cff1=-cff1*cff1*tl_dc(0)
853 tl_cff2=tl_cf(0)*cff1+cf(0)*tl_cff1
854# ifdef MASKING
855 tl_cff2=tl_cff2*vmask(i,j)
856# endif
857 boundary(ng)%tl_vbar_north(i)=tl_cff2
858 END DO
859 END IF
860!
861 RETURN
862 END SUBROUTINE tl_obc2d_adjust_tile
863# endif
864#endif
865 END MODULE tl_obc_adjust_mod
type(t_boundary), dimension(:), allocatable boundary
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
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_bounds), dimension(:), allocatable bounds
Definition mod_param.F:232
type(t_lbc), dimension(:,:,:), allocatable tl_lbc
Definition mod_param.F:379
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, parameter itlm
Definition mod_param.F:663
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer, parameter r2dvar
Definition mod_param.F:717
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
integer, dimension(:), allocatable nobc
logical, dimension(:,:,:), allocatable lobc
integer, parameter iwest
real(dp), dimension(:,:), allocatable obc_time
integer, parameter isouth
integer, parameter ieast
real(dp), dimension(:), allocatable time
integer, dimension(:), allocatable ntstart
integer, parameter inorth
integer, dimension(:), allocatable nbrec
subroutine, public tl_obc2d_adjust(ng, tile, linp)
subroutine, public tl_obc_adjust(ng, tile, linp)
subroutine tl_obc2d_adjust_tile(ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, imins, imaxs, jmins, jmaxs, linp, umask, vmask, hz, hz_bry, tl_hz, tl_hz_bry)
subroutine tl_obc_adjust_tile(ng, tile, lbi, ubi, lbj, ubj, lbij, ubij, imins, imaxs, jmins, jmaxs, linp)
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