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

Functions/Subroutines

subroutine, public tl_t3dbc (ng, tile, nout, itrc, ic)
 
subroutine, public tl_t3dbc_tile (ng, tile, itrc, ic, lbi, ubi, lbj, ubj, ubk, ubt, imins, imaxs, jmins, jmaxs, nstp, nout, tl_t)
 

Function/Subroutine Documentation

◆ tl_t3dbc()

subroutine, public tl_t3dbc_mod::tl_t3dbc ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) nout,
integer, intent(in) itrc,
integer, intent(in) ic )

Definition at line 28 of file tl_t3dbc_im.F.

29!***********************************************************************
30!
31 USE mod_param
32 USE mod_ocean
33 USE mod_stepping
34!
35! Imported variable declarations.
36!
37 integer, intent(in) :: ng, tile, nout, itrc, ic
38!
39! Local variable declarations.
40!
41# include "tile.h"
42!
43 CALL tl_t3dbc_tile (ng, tile, itrc, ic, &
44 & lbi, ubi, lbj, ubj, n(ng), nt(ng), &
45 & imins, imaxs, jmins, jmaxs, &
46 & nstp(ng), nout, &
47 & ocean(ng)% tl_t)
48
49 RETURN
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer, dimension(:), allocatable nstp

References mod_param::n, mod_stepping::nstp, mod_param::nt, mod_ocean::ocean, and tl_t3dbc_tile().

Here is the call graph for this function:

◆ tl_t3dbc_tile()

subroutine, public tl_t3dbc_mod::tl_t3dbc_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) itrc,
integer, intent(in) ic,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) ubk,
integer, intent(in) ubt,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) nstp,
integer, intent(in) nout,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) tl_t )

Definition at line 53 of file tl_t3dbc_im.F.

