ROMS
Loading...
Searching...
No Matches
t3dbc_im.F
Go to the documentation of this file.
1#include "cppdefs.h"
2 MODULE t3dbc_mod
3#ifdef SOLVE3D
4!
5!git $Id$
6!=======================================================================
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md Hernan G. Arango !
10!========================================== Alexander F. Shchepetkin ===
11! !
12! This subroutine sets lateral boundary conditions for the ITRC-th !
13! tracer field. !
14! !
15!=======================================================================
16!
17 implicit none
18!
19 PRIVATE
20 PUBLIC :: t3dbc_tile
21!
22 CONTAINS
23!
24!***********************************************************************
25 SUBROUTINE t3dbc (ng, tile, nout, itrc, ic)
26!***********************************************************************
27!
28 USE mod_param
29 USE mod_ocean
30 USE mod_stepping
31!
32! Imported variable declarations.
33!
34 integer, intent(in) :: ng, tile, nout, itrc, ic
35!
36! Local variable declarations.
37!
38# include "tile.h"
39!
40 CALL t3dbc_tile (ng, tile, itrc, ic, &
41 & lbi, ubi, lbj, ubj, n(ng), nt(ng), &
42 & imins, imaxs, jmins, jmaxs, &
43 & nstp(ng), nout, &
44 & ocean(ng)% t)
45
46 RETURN
47 END SUBROUTINE t3dbc
48!
49!***********************************************************************
50 SUBROUTINE t3dbc_tile (ng, tile, itrc, ic, &
51 & LBi, UBi, LBj, UBj, UBk, UBt, &
52 & IminS, ImaxS, JminS, JmaxS, &
53 & nstp, nout, &
54 & t)
55!***********************************************************************
56!
57 USE mod_param
58 USE mod_boundary
59 USE mod_clima
60 USE mod_grid
61 USE mod_ncparam
62 USE mod_scalars
63!
64! Imported variable declarations.
65!
66 integer, intent(in) :: ng, tile, itrc, ic
67 integer, intent(in) :: lbi, ubi, lbj, ubj, ubk, ubt
68 integer, intent(in) :: imins, imaxs, jmins, jmaxs
69 integer, intent(in) :: nstp, nout
70!
71# ifdef ASSUMED_SHAPE
72 real(r8), intent(inout) :: t(lbi:,lbj:,:,:,:)
73# else
74 real(r8), intent(inout) :: t(lbi:ubi,lbj:ubj,ubk,3,ubt)
75# endif
76!
77! Local variable declarations.
78!
79 integer :: i, j, k
80
81 real(r8), parameter :: eps =1.0e-20_r8
82
83 real(r8) :: ce, cx, cff, dtde, dtdt, dtdx
84 real(r8) :: obc_in, obc_out, tau
85
86 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: grad
87
88# include "set_bounds.h"
89!
90!-----------------------------------------------------------------------
91! Lateral boundary conditions at the western edge.
92!-----------------------------------------------------------------------
93!
94 IF (domain(ng)%Western_Edge(tile)) THEN
95!
96! Western edge, implicit upstream radiation condition.
97!
98 IF (lbc(iwest,istvar(itrc),ng)%radiation) THEN
99 DO k=1,n(ng)
100 DO j=jstr,jend+1
101 grad(istr-1,j)=t(istr-1,j ,k,nstp,itrc)- &
102 & t(istr-1,j-1,k,nstp,itrc)
103# ifdef MASKING
104 grad(istr-1,j)=grad(istr-1,j)* &
105 & grid(ng)%vmask(istr-1,j)
106# endif
107 grad(istr ,j)=t(istr ,j ,k,nstp,itrc)- &
108 & t(istr ,j-1,k,nstp,itrc)
109# ifdef MASKING
110 grad(istr ,j)=grad(istr ,j)* &
111 & grid(ng)%vmask(istr ,j)
112# endif
113 END DO
114 DO j=jstr,jend
115 IF (lbc_apply(ng)%west(j)) THEN
116 dtdt=t(istr,j,k,nstp,itrc)-t(istr ,j,k,nout,itrc)
117 dtdx=t(istr,j,k,nout,itrc)-t(istr+1,j,k,nout,itrc)
118
119 IF (lbc(iwest,istvar(itrc),ng)%nudging) THEN
120 IF (lnudgetclm(itrc,ng)) THEN
121 obc_out=clima(ng)%Tnudgcof(istr-1,j,k,ic)
122 obc_in =obcfac(ng)*obc_out
123 ELSE
124 obc_out=tobc_out(itrc,ng,iwest)
125 obc_in =tobc_in(itrc,ng,iwest)
126 END IF
127 IF ((dtdt*dtdx).lt.0.0_r8) THEN
128 tau=obc_in
129 ELSE
130 tau=obc_out
131 END IF
132 tau=tau*dt(ng)
133 END IF
134
135 IF ((dtdt*dtdx).lt.0.0_r8) dtdt=0.0_r8
136 IF ((dtdt*(grad(istr,j )+ &
137 & grad(istr,j+1))).gt.0.0_r8) THEN
138 dtde=grad(istr,j )
139 ELSE
140 dtde=grad(istr,j+1)
141 END IF
142 cff=max(dtdx*dtdx+dtde*dtde,eps)
143 cx=dtdt*dtdx
144# ifdef RADIATION_2D
145 ce=min(cff,max(dtdt*dtde,-cff))
146# else
147 ce=0.0_r8
148# endif
149# if defined CELERITY_WRITE && defined FORWARD_WRITE
150 boundary(ng)%t_west_Cx(j,k,itrc)=cx
151 boundary(ng)%t_west_Ce(j,k,itrc)=ce
152 boundary(ng)%t_west_C2(j,k,itrc)=cff
153# endif
154 t(istr-1,j,k,nout,itrc)=(cff*t(istr-1,j,k,nstp,itrc)+ &
155 & cx *t(istr ,j,k,nout,itrc)- &
156 & max(ce,0.0_r8)* &
157 & grad(istr-1,j )- &
158 & min(ce,0.0_r8)* &
159 & grad(istr-1,j+1))/ &
160 & (cff+cx)
161
162 IF (lbc(iwest,istvar(itrc),ng)%nudging) THEN
163 t(istr-1,j,k,nout,itrc)=t(istr-1,j,k,nout,itrc)+ &
164 & tau* &
165 & (boundary(ng)%t_west(j,k, &
166 & itrc)- &
167 & t(istr-1,j,k,nstp,itrc))
168 END IF
169# ifdef MASKING
170 t(istr-1,j,k,nout,itrc)=t(istr-1,j,k,nout,itrc)* &
171 & grid(ng)%rmask(istr-1,j)
172# endif
173 END IF
174 END DO
175 END DO
176!
177! Western edge, clamped boundary condition.
178!
179 ELSE IF (lbc(iwest,istvar(itrc),ng)%clamped) THEN
180 DO k=1,n(ng)
181 DO j=jstr,jend
182 IF (lbc_apply(ng)%west(j)) THEN
183 t(istr-1,j,k,nout,itrc)=boundary(ng)%t_west(j,k,itrc)
184# ifdef MASKING
185 t(istr-1,j,k,nout,itrc)=t(istr-1,j,k,nout,itrc)* &
186 & grid(ng)%rmask(istr-1,j)
187# endif
188 END IF
189 END DO
190 END DO
191!
192! Western edge, gradient boundary condition.
193!
194 ELSE IF (lbc(iwest,istvar(itrc),ng)%gradient) THEN
195 DO k=1,n(ng)
196 DO j=jstr,jend
197 IF (lbc_apply(ng)%west(j)) THEN
198 t(istr-1,j,k,nout,itrc)=t(istr,j,k,nout,itrc)
199# ifdef MASKING
200 t(istr-1,j,k,nout,itrc)=t(istr-1,j,k,nout,itrc)* &
201 & grid(ng)%rmask(istr-1,j)
202# endif
203 END IF
204 END DO
205 END DO
206!
207! Western edge, closed boundary condition.
208!
209 ELSE IF (lbc(iwest,istvar(itrc),ng)%closed) THEN
210 DO k=1,n(ng)
211 DO j=jstr,jend
212 IF (lbc_apply(ng)%west(j)) THEN
213 t(istr-1,j,k,nout,itrc)=t(istr,j,k,nout,itrc)
214# ifdef MASKING
215 t(istr-1,j,k,nout,itrc)=t(istr-1,j,k,nout,itrc)* &
216 & grid(ng)%rmask(istr-1,j)
217# endif
218 END IF
219 END DO
220 END DO
221 END IF
222 END IF
223!
224!-----------------------------------------------------------------------
225! Lateral boundary conditions at the eastern edge.
226!-----------------------------------------------------------------------
227!
228 IF (domain(ng)%Eastern_Edge(tile)) THEN
229!
230! Eastern edge, implicit upstream radiation condition.
231!
232 IF (lbc(ieast,istvar(itrc),ng)%radiation) THEN
233 DO k=1,n(ng)
234 DO j=jstr,jend+1
235 grad(iend ,j)=t(iend ,j ,k,nstp,itrc)- &
236 & t(iend ,j-1,k,nstp,itrc)
237# ifdef MASKING
238 grad(iend ,j)=grad(iend ,j)* &
239 & grid(ng)%vmask(iend ,j)
240# endif
241 grad(iend+1,j)=t(iend+1,j ,k,nstp,itrc)- &
242 & t(iend+1,j-1,k,nstp,itrc)
243# ifdef MASKING
244 grad(iend+1,j)=grad(iend+1,j)* &
245 & grid(ng)%vmask(iend+1,j)
246# endif
247 END DO
248 DO j=jstr,jend
249 IF (lbc_apply(ng)%east(j)) THEN
250 dtdt=t(iend,j,k,nstp,itrc)-t(iend ,j,k,nout,itrc)
251 dtdx=t(iend,j,k,nout,itrc)-t(iend-1,j,k,nout,itrc)
252
253 IF (lbc(ieast,istvar(itrc),ng)%nudging) THEN
254 IF (lnudgetclm(itrc,ng)) THEN
255 obc_out=clima(ng)%Tnudgcof(iend+1,j,k,ic)
256 obc_in =obcfac(ng)*obc_out
257 ELSE
258 obc_out=tobc_out(itrc,ng,ieast)
259 obc_in =tobc_in(itrc,ng,ieast)
260 END IF
261 IF ((dtdt*dtdx).lt.0.0_r8) THEN
262 tau=obc_in
263 ELSE
264 tau=obc_out
265 END IF
266 tau=tau*dt(ng)
267 END IF
268
269 IF ((dtdt*dtdx).lt.0.0_r8) dtdt=0.0_r8
270 IF ((dtdt*(grad(iend,j )+ &
271 & grad(iend,j+1))).gt.0.0_r8) THEN
272 dtde=grad(iend,j )
273 ELSE
274 dtde=grad(iend,j+1)
275 END IF
276 cff=max(dtdx*dtdx+dtde*dtde,eps)
277 cx=dtdt*dtdx
278# ifdef RADIATION_2D
279 ce=min(cff,max(dtdt*dtde,-cff))
280# else
281 ce=0.0_r8
282# endif
283# if defined CELERITY_WRITE && defined FORWARD_WRITE
284 boundary(ng)%t_east_Cx(j,k,itrc)=cx
285 boundary(ng)%t_east_Ce(j,k,itrc)=ce
286 boundary(ng)%t_east_C2(j,k,itrc)=cff
287# endif
288 t(iend+1,j,k,nout,itrc)=(cff*t(iend+1,j,k,nstp,itrc)+ &
289 & cx *t(iend ,j,k,nout,itrc)- &
290 & max(ce,0.0_r8)* &
291 & grad(iend+1,j )- &
292 & min(ce,0.0_r8)* &
293 & grad(iend+1,j+1))/ &
294 & (cff+cx)
295
296 IF (lbc(ieast,istvar(itrc),ng)%nudging) THEN
297 t(iend+1,j,k,nout,itrc)=t(iend+1,j,k,nout,itrc)+ &
298 & tau* &
299 & (boundary(ng)%t_east(j,k, &
300 & itrc)- &
301 & t(iend+1,j,k,nstp,itrc))
302 END IF
303# ifdef MASKING
304 t(iend+1,j,k,nout,itrc)=t(iend+1,j,k,nout,itrc)* &
305 & grid(ng)%rmask(iend+1,j)
306# endif
307 END IF
308 END DO
309 END DO
310!
311! Eastern edge, clamped boundary condition.
312!
313 ELSE IF (lbc(ieast,istvar(itrc),ng)%clamped) THEN
314 DO k=1,n(ng)
315 DO j=jstr,jend
316 IF (lbc_apply(ng)%east(j)) THEN
317 t(iend+1,j,k,nout,itrc)=boundary(ng)%t_east(j,k,itrc)
318# ifdef MASKING
319 t(iend+1,j,k,nout,itrc)=t(iend+1,j,k,nout,itrc)* &
320 & grid(ng)%rmask(iend+1,j)
321# endif
322 END IF
323 END DO
324 END DO
325!
326! Eastern edge, gradient boundary condition.
327!
328 ELSE IF (lbc(ieast,istvar(itrc),ng)%gradient) THEN
329 DO k=1,n(ng)
330 DO j=jstr,jend
331 IF (lbc_apply(ng)%east(j)) THEN
332 t(iend+1,j,k,nout,itrc)=t(iend,j,k,nout,itrc)
333# ifdef MASKING
334 t(iend+1,j,k,nout,itrc)=t(iend+1,j,k,nout,itrc)* &
335 & grid(ng)%rmask(iend+1,j)
336# endif
337 END IF
338 END DO
339 END DO
340!
341! Eastern edge, closed boundary condition.
342!
343 ELSE IF (lbc(ieast,istvar(itrc),ng)%closed) THEN
344 DO k=1,n(ng)
345 DO j=jstr,jend
346 IF (lbc_apply(ng)%east(j)) THEN
347 t(iend+1,j,k,nout,itrc)=t(iend,j,k,nout,itrc)
348# ifdef MASKING
349 t(iend+1,j,k,nout,itrc)=t(iend+1,j,k,nout,itrc)* &
350 & grid(ng)%rmask(iend+1,j)
351# endif
352 END IF
353 END DO
354 END DO
355 END IF
356 END IF
357!
358!-----------------------------------------------------------------------
359! Lateral boundary conditions at the southern edge.
360!-----------------------------------------------------------------------
361!
362 IF (domain(ng)%Southern_Edge(tile)) THEN
363!
364! Southern edge, implicit upstream radiation condition.
365!
366 IF (lbc(isouth,istvar(itrc),ng)%radiation) THEN
367 DO k=1,n(ng)
368 DO i=istr,iend+1
369 grad(i,jstr )=t(i ,jstr ,k,nstp,itrc)- &
370 & t(i-1,jstr ,k,nstp,itrc)
371# ifdef MASKING
372 grad(i,jstr )=grad(i,jstr )* &
373 & grid(ng)%umask(i,jstr )
374# endif
375 grad(i,jstr-1)=t(i ,jstr-1,k,nstp,itrc)- &
376 & t(i-1,jstr-1,k,nstp,itrc)
377# ifdef MASKING
378 grad(i,jstr-1)=grad(i,jstr-1)* &
379 & grid(ng)%umask(i,jstr-1)
380# endif
381 END DO
382 DO i=istr,iend
383 IF (lbc_apply(ng)%south(i)) THEN
384 dtdt=t(i,jstr,k,nstp,itrc)-t(i,jstr ,k,nout,itrc)
385 dtde=t(i,jstr,k,nout,itrc)-t(i,jstr+1,k,nout,itrc)
386
387 IF (lbc(isouth,istvar(itrc),ng)%nudging) THEN
388 IF (lnudgetclm(itrc,ng)) THEN
389 obc_out=clima(ng)%Tnudgcof(i,jstr-1,k,ic)
390 obc_in =obcfac(ng)*obc_out
391 ELSE
392 obc_out=tobc_out(itrc,ng,isouth)
393 obc_in =tobc_in(itrc,ng,isouth)
394 END IF
395 IF ((dtdt*dtde).lt.0.0_r8) THEN
396 tau=obc_in
397 ELSE
398 tau=obc_out
399 END IF
400 tau=tau*dt(ng)
401 END IF
402
403 IF ((dtdt*dtde).lt.0.0_r8) dtdt=0.0_r8
404 IF ((dtdt*(grad(i ,jstr)+ &
405 & grad(i+1,jstr))).gt.0.0_r8) THEN
406 dtdx=grad(i ,jstr)
407 ELSE
408 dtdx=grad(i+1,jstr)
409 END IF
410 cff=max(dtdx*dtdx+dtde*dtde,eps)
411# ifdef RADIATION_2D
412 cx=min(cff,max(dtdt*dtdx,-cff))
413# else
414 cx=0.0_r8
415# endif
416 ce=dtdt*dtde
417# if defined CELERITY_WRITE && defined FORWARD_WRITE
418 boundary(ng)%t_south_Cx(i,k,itrc)=cx
419 boundary(ng)%t_south_Ce(i,k,itrc)=ce
420 boundary(ng)%t_south_C2(i,k,itrc)=cff
421# endif
422 t(i,jstr-1,k,nout,itrc)=(cff*t(i,jstr-1,k,nstp,itrc)+ &
423 & ce *t(i,jstr ,k,nout,itrc )- &
424 & max(cx,0.0_r8)* &
425 & grad(i ,jstr-1)- &
426 & min(cx,0.0_r8)* &
427 & grad(i+1,jstr-1))/ &
428 & (cff+ce)
429
430 IF (lbc(isouth,istvar(itrc),ng)%nudging) THEN
431 t(i,jstr-1,k,nout,itrc)=t(i,jstr-1,k,nout,itrc)+ &
432 & tau* &
433 & (boundary(ng)%t_south(i,k, &
434 & itrc)- &
435 & t(i,jstr-1,k,nstp,itrc))
436 END IF
437# ifdef MASKING
438 t(i,jstr-1,k,nout,itrc)=t(i,jstr-1,k,nout,itrc)* &
439 & grid(ng)%rmask(i,jstr-1)
440# endif
441 END IF
442 END DO
443 END DO
444!
445! Southern edge, clamped boundary condition.
446!
447 ELSE IF (lbc(isouth,istvar(itrc),ng)%clamped) THEN
448 DO k=1,n(ng)
449 DO i=istr,iend
450 IF (lbc_apply(ng)%south(i)) THEN
451 t(i,jstr-1,k,nout,itrc)=boundary(ng)%t_south(i,k,itrc)
452# ifdef MASKING
453 t(i,jstr-1,k,nout,itrc)=t(i,jstr-1,k,nout,itrc)* &
454 & grid(ng)%rmask(i,jstr-1)
455# endif
456 END IF
457 END DO
458 END DO
459!
460! Southern edge, gradient boundary condition.
461!
462 ELSE IF (lbc(isouth,istvar(itrc),ng)%gradient) THEN
463 DO k=1,n(ng)
464 DO i=istr,iend
465 IF (lbc_apply(ng)%south(i)) THEN
466 t(i,jstr-1,k,nout,itrc)=t(i,jstr,k,nout,itrc)
467# ifdef MASKING
468 t(i,jstr-1,k,nout,itrc)=t(i,jstr-1,k,nout,itrc)* &
469 & grid(ng)%rmask(i,jstr-1)
470# endif
471 END IF
472 END DO
473 END DO
474!
475! Southern edge, closed boundary condition.
476!
477 ELSE IF (lbc(isouth,istvar(itrc),ng)%closed) THEN
478 DO k=1,n(ng)
479 DO i=istr,iend
480 IF (lbc_apply(ng)%south(i)) THEN
481 t(i,jstr-1,k,nout,itrc)=t(i,jstr,k,nout,itrc)
482# ifdef MASKING
483 t(i,jstr-1,k,nout,itrc)=t(i,jstr-1,k,nout,itrc)* &
484 & grid(ng)%rmask(i,jstr-1)
485# endif
486 END IF
487 END DO
488 END DO
489 END IF
490 END IF
491!
492!-----------------------------------------------------------------------
493! Lateral boundary conditions at the northern edge.
494!-----------------------------------------------------------------------
495!
496 IF (domain(ng)%Northern_Edge(tile)) THEN
497!
498! Northern edge, implicit upstream radiation condition.
499!
500 IF (lbc(inorth,istvar(itrc),ng)%radiation) THEN
501 DO k=1,n(ng)
502 DO i=istr,iend+1
503 grad(i,jend )=t(i ,jend ,k,nstp,itrc)- &
504 & t(i-1,jend ,k,nstp,itrc)
505# ifdef MASKING
506 grad(i,jend )=grad(i,jend )* &
507 & grid(ng)%umask(i,jend )
508# endif
509 grad(i,jend+1)=t(i ,jend+1,k,nstp,itrc)- &
510 & t(i-1,jend+1,k,nstp,itrc)
511# ifdef MASKING
512 grad(i,jend+1)=grad(i,jend+1)* &
513 & grid(ng)%umask(i,jend+1)
514# endif
515 END DO
516 DO i=istr,iend
517 IF (lbc_apply(ng)%north(i)) THEN
518 dtdt=t(i,jend,k,nstp,itrc)-t(i,jend ,k,nout,itrc)
519 dtde=t(i,jend,k,nout,itrc)-t(i,jend-1,k,nout,itrc)
520
521 IF (lbc(inorth,istvar(itrc),ng)%nudging) THEN
522 IF (lnudgetclm(itrc,ng)) THEN
523 obc_out=clima(ng)%Tnudgcof(i,jend+1,k,ic)
524 obc_in =obcfac(ng)*obc_out
525 ELSE
526 obc_out=tobc_out(itrc,ng,inorth)
527 obc_in =tobc_in(itrc,ng,inorth)
528 END IF
529 IF ((dtdt*dtde).lt.0.0_r8) THEN
530 tau=obc_in
531 ELSE
532 tau=obc_out
533 END IF
534 tau=tau*dt(ng)
535 END IF
536
537 IF ((dtdt*dtde).lt.0.0_r8) dtdt=0.0_r8
538 IF ((dtdt*(grad(i ,jend)+ &
539 & grad(i+1,jend))).gt.0.0_r8) THEN
540 dtdx=grad(i ,jend)
541 ELSE
542 dtdx=grad(i+1,jend)
543 END IF
544 cff=max(dtdx*dtdx+dtde*dtde,eps)
545# ifdef RADIATION_2D
546 cx=min(cff,max(dtdt*dtdx,-cff))
547# else
548 cx=0.0_r8
549# endif
550 ce=dtdt*dtde
551# if defined CELERITY_WRITE && defined FORWARD_WRITE
552 boundary(ng)%t_north_Cx(i,k,itrc)=cx
553 boundary(ng)%t_north_Ce(i,k,itrc)=ce
554 boundary(ng)%t_north_C2(i,k,itrc)=cff
555# endif
556 t(i,jend+1,k,nout,itrc)=(cff*t(i,jend+1,k,nstp,itrc)+ &
557 & ce *t(i,jend ,k,nout,itrc)- &
558 & max(cx,0.0_r8)* &
559 & grad(i ,jend+1)- &
560 & min(cx,0.0_r8)* &
561 & grad(i+1,jend+1))/ &
562 & (cff+ce)
563
564 IF (lbc(inorth,istvar(itrc),ng)%nudging) THEN
565 t(i,jend+1,k,nout,itrc)=t(i,jend+1,k,nout,itrc)+ &
566 & tau* &
567 & (boundary(ng)%t_north(i,k, &
568 & itrc)- &
569 & t(i,jend+1,k,nstp,itrc))
570 END IF
571# ifdef MASKING
572 t(i,jend+1,k,nout,itrc)=t(i,jend+1,k,nout,itrc)* &
573 & grid(ng)%rmask(i,jend+1)
574# endif
575 END IF
576 END DO
577 END DO
578!
579! Northern edge, clamped boundary condition.
580!
581 ELSE IF (lbc(inorth,istvar(itrc),ng)%clamped) THEN
582 DO k=1,n(ng)
583 DO i=istr,iend
584 IF (lbc_apply(ng)%north(i)) THEN
585 t(i,jend+1,k,nout,itrc)=boundary(ng)%t_north(i,k,itrc)
586# ifdef MASKING
587 t(i,jend+1,k,nout,itrc)=t(i,jend+1,k,nout,itrc)* &
588 & grid(ng)%rmask(i,jend+1)
589# endif
590 END IF
591 END DO
592 END DO
593!
594! Northern edge, gradient boundary condition.
595!
596 ELSE IF (lbc(inorth,istvar(itrc),ng)%gradient) THEN
597 DO k=1,n(ng)
598 DO i=istr,iend
599 IF (lbc_apply(ng)%north(i)) THEN
600 t(i,jend+1,k,nout,itrc)=t(i,jend,k,nout,itrc)
601# ifdef MASKING
602 t(i,jend+1,k,nout,itrc)=t(i,jend+1,k,nout,itrc)* &
603 & grid(ng)%rmask(i,jend+1)
604# endif
605 END IF
606 END DO
607 END DO
608!
609! Northern edge, closed boundary condition.
610!
611 ELSE IF (lbc(inorth,istvar(itrc),ng)%closed) THEN
612 DO k=1,n(ng)
613 DO i=istr,iend
614 IF (lbc_apply(ng)%north(i)) THEN
615 t(i,jend+1,k,nout,itrc)=t(i,jend,k,nout,itrc)
616# ifdef MASKING
617 t(i,jend+1,k,nout,itrc)=t(i,jend+1,k,nout,itrc)* &
618 & grid(ng)%rmask(i,jend+1)
619# endif
620 END IF
621 END DO
622 END DO
623 END IF
624 END IF
625!
626!-----------------------------------------------------------------------
627! Boundary corners.
628!-----------------------------------------------------------------------
629!
630 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
631 IF (domain(ng)%SouthWest_Corner(tile)) THEN
632 IF (lbc_apply(ng)%south(istr-1).and. &
633 & lbc_apply(ng)%west (jstr-1)) THEN
634 DO k=1,n(ng)
635 t(istr-1,jstr-1,k,nout,itrc)=0.5_r8* &
636 & (t(istr ,jstr-1,k,nout, &
637 & itrc)+ &
638 & t(istr-1,jstr ,k,nout, &
639 & itrc))
640 END DO
641 END IF
642 END IF
643 IF (domain(ng)%SouthEast_Corner(tile)) THEN
644 IF (lbc_apply(ng)%south(iend+1).and. &
645 & lbc_apply(ng)%east (jstr-1)) THEN
646 DO k=1,n(ng)
647 t(iend+1,jstr-1,k,nout,itrc)=0.5_r8* &
648 & (t(iend ,jstr-1,k,nout, &
649 & itrc)+ &
650 & t(iend+1,jstr ,k,nout, &
651 & itrc))
652 END DO
653 END IF
654 END IF
655 IF (domain(ng)%NorthWest_Corner(tile)) THEN
656 IF (lbc_apply(ng)%north(istr-1).and. &
657 & lbc_apply(ng)%west (jend+1)) THEN
658 DO k=1,n(ng)
659 t(istr-1,jend+1,k,nout,itrc)=0.5_r8* &
660 & (t(istr-1,jend ,k,nout, &
661 & itrc)+ &
662 & t(istr ,jend+1,k,nout, &
663 & itrc))
664 END DO
665 END IF
666 END IF
667 IF (domain(ng)%NorthEast_Corner(tile)) THEN
668 IF (lbc_apply(ng)%north(iend+1).and. &
669 & lbc_apply(ng)%east (jend+1)) THEN
670 DO k=1,n(ng)
671 t(iend+1,jend+1,k,nout,itrc)=0.5_r8* &
672 & (t(iend+1,jend ,k,nout, &
673 & itrc)+ &
674 & t(iend ,jend+1,k,nout, &
675 & itrc))
676 END DO
677 END IF
678 END IF
679 END IF
680
681 RETURN
682 END SUBROUTINE t3dbc_tile
683#endif
684 END MODULE t3dbc_mod
type(t_boundary), dimension(:), allocatable boundary
type(t_apply), dimension(:), allocatable lbc_apply
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer, dimension(:), allocatable istvar
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
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
real(dp), dimension(:), allocatable dt
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
real(dp), dimension(:,:,:), allocatable tobc_out
real(dp), dimension(:,:,:), allocatable tobc_in
real(dp), dimension(:), allocatable obcfac
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth
logical, dimension(:,:), allocatable lnudgetclm
integer, dimension(:), allocatable nstp
subroutine t3dbc(ng, tile, nout, itrc, ic)
Definition t3dbc_im.F:26
subroutine, public t3dbc_tile(ng, tile, itrc, ic, lbi, ubi, lbj, ubj, ubk, ubt, imins, imaxs, jmins, jmaxs, nstp, nout, t)
Definition t3dbc_im.F:55