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