58!***********************************************************************
59!
60 USE mod_param
61 USE mod_boundary
62 USE mod_clima
63 USE mod_grid
64 USE mod_ncparam
65 USE mod_scalars
66!
67! Imported variable declarations.
68!
69 integer, intent(in) :: ng, tile, itrc, ic
70 integer, intent(in) :: LBi, UBi, LBj, UBj, UBk, UBt
71 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
72 integer, intent(in) :: nstp, nout
73!
74# ifdef ASSUMED_SHAPE
75 real(r8), intent(inout) :: tl_t(LBi:,LBj:,:,:,:)
76# else
77 real(r8), intent(inout) :: tl_t(LBi:UBi,LBj:UBj,UBk,3,UBt)
78# endif
79!
80! Local variable declarations.
81!
82 integer :: i, j, k
83
84 real(r8) :: Ce, Cx, cff
85 real(r8) :: obc_in, obc_out, tau
86
87 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_grad
88
89# include "set_bounds.h"
90!
91!-----------------------------------------------------------------------
92! Lateral boundary conditions at the western edge.
93!-----------------------------------------------------------------------
94!
95 IF (domain(ng)%Western_Edge(tile)) THEN
96!
97! Western edge, implicit upstream radiation condition.
98!
99 IF (tl_lbc(iwest,istvar(itrc),ng)%radiation) THEN
100 IF (iic(ng).ne.0) THEN
101 DO k=1,n(ng)
102 DO j=jstr,jend+1
103!^ grad(Istr-1,j)=t(Istr-1,j ,k,nstp,itrc)- &
104!^ & t(Istr-1,j-1,k,nstp,itrc)
105!^
106 tl_grad(istr-1,j)=0.0_r8
107 END DO
108 DO j=jstr,jend
109 IF (lbc_apply(ng)%west(j)) THEN
110# if defined CELERITY_READ && defined FORWARD_READ
111 IF (tl_lbc(iwest,istvar(itrc),ng)%nudging) THEN
112 IF (lnudgetclm(itrc,ng)) THEN
113 obc_out=clima(ng)%Tnudgcof(istr-1,j,k,ic)
114 obc_in =obcfac(ng)*obc_out
115 ELSE
116 obc_out=tobc_out(itrc,ng,iwest)
117 obc_in =tobc_in(itrc,ng,iwest)
118 END IF
119 IF (boundary(ng)%t_west_Cx(j,k,itrc).lt. &
120 & 0.0_r8) THEN
121 tau=obc_in
122 ELSE
123 tau=obc_out
124 END IF
125 tau=tau*dt(ng)
126 END IF
127 cx=boundary(ng)%t_west_Cx(j,k,itrc)
128# ifdef RADIATION_2D
129 ce=boundary(ng)%t_west_Ce(j,k,itrc)
130# else
131 ce=0.0_r8
132# endif
133 cff=boundary(ng)%t_west_C2(j,k,itrc)
134# endif
135!^ t(Istr-1,j,k,nout,itrc)=(cff*t(Istr-1,j,k,nstp,itrc)+ &
136!^ & Cx *t(Istr ,j,k,nout,itrc)- &
137!^ & MAX(Ce,0.0_r8)* &
138!^ & grad(Istr-1,j )- &
139!^ & MIN(Ce,0.0_r8)* &
140!^ & grad(Istr-1,j+1))/ &
141!^ & (cff+Cx)
142!^
143 tl_t(istr-1,j,k,nout,itrc)=(cff* &
144 & tl_t(istr-1,j,k,nstp, &
145 & itrc)+ &
146 & cx* &
147 & tl_t(istr ,j,k,nout, &
148 & itrc)- &
149 & max(ce,0.0_r8)* &
150 & tl_grad(istr-1,j )- &
151 & min(ce,0.0_r8)* &
152 & tl_grad(istr-1,j+1))/ &
153 & (cff+cx)
154
155 IF (tl_lbc(iwest,istvar(itrc),ng)%nudging) THEN
156!^ t(Istr-1,j,k,nout,itrc)=t(Istr-1,j,k,nout,itrc)+ &
157!^ & tau* &
158!^ & (BOUNDARY(ng)% &
159!^ & t_west(j,k,itrc)- &
160!^ & t(Istr-1,j,k,nstp,itrc))
161!^
162 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr-1,j,k,nout, &
163 & itrc)- &
164 & tau* &
165 & tl_t(istr-1,j,k,nstp, &
166 & itrc)
167 END IF
168# ifdef MASKING
169!^ t(Istr-1,j,k,nout,itrc)=t(Istr-1,j,k,nout,itrc)* &
170!^ & GRID(ng)%rmask(Istr-1,j)
171!^
172 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr-1,j,k,nout, &
173 & itrc)* &
174 & grid(ng)%rmask(istr-1,j)
175# endif
176 END IF
177 END DO
178 END DO
179 END IF
180!
181! Western edge, clamped boundary condition.
182!
183 ELSE IF (tl_lbc(iwest,istvar(itrc),ng)%clamped) THEN
184 DO k=1,n(ng)
185 DO j=jstr,jend
186 IF (lbc_apply(ng)%west(j)) THEN
187!^ t(Istr-1,j,k,nout,itrc)=BOUNDARY(ng)%t_west(j,k,itrc)
188!^
189# ifdef ADJUST_BOUNDARY
190 IF (lobc(iwest,istvar(itrc),ng)) THEN
191 tl_t(istr-1,j,k,nout,itrc)=boundary(ng)% &
192 & tl_t_west(j,k,itrc)
193 ELSE
194 tl_t(istr-1,j,k,nout,itrc)=0.0_r8
195 END IF
196# else
197 tl_t(istr-1,j,k,nout,itrc)=0.0_r8
198# endif
199
200# ifdef MASKING
201!^ t(Istr-1,j,k,nout,itrc)=t(Istr-1,j,k,nout,itrc)* &
202!^ & GRID(ng)%rmask(Istr-1,j)
203!^
204 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr-1,j,k,nout,itrc)* &
205 & grid(ng)%rmask(istr-1,j)
206# endif
207 END IF
208 END DO
209 END DO
210!
211! Western edge, gradient boundary condition.
212!
213 ELSE IF (tl_lbc(iwest,istvar(itrc),ng)%gradient) THEN
214 DO k=1,n(ng)
215 DO j=jstr,jend
216 IF (lbc_apply(ng)%west(j)) THEN
217!^ t(Istr-1,j,k,nout,itrc)=t(Istr,j,k,nout,itrc)
218!^
219 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr,j,k,nout,itrc)
220# ifdef MASKING
221!^ t(Istr-1,j,k,nout,itrc)=t(Istr-1,j,k,nout,itrc)* &
222!^ & GRID(ng)%rmask(Istr-1,j)
223!^
224 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr-1,j,k,nout,itrc)* &
225 & grid(ng)%rmask(istr-1,j)
226# endif
227 END IF
228 END DO
229 END DO
230!
231! Western edge, closed boundary condition.
232!
233 ELSE IF (tl_lbc(iwest,istvar(itrc),ng)%closed) THEN
234 DO k=1,n(ng)
235 DO j=jstr,jend
236 IF (lbc_apply(ng)%west(j)) THEN
237!^ t(Istr-1,j,k,nout,itrc)=t(Istr,j,k,nout,itrc)
238!^
239 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr,j,k,nout,itrc)
240# ifdef MASKING
241!^ t(Istr-1,j,k,nout,itrc)=t(Istr-1,j,k,nout,itrc)* &
242!^ & GRID(ng)%rmask(Istr-1,j)
243!^
244 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr-1,j,k,nout,itrc)* &
245 & grid(ng)%rmask(istr-1,j)
246# endif
247 END IF
248 END DO
249 END DO
250 END IF
251 END IF
252!
253!-----------------------------------------------------------------------
254! Lateral boundary conditions at the eastern edge.
255!-----------------------------------------------------------------------
256!
257 IF (domain(ng)%Eastern_Edge(tile)) THEN
258!
259! Eastern edge, implicit upstream radiation condition.
260!
261 IF (tl_lbc(ieast,istvar(itrc),ng)%radiation) THEN
262 IF (iic(ng).ne.0) THEN
263 DO k=1,n(ng)
264 DO j=jstr,jend+1
265!^ grad(Iend+1,j)=t(Iend+1,j ,k,nstp,itrc)- &
266!^ & t(Iend+1,j-1,k,nstp,itrc)
267!^
268 tl_grad(iend+1,j)=0.0_r8
269 END DO
270 DO j=jstr,jend
271 IF (lbc_apply(ng)%east(j)) THEN
272# if defined CELERITY_READ && defined FORWARD_READ
273 IF (tl_lbc(ieast,istvar(itrc),ng)%nudging) THEN
274 IF (lnudgetclm(itrc,ng)) THEN
275 obc_out=clima(ng)%Tnudgcof(iend+1,j,k,ic)
276 obc_in =obcfac(ng)*obc_out
277 ELSE
278 obc_out=tobc_out(itrc,ng,ieast)
279 obc_in =tobc_in(itrc,ng,ieast)
280 END IF
281 IF (boundary(ng)%t_east_Cx(j,k,itrc).lt. &
282 & 0.0_r8) THEN
283 tau=obc_in
284 ELSE
285 tau=obc_out
286 END IF
287 tau=tau*dt(ng)
288 END IF
289 cx=boundary(ng)%t_east_Cx(j,k,itrc)
290# ifdef RADIATION_2D
291 ce=boundary(ng)%t_east_Ce(j,k,itrc)
292# else
293 ce=0.0_r8
294# endif
295 cff=boundary(ng)%t_east_C2(j,k,itrc)
296# endif
297!^ t(Iend+1,j,k,nout,itrc)=(cff*t(Iend+1,j,k,nstp,itrc)+ &
298!^ & Cx *t(Iend ,j,k,nout,itrc)- &
299!^ & MAX(Ce,0.0_r8)* &
300!^ & grad(Iend+1,j )- &
301!^ & MIN(Ce,0.0_r8)* &
302!^ & grad(Iend+1,j+1))/ &
303!^ & (cff+Cx)
304!^
305 tl_t(iend+1,j,k,nout,itrc)=(cff* &
306 & tl_t(iend+1,j,k,nstp, &
307 & itrc)+ &
308 & cx* &
309 & tl_t(iend ,j,k,nout, &
310 & itrc)- &
311 & max(ce,0.0_r8)* &
312 & tl_grad(iend+1,j )- &
313 & min(ce,0.0_r8)* &
314 & tl_grad(iend+1,j+1))/ &
315 & (cff+cx)
316
317 IF (tl_lbc(ieast,istvar(itrc),ng)%nudging) THEN
318!^ t(Iend+1,j,k,nout,itrc)=t(Iend+1,j,k,nout,itrc)+ &
319!^ & tau* &
320!^ & (BOUNDARY(ng)% &
321!^ & t_east(j,k,itrc)- &
322!^ & t(Iend+1,j,k,nstp,itrc))
323!^
324 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend+1,j,k,nout, &
325 & itrc)- &
326 & tau* &
327 & tl_t(iend+1,j,k,nstp, &
328 & itrc)
329 END IF
330# ifdef MASKING
331!^ t(Iend+1,j,k,nout,itrc)=t(Iend+1,j,k,nout,itrc)* &
332!^ & GRID(ng)%rmask(Iend+1,j)
333!^
334 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend+1,j,k,nout, &
335 & itrc)* &
336 & grid(ng)%rmask(iend+1,j)
337# endif
338 END IF
339 END DO
340 END DO
341 END IF
342!
343! Eastern edge, clamped boundary condition.
344!
345 ELSE IF (tl_lbc(ieast,istvar(itrc),ng)%clamped) THEN
346 DO k=1,n(ng)
347 DO j=jstr,jend
348 IF (lbc_apply(ng)%east(j)) THEN
349!^ t(Iend+1,j,k,nout,itrc)=BOUNDARY(ng)%t_east(j,k,itrc)
350!^
351# ifdef ADJUST_BOUNDARY
352 IF (lobc(ieast,istvar(itrc),ng)) THEN
353 tl_t(iend+1,j,k,nout,itrc)=boundary(ng)% &
354 & tl_t_east(j,k,itrc)
355 ELSE
356 tl_t(iend+1,j,k,nout,itrc)=0.0_r8
357 END IF
358# else
359 tl_t(iend+1,j,k,nout,itrc)=0.0_r8
360# endif
361# ifdef MASKING
362!^ t(Iend+1,j,k,nout,itrc)=t(Iend+1,j,k,nout,itrc)* &
363!^ & GRID(ng)%rmask(Iend+1,j)
364!^
365 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend+1,j,k,nout,itrc)* &
366 & grid(ng)%rmask(iend+1,j)
367# endif
368 END IF
369 END DO
370 END DO
371!
372! Eastern edge, gradient boundary condition.
373!
374 ELSE IF (tl_lbc(ieast,istvar(itrc),ng)%gradient) THEN
375 DO k=1,n(ng)
376 DO j=jstr,jend
377 IF (lbc_apply(ng)%east(j)) THEN
378!^ t(Iend+1,j,k,nout,itrc)=t(Iend,j,k,nout,itrc)
379!^
380 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend,j,k,nout,itrc)
381# ifdef MASKING
382!^ t(Iend+1,j,k,nout,itrc)=t(Iend+1,j,k,nout,itrc)* &
383!^ & GRID(ng)%rmask(Iend+1,j)
384!^
385 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend+1,j,k,nout,itrc)* &
386 & grid(ng)%rmask(iend+1,j)
387# endif
388 END IF
389 END DO
390 END DO
391!
392! Eastern edge, closed boundary condition.
393!
394 ELSE IF (tl_lbc(ieast,istvar(itrc),ng)%closed) THEN
395 DO k=1,n(ng)
396 DO j=jstr,jend
397 IF (lbc_apply(ng)%east(j)) THEN
398!^ t(Iend+1,j,k,nout,itrc)=t(Iend,j,k,nout,itrc)
399!^
400 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend,j,k,nout,itrc)
401# ifdef MASKING
402!^ t(Iend+1,j,k,nout,itrc)=t(Iend+1,j,k,nout,itrc)* &
403!^ & GRID(ng)%rmask(Iend+1,j)
404!^
405 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend+1,j,k,nout,itrc)* &
406 & grid(ng)%rmask(iend+1,j)
407# endif
408 END IF
409 END DO
410 END DO
411 END IF
412 END IF
413!
414!-----------------------------------------------------------------------
415! Lateral boundary conditions at the southern edge.
416!-----------------------------------------------------------------------
417!
418 IF (domain(ng)%Southern_Edge(tile)) THEN
419!
420! Southern edge, implicit upstream radiation condition.
421!
422 IF (tl_lbc(isouth,istvar(itrc),ng)%radiation) THEN
423 IF (iic(ng).ne.0) THEN
424 DO k=1,n(ng)
425 DO i=istr,iend+1
426!^ grad(i,Jstr-1)=t(i ,Jstr-1,k,nstp,itrc)- &
427!^ & t(i-1,Jstr-1,k,nstp,itrc)
428!^
429 tl_grad(i,jstr-1)=0.0_r8
430 END DO
431 DO i=istr,iend
432 IF (lbc_apply(ng)%south(i)) THEN
433# if defined CELERITY_READ && defined FORWARD_READ
434 IF (tl_lbc(isouth,istvar(itrc),ng)%nudging) THEN
435 IF (lnudgetclm(itrc,ng)) THEN
436 obc_out=clima(ng)%Tnudgcof(i,jstr-1,k,ic)
437 obc_in =obcfac(ng)*obc_out
438 ELSE
439 obc_out=tobc_out(itrc,ng,isouth)
440 obc_in =tobc_in(itrc,ng,isouth)
441 END IF
442 IF (boundary(ng)%t_south_Ce(i,k,itrc).lt. &
443 & 0.0_r8) THEN
444 tau=obc_in
445 ELSE
446 tau=obc_out
447 END IF
448 tau=tau*dt(ng)
449 END IF
450# ifdef RADIATION_2D
451 cx=boundary(ng)%t_south_Cx(i,k,itrc)
452# else
453 cx=0.0_r8
454# endif
455 ce=boundary(ng)%t_south_Ce(i,k,itrc)
456 cff=boundary(ng)%t_south_C2(i,k,itrc)
457# endif
458!^ t(i,Jstr-1,k,nout,itrc)=(cff*t(i,Jstr-1,k,nstp,itrc)+ &
459!^ & Ce *t(i,Jstr ,k,nout,itrc )-&
460!^ & MAX(Cx,0.0_r8)* &
461!^ & grad(i ,Jstr-1)- &
462!^ & MIN(Cx,0.0_r8)* &
463!^ & grad(i+1,Jstr-1))/ &
464!^ & (cff+Ce)
465!^
466 tl_t(i,jstr-1,k,nout,itrc)=(cff* &
467 & tl_t(i,jstr-1,k,nstp, &
468 & itrc)+ &
469 & ce* &
470 & tl_t(i,jstr ,k,nout, &
471 & itrc)- &
472 & max(cx,0.0_r8)* &
473 & tl_grad(i ,jstr-1)- &
474 & min(cx,0.0_r8)* &
475 & tl_grad(i+1,jstr-1))/ &
476 & (cff+ce)
477
478 IF (tl_lbc(isouth,istvar(itrc),ng)%nudging) THEN
479!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr-1,k,nout,itrc)+ &
480!^ & tau* &
481!^ & (BOUNDARY(ng)% &
482!^ & t_south(i,k,itrc)- &
483!^ & t(i,Jstr-1,k,nstp,itrc))
484!^
485 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr-1,k,nout, &
486 & itrc)- &
487 & tau* &
488 & tl_t(i,jstr-1,k,nstp, &
489 & itrc)
490 END IF
491# ifdef MASKING
492!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr-1,k,nout,itrc)* &
493!^ & GRID(ng)%rmask(i,Jstr-1)
494!^
495 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr-1,k,nout, &
496 & itrc)* &
497 & grid(ng)%rmask(i,jstr-1)
498# endif
499 END IF
500 END DO
501 END DO
502 END IF
503!
504! Southern edge, clamped boundary condition.
505!
506 ELSE IF (tl_lbc(isouth,istvar(itrc),ng)%clamped) THEN
507 DO k=1,n(ng)
508 DO i=istr,iend
509 IF (lbc_apply(ng)%south(i)) THEN
510!^ t(i,Jstr-1,k,nout,itrc)=BOUNDARY(ng)%t_south(i,k,itrc)
511!^
512# ifdef ADJUST_BOUNDARY
513 IF (lobc(isouth,istvar(itrc),ng)) THEN
514 tl_t(i,jstr-1,k,nout,itrc)=boundary(ng)% &
515 & tl_t_south(i,k,itrc)
516 ELSE
517 tl_t(i,jstr-1,k,nout,itrc)=0.0_r8
518 END IF
519# else
520 tl_t(i,jstr-1,k,nout,itrc)=0.0_r8
521# endif
522# ifdef MASKING
523!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr-1,k,nout,itrc)* &
524!^ & GRID(ng)%rmask(i,Jstr-1)
525!^
526 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr-1,k,nout,itrc)* &
527 & grid(ng)%rmask(i,jstr-1)
528# endif
529 END IF
530 END DO
531 END DO
532!
533! Southern edge, gradient boundary condition.
534!
535 ELSE IF (tl_lbc(isouth,istvar(itrc),ng)%gradient) THEN
536 DO k=1,n(ng)
537 DO i=istr,iend
538 IF (lbc_apply(ng)%south(i)) THEN
539!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr,k,nout,itrc)
540!^
541 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr,k,nout,itrc)
542# ifdef MASKING
543!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr-1,k,nout,itrc)* &
544!^ & GRID(ng)%rmask(i,Jstr-1)
545!^
546 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr-1,k,nout,itrc)* &
547 & grid(ng)%rmask(i,jstr-1)
548# endif
549 END IF
550 END DO
551 END DO
552!
553! Southern edge, closed boundary condition.
554!
555 ELSE IF (tl_lbc(isouth,istvar(itrc),ng)%closed) THEN
556 DO k=1,n(ng)
557 DO i=istr,iend
558 IF (lbc_apply(ng)%south(i)) THEN
559!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr,k,nout,itrc)
560!^
561 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr,k,nout,itrc)
562# ifdef MASKING
563!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr-1,k,nout,itrc)* &
564!^ & GRID(ng)%rmask(i,Jstr-1)
565!^
566 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr-1,k,nout,itrc)* &
567 & grid(ng)%rmask(i,jstr-1)
568# endif
569 END IF
570 END DO
571 END DO
572 END IF
573 END IF
574!
575!-----------------------------------------------------------------------
576! Lateral boundary conditions at the northern edge.
577!-----------------------------------------------------------------------
578!
579 IF (domain(ng)%Northern_Edge(tile)) THEN
580!
581! Northern edge, implicit upstream radiation condition.
582!
583 IF (tl_lbc(inorth,istvar(itrc),ng)%radiation) THEN
584 IF (iic(ng).ne.0) THEN
585 DO k=1,n(ng)
586 DO i=istr,iend+1
587!^ grad(i,Jend+1)=t(i ,Jend+1,k,nstp,itrc)- &
588!^ & t(i-1,Jend+1,k,nstp,itrc)
589!^
590 tl_grad(i,jend+1)=0.0_r8
591 END DO
592 DO i=istr,iend
593 IF (lbc_apply(ng)%north(i)) THEN
594# if defined CELERITY_READ && defined FORWARD_READ
595 IF (tl_lbc(inorth,istvar(itrc),ng)%nudging) THEN
596 IF (lnudgetclm(itrc,ng)) THEN
597 obc_out=clima(ng)%Tnudgcof(i,jend+1,k,ic)
598 obc_in =obcfac(ng)*obc_out
599 ELSE
600 obc_out=tobc_out(itrc,ng,inorth)
601 obc_in =tobc_in(itrc,ng,inorth)
602 END IF
603 IF (boundary(ng)%t_north_Ce(i,k,itrc).lt. &
604 & 0.0_r8) THEN
605 tau=obc_in
606 ELSE
607 tau=obc_out
608 END IF
609 tau=tau*dt(ng)
610 END IF
611# ifdef RADIATION_2D
612 cx=boundary(ng)%t_north_Cx(i,k,itrc)
613# else
614 cx=0.0_r8
615# endif
616 ce=boundary(ng)%t_north_Ce(i,k,itrc)
617 cff=boundary(ng)%t_north_C2(i,k,itrc)
618# endif
619!^ t(i,Jend+1,k,nout,itrc)=(cff*t(i,Jend+1,k,nstp,itrc)+ &
620!^ & Ce *t(i,Jend ,k,nout,itrc)- &
621!^ & MAX(Cx,0.0_r8)* &
622!^ & grad(i ,Jend+1)- &
623!^ & MIN(Cx,0.0_r8)* &
624!^ & grad(i+1,Jend+1))/ &
625!^ & (cff+Ce)
626!^
627 tl_t(i,jend+1,k,nout,itrc)=(cff* &
628 & tl_t(i,jend+1,k,nstp, &
629 & itrc)+ &
630 & ce* &
631 & tl_t(i,jend ,k,nout, &
632 & itrc)- &
633 & max(cx,0.0_r8)* &
634 & tl_grad(i ,jend+1)- &
635 & min(cx,0.0_r8)* &
636 & tl_grad(i+1,jend+1))/ &
637 & (cff+ce)
638
639 IF (tl_lbc(inorth,istvar(itrc),ng)%nudging) THEN
640!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend+1,k,nout,itrc)+ &
641!^ & tau* &
642!^ & (BOUNDARY(ng)% &
643!^ & t_north(i,k,itrc)- &
644!^ & t(i,Jend+1,k,nstp,itrc))
645!^
646 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend+1,k,nout, &
647 & itrc)- &
648 & tau* &
649 & tl_t(i,jend+1,k,nstp, &
650 & itrc)
651 END IF
652# ifdef MASKING
653!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend+1,k,nout,itrc)* &
654!^ & GRID(ng)%rmask(i,Jend+1)
655!^
656 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend+1,k,nout, &
657 & itrc)* &
658 & grid(ng)%rmask(i,jend+1)
659# endif
660 END IF
661 END DO
662 END DO
663 END IF
664!
665! Northern edge, clamped boundary condition.
666!
667 ELSE IF (tl_lbc(inorth,istvar(itrc),ng)%clamped) THEN
668 DO k=1,n(ng)
669 DO i=istr,iend
670 IF (lbc_apply(ng)%north(i)) THEN
671!^ t(i,Jend+1,k,nout,itrc)=BOUNDARY(ng)%t_north(i,k,itrc)
672!^
673# ifdef ADJUST_BOUNDARY
674 IF (lobc(inorth,istvar(itrc),ng)) THEN
675 tl_t(i,jend+1,k,nout,itrc)=boundary(ng)% &
676 & tl_t_north(i,k,itrc)
677 ELSE
678 tl_t(i,jend+1,k,nout,itrc)=0.0_r8
679 END IF
680# else
681 tl_t(i,jend+1,k,nout,itrc)=0.0_r8
682# endif
683# ifdef MASKING
684!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend+1,k,nout,itrc)* &
685!^ & GRID(ng)%rmask(i,Jend+1)
686!^
687 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend+1,k,nout,itrc)* &
688 & grid(ng)%rmask(i,jend+1)
689# endif
690 END IF
691 END DO
692 END DO
693!
694! Northern edge, gradient boundary condition.
695!
696 ELSE IF (tl_lbc(inorth,istvar(itrc),ng)%gradient) THEN
697 DO k=1,n(ng)
698 DO i=istr,iend
699 IF (lbc_apply(ng)%north(i)) THEN
700!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend,k,nout,itrc)
701!^
702 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend,k,nout,itrc)
703# ifdef MASKING
704!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend+1,k,nout,itrc)* &
705!^ & GRID(ng)%rmask(i,Jend+1)
706!^
707 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend+1,k,nout,itrc)* &
708 & grid(ng)%rmask(i,jend+1)
709# endif
710 END IF
711 END DO
712 END DO
713!
714! Northern edge, closed boundary condition.
715!
716 ELSE IF (tl_lbc(inorth,istvar(itrc),ng)%closed) THEN
717 DO k=1,n(ng)
718 DO i=istr,iend
719 IF (lbc_apply(ng)%north(i)) THEN
720!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend,k,nout,itrc)
721!^
722 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend,k,nout,itrc)
723# ifdef MASKING
724!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend+1,k,nout,itrc)* &
725!^ & GRID(ng)%rmask(i,Jend+1)
726!^
727 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend+1,k,nout,itrc)* &
728 & grid(ng)%rmask(i,jend+1)
729# endif
730 END IF
731 END DO
732 END DO
733 END IF
734 END IF
735!
736!-----------------------------------------------------------------------
737! Boundary corners.
738!-----------------------------------------------------------------------
739!
740 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
741 IF (domain(ng)%SouthWest_Corner(tile)) THEN
742 IF (lbc_apply(ng)%south(istr-1).and. &
743 & lbc_apply(ng)%west (jstr-1)) THEN
744 DO k=1,n(ng)
745!^ t(Istr-1,Jstr-1,k,nout,itrc)=0.5_r8* &
746!^ & (t(Istr ,Jstr-1,k,nout, &
747!^ & itrc)+ &
748!^ & t(Istr-1,Jstr ,k,nout, &
749!^ & itrc))
750!^
751 tl_t(istr-1,jstr-1,k,nout,itrc)=0.5_r8* &
752 & (tl_t(istr ,jstr-1,k, &
753 & nout,itrc)+ &
754 & tl_t(istr-1,jstr ,k, &
755 & nout,itrc))
756 END DO
757 END IF
758 END IF
759 IF (domain(ng)%SouthEast_Corner(tile)) THEN
760 IF (lbc_apply(ng)%south(iend+1).and. &
761 & lbc_apply(ng)%east (jstr-1)) THEN
762 DO k=1,n(ng)
763!^ t(Iend+1,Jstr-1,k,nout,itrc)=0.5_r8* &
764!^ & (t(Iend ,Jstr-1,k,nout, &
765!^ & itrc)+ &
766!^ & t(Iend+1,Jstr ,k,nout, &
767!^ & itrc))
768!^
769 tl_t(iend+1,jstr-1,k,nout,itrc)=0.5_r8* &
770 & (tl_t(iend ,jstr-1,k, &
771 & nout,itrc)+ &
772 & tl_t(iend+1,jstr ,k, &
773 & nout,itrc))
774 END DO
775 END IF
776 END IF
777 IF (domain(ng)%NorthWest_Corner(tile)) THEN
778 IF (lbc_apply(ng)%north(istr-1).and. &
779 & lbc_apply(ng)%west (jend+1)) THEN
780 DO k=1,n(ng)
781!^ t(Istr-1,Jend+1,k,nout,itrc)=0.5_r8* &
782!^ & (t(Istr-1,Jend ,k,nout, &
783!^ & itrc)+ &
784!^ & t(Istr ,Jend+1,k,nout, &
785!^ & itrc))
786!^
787 tl_t(istr-1,jend+1,k,nout,itrc)=0.5_r8* &
788 & (tl_t(istr-1,jend ,k, &
789 & nout,itrc)+ &
790 & tl_t(istr ,jend+1,k, &
791 & nout,itrc))
792 END DO
793 END IF
794 END IF
795 IF (domain(ng)%NorthEast_Corner(tile)) THEN
796 IF (lbc_apply(ng)%north(iend+1).and. &
797 & lbc_apply(ng)%east (jend+1)) THEN
798 DO k=1,n(ng)
799!^ t(Iend+1,Jend+1,k,nout,itrc)=0.5_r8* &
800!^ & (t(Iend+1,Jend ,k,nout, &
801!^ & itrc)+ &
802!^ & t(Iend ,Jend+1,k,nout, &
803!^ & itrc))
804!^
805 tl_t(iend+1,jend+1,k,nout,itrc)=0.5_r8* &
806 & (tl_t(iend+1,jend ,k, &
807 & nout,itrc)+ &
808 & tl_t(iend ,jend+1,k, &
809 & nout,itrc))
810 END DO
811 END IF
812 END IF
813 END IF
814
815 RETURN
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_lbc), dimension(:,:,:), allocatable tl_lbc
Definition mod_param.F:379
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
logical, dimension(:,:,:), allocatable lobc
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

References mod_boundary::boundary, mod_clima::clima, mod_param::domain, mod_scalars::dt, mod_scalars::ewperiodic, mod_grid::grid, mod_scalars::ieast, mod_scalars::iic, mod_scalars::inorth, mod_scalars::isouth, mod_ncparam::istvar, mod_scalars::iwest, mod_boundary::lbc_apply, mod_scalars::lnudgetclm, mod_scalars::lobc, mod_param::n, mod_scalars::nsperiodic, mod_scalars::obcfac, mod_param::tl_lbc, mod_scalars::tobc_in, and mod_scalars::tobc_out.

Referenced by tl_ini_fields_mod::tl_ini_fields_tile(), ini_adjust_mod::tl_ini_perturb_tile(), tl_pre_step3d_mod::tl_pre_step3d_tile(), tl_step3d_t_mod::tl_step3d_t_tile(), and tl_t3dbc().

Here is the caller graph for this function: