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

Functions/Subroutines

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

Function/Subroutine Documentation

◆ rp_t3dbc()

subroutine, public rp_t3dbc_mod::rp_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 rp_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 rp_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 rp_t3dbc_tile().

Here is the call graph for this function:

◆ rp_t3dbc_tile()

subroutine, public rp_t3dbc_mod::rp_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 rp_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 tl_t(istr-1,j,k,nout,itrc)=boundary(ng)% &
190 & tl_t_west(j,k,itrc)
191# ifdef MASKING
192!^ t(Istr-1,j,k,nout,itrc)=t(Istr-1,j,k,nout,itrc)* &
193!^ & GRID(ng)%rmask(Istr-1,j)
194!^
195 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr-1,j,k,nout,itrc)* &
196 & grid(ng)%rmask(istr-1,j)
197# endif
198 END IF
199 END DO
200 END DO
201!
202! Western edge, gradient boundary condition.
203!
204 ELSE IF (tl_lbc(iwest,istvar(itrc),ng)%gradient) THEN
205 DO k=1,n(ng)
206 DO j=jstr,jend
207 IF (lbc_apply(ng)%west(j)) THEN
208!^ t(Istr-1,j,k,nout,itrc)=t(Istr,j,k,nout,itrc)
209!^
210 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr,j,k,nout,itrc)
211# ifdef MASKING
212!^ t(Istr-1,j,k,nout,itrc)=t(Istr-1,j,k,nout,itrc)* &
213!^ & GRID(ng)%rmask(Istr-1,j)
214!^
215 tl_t(istr-1,j,k,nout,itrc)=tl_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!
222! Western edge, closed boundary condition.
223!
224 ELSE IF (tl_lbc(iwest,istvar(itrc),ng)%closed) THEN
225 DO k=1,n(ng)
226 DO j=jstr,jend
227 IF (lbc_apply(ng)%west(j)) THEN
228!^ t(Istr-1,j,k,nout,itrc)=t(Istr,j,k,nout,itrc)
229!^
230 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr,j,k,nout,itrc)
231# ifdef MASKING
232!^ t(Istr-1,j,k,nout,itrc)=t(Istr-1,j,k,nout,itrc)* &
233!^ & GRID(ng)%rmask(Istr-1,j)
234!^
235 tl_t(istr-1,j,k,nout,itrc)=tl_t(istr-1,j,k,nout,itrc)* &
236 & grid(ng)%rmask(istr-1,j)
237# endif
238 END IF
239 END DO
240 END DO
241 END IF
242 END IF
243!
244!-----------------------------------------------------------------------
245! Lateral boundary conditions at the eastern edge.
246!-----------------------------------------------------------------------
247!
248 IF (domain(ng)%Eastern_Edge(tile)) THEN
249!
250! Eastern edge, implicit upstream radiation condition.
251!
252 IF (tl_lbc(ieast,istvar(itrc),ng)%radiation) THEN
253 IF (iic(ng).ne.0) THEN
254 DO k=1,n(ng)
255 DO j=jstr,jend+1
256!^ grad(Iend+1,j)=t(Iend+1,j ,k,nstp,itrc)- &
257!^ & t(Iend+1,j-1,k,nstp,itrc)
258!^
259 tl_grad(iend+1,j)=0.0_r8
260 END DO
261 DO j=jstr,jend
262 IF (lbc_apply(ng)%east(j)) THEN
263# if defined CELERITY_READ && defined FORWARD_READ
264 IF (tl_lbc(ieast,istvar(itrc),ng)%nudging) THEN
265 IF (lnudgetclm(itrc,ng)) THEN
266 obc_out=clima(ng)%Tnudgcof(iend+1,j,k,ic)
267 obc_in =obcfac(ng)*obc_out
268 ELSE
269 obc_out=tobc_out(itrc,ng,ieast)
270 obc_in =tobc_in(itrc,ng,ieast)
271 END IF
272 IF (boundary(ng)%t_east_Cx(j,k,itrc).lt. &
273 & 0.0_r8) THEN
274 tau=obc_in
275 ELSE
276 tau=obc_out
277 END IF
278 tau=tau*dt(ng)
279 END IF
280 cx=boundary(ng)%t_east_Cx(j,k,itrc)
281# ifdef RADIATION_2D
282 ce=boundary(ng)%t_east_Ce(j,k,itrc)
283# else
284 ce=0.0_r8
285# endif
286 cff=boundary(ng)%t_east_C2(j,k,itrc)
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 tl_t(iend+1,j,k,nout,itrc)=(cff* &
297 & tl_t(iend+1,j,k,nstp, &
298 & itrc)+ &
299 & cx* &
300 & tl_t(iend ,j,k,nout, &
301 & itrc)- &
302 & max(ce,0.0_r8)* &
303 & tl_grad(iend+1,j )- &
304 & min(ce,0.0_r8)* &
305 & tl_grad(iend+1,j+1))/ &
306 & (cff+cx)
307
308 IF (tl_lbc(ieast,istvar(itrc),ng)%nudging) THEN
309!^ t(Iend+1,j,k,nout,itrc)=t(Iend+1,j,k,nout,itrc)+ &
310!^ & tau* &
311!^ & (BOUNDARY(ng)% &
312!^ & t_east(j,k,itrc)- &
313!^ & t(Iend+1,j,k,nstp,itrc))
314!^
315 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend+1,j,k,nout, &
316 & itrc)- &
317 & tau* &
318 & tl_t(iend+1,j,k,nstp, &
319 & itrc)
320 END IF
321# ifdef MASKING
322!^ t(Iend+1,j,k,nout,itrc)=t(Iend+1,j,k,nout,itrc)* &
323!^ & GRID(ng)%rmask(Iend+1,j)
324!^
325 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend+1,j,k,nout, &
326 & itrc)* &
327 & grid(ng)%rmask(iend+1,j)
328# endif
329 END IF
330 END DO
331 END DO
332 END IF
333!
334! Eastern edge, clamped boundary condition.
335!
336 ELSE IF (tl_lbc(ieast,istvar(itrc),ng)%clamped) THEN
337 DO k=1,n(ng)
338 DO j=jstr,jend
339 IF (lbc_apply(ng)%east(j)) THEN
340!^ t(Iend+1,j,k,nout,itrc)=BOUNDARY(ng)%t_east(j,k,itrc)
341!^
342 tl_t(iend+1,j,k,nout,itrc)=boundary(ng)% &
343 & tl_t_east(j,k,itrc)
344# ifdef MASKING
345!^ t(Iend+1,j,k,nout,itrc)=t(Iend+1,j,k,nout,itrc)* &
346!^ & GRID(ng)%rmask(Iend+1,j)
347!^
348 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend+1,j,k,nout,itrc)* &
349 & grid(ng)%rmask(iend+1,j)
350# endif
351 END IF
352 END DO
353 END DO
354!
355! Eastern edge, gradient boundary condition.
356!
357 ELSE IF (tl_lbc(ieast,istvar(itrc),ng)%gradient) THEN
358 DO k=1,n(ng)
359 DO j=jstr,jend
360 IF (lbc_apply(ng)%east(j)) THEN
361!^ t(Iend+1,j,k,nout,itrc)=t(Iend,j,k,nout,itrc)
362!^
363 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend,j,k,nout,itrc)
364# ifdef MASKING
365!^ t(Iend+1,j,k,nout,itrc)=t(Iend+1,j,k,nout,itrc)* &
366!^ & GRID(ng)%rmask(Iend+1,j)
367!^
368 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend+1,j,k,nout,itrc)* &
369 & grid(ng)%rmask(iend+1,j)
370# endif
371 END IF
372 END DO
373 END DO
374!
375! Eastern edge, closed boundary condition.
376!
377 ELSE IF (tl_lbc(ieast,istvar(itrc),ng)%closed) THEN
378 DO k=1,n(ng)
379 DO j=jstr,jend
380 IF (lbc_apply(ng)%east(j)) THEN
381!^ t(Iend+1,j,k,nout,itrc)=t(Iend,j,k,nout,itrc)
382!^
383 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend,j,k,nout,itrc)
384# ifdef MASKING
385!^ t(Iend+1,j,k,nout,itrc)=t(Iend+1,j,k,nout,itrc)* &
386!^ & GRID(ng)%rmask(Iend+1,j)
387!^
388 tl_t(iend+1,j,k,nout,itrc)=tl_t(iend+1,j,k,nout,itrc)* &
389 & grid(ng)%rmask(iend+1,j)
390# endif
391 END IF
392 END DO
393 END DO
394 END IF
395 END IF
396!
397!-----------------------------------------------------------------------
398! Lateral boundary conditions at the southern edge.
399!-----------------------------------------------------------------------
400!
401 IF (domain(ng)%Southern_Edge(tile)) THEN
402!
403! Southern edge, implicit upstream radiation condition.
404!
405 IF (tl_lbc(isouth,istvar(itrc),ng)%radiation) THEN
406 IF (iic(ng).ne.0) THEN
407 DO k=1,n(ng)
408 DO i=istr,iend+1
409!^ grad(i,Jstr-1)=t(i ,Jstr-1,k,nstp,itrc)- &
410!^ & t(i-1,Jstr-1,k,nstp,itrc)
411!^
412 tl_grad(i,jstr-1)=0.0_r8
413 END DO
414 DO i=istr,iend
415 IF (lbc_apply(ng)%south(i)) THEN
416# if defined CELERITY_READ && defined FORWARD_READ
417 IF (tl_lbc(isouth,istvar(itrc),ng)%nudging) THEN
418 IF (lnudgetclm(itrc,ng)) THEN
419 obc_out=clima(ng)%Tnudgcof(i,jstr-1,k,ic)
420 obc_in =obcfac(ng)*obc_out
421 ELSE
422 obc_out=tobc_out(itrc,ng,isouth)
423 obc_in =tobc_in(itrc,ng,isouth)
424 END IF
425 IF (boundary(ng)%t_south_Ce(i,k,itrc).lt. &
426 & 0.0_r8) THEN
427 tau=obc_in
428 ELSE
429 tau=obc_out
430 END IF
431 tau=tau*dt(ng)
432 END IF
433# ifdef RADIATION_2D
434 cx=boundary(ng)%t_south_Cx(i,k,itrc)
435# else
436 cx=0.0_r8
437# endif
438 ce=boundary(ng)%t_south_Ce(i,k,itrc)
439 cff=boundary(ng)%t_south_C2(i,k,itrc)
440# endif
441!^ t(i,Jstr-1,k,nout,itrc)=(cff*t(i,Jstr-1,k,nstp,itrc)+ &
442!^ & Ce *t(i,Jstr ,k,nout,itrc )-&
443!^ & MAX(Cx,0.0_r8)* &
444!^ & grad(i ,Jstr-1)- &
445!^ & MIN(Cx,0.0_r8)* &
446!^ & grad(i+1,Jstr-1))/ &
447!^ & (cff+Ce)
448!^
449 tl_t(i,jstr-1,k,nout,itrc)=(cff* &
450 & tl_t(i,jstr-1,k,nstp, &
451 & itrc)+ &
452 & ce* &
453 & tl_t(i,jstr ,k,nout, &
454 & itrc)- &
455 & max(cx,0.0_r8)* &
456 & tl_grad(i ,jstr-1)- &
457 & min(cx,0.0_r8)* &
458 & tl_grad(i+1,jstr-1))/ &
459 & (cff+ce)
460
461 IF (tl_lbc(isouth,istvar(itrc),ng)%nudging) THEN
462!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr-1,k,nout,itrc)+ &
463!^ & tau* &
464!^ & (BOUNDARY(ng)% &
465!^ & t_south(i,k,itrc)- &
466!^ & t(i,Jstr-1,k,nstp,itrc))
467!^
468 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr-1,k,nout, &
469 & itrc)- &
470 & tau* &
471 & tl_t(i,jstr-1,k,nstp, &
472 & itrc)
473 END IF
474# ifdef MASKING
475!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr-1,k,nout,itrc)* &
476!^ & GRID(ng)%rmask(i,Jstr-1)
477!^
478 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr-1,k,nout, &
479 & itrc)* &
480 & grid(ng)%rmask(i,jstr-1)
481# endif
482 END IF
483 END DO
484 END DO
485 END IF
486!
487! Southern edge, clamped boundary condition.
488!
489 ELSE IF (tl_lbc(isouth,istvar(itrc),ng)%clamped) THEN
490 DO k=1,n(ng)
491 DO i=istr,iend
492 IF (lbc_apply(ng)%south(i)) THEN
493!^ t(i,Jstr-1,k,nout,itrc)=BOUNDARY(ng)%t_south(i,k,itrc)
494!^
495 tl_t(i,jstr-1,k,nout,itrc)=boundary(ng)% &
496 & tl_t_south(i,k,itrc)
497# ifdef MASKING
498!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr-1,k,nout,itrc)* &
499!^ & GRID(ng)%rmask(i,Jstr-1)
500!^
501 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr-1,k,nout,itrc)* &
502 & grid(ng)%rmask(i,jstr-1)
503# endif
504 END IF
505 END DO
506 END DO
507!
508! Southern edge, gradient boundary condition.
509!
510 ELSE IF (tl_lbc(isouth,istvar(itrc),ng)%gradient) THEN
511 DO k=1,n(ng)
512 DO i=istr,iend
513 IF (lbc_apply(ng)%south(i)) THEN
514!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr,k,nout,itrc)
515!^
516 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr,k,nout,itrc)
517# ifdef MASKING
518!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr-1,k,nout,itrc)* &
519!^ & GRID(ng)%rmask(i,Jstr-1)
520!^
521 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr-1,k,nout,itrc)* &
522 & grid(ng)%rmask(i,jstr-1)
523# endif
524 END IF
525 END DO
526 END DO
527!
528! Southern edge, closed boundary condition.
529!
530 ELSE IF (tl_lbc(isouth,istvar(itrc),ng)%closed) THEN
531 DO k=1,n(ng)
532 DO i=istr,iend
533 IF (lbc_apply(ng)%south(i)) THEN
534!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr,k,nout,itrc)
535!^
536 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr,k,nout,itrc)
537# ifdef MASKING
538!^ t(i,Jstr-1,k,nout,itrc)=t(i,Jstr-1,k,nout,itrc)* &
539!^ & GRID(ng)%rmask(i,Jstr-1)
540!^
541 tl_t(i,jstr-1,k,nout,itrc)=tl_t(i,jstr-1,k,nout,itrc)* &
542 & grid(ng)%rmask(i,jstr-1)
543# endif
544 END IF
545 END DO
546 END DO
547 END IF
548 END IF
549!
550!-----------------------------------------------------------------------
551! Lateral boundary conditions at the northern edge.
552!-----------------------------------------------------------------------
553!
554 IF (domain(ng)%Northern_Edge(tile)) THEN
555!
556! Northern edge, implicit upstream radiation condition.
557!
558 IF (tl_lbc(inorth,istvar(itrc),ng)%radiation) THEN
559 IF (iic(ng).ne.0) THEN
560 DO k=1,n(ng)
561 DO i=istr,iend+1
562!^ grad(i,Jend+1)=t(i ,Jend+1,k,nstp,itrc)- &
563!^ & t(i-1,Jend+1,k,nstp,itrc)
564!^
565 tl_grad(i,jend+1)=0.0_r8
566 END DO
567 DO i=istr,iend
568 IF (lbc_apply(ng)%north(i)) THEN
569# if defined CELERITY_READ && defined FORWARD_READ
570 IF (tl_lbc(inorth,istvar(itrc),ng)%nudging) THEN
571 IF (lnudgetclm(itrc,ng)) THEN
572 obc_out=clima(ng)%Tnudgcof(i,jend+1,k,ic)
573 obc_in =obcfac(ng)*obc_out
574 ELSE
575 obc_out=tobc_out(itrc,ng,inorth)
576 obc_in =tobc_in(itrc,ng,inorth)
577 END IF
578 IF (boundary(ng)%t_north_Ce(i,k,itrc).lt. &
579 & 0.0_r8) THEN
580 tau=obc_in
581 ELSE
582 tau=obc_out
583 END IF
584 tau=tau*dt(ng)
585 END IF
586# ifdef RADIATION_2D
587 cx=boundary(ng)%t_north_Cx(i,k,itrc)
588# else
589 cx=0.0_r8
590# endif
591 ce=boundary(ng)%t_north_Ce(i,k,itrc)
592 cff=boundary(ng)%t_north_C2(i,k,itrc)
593# endif
594!^ t(i,Jend+1,k,nout,itrc)=(cff*t(i,Jend+1,k,nstp,itrc)+ &
595!^ & Ce *t(i,Jend ,k,nout,itrc)- &
596!^ & MAX(Cx,0.0_r8)* &
597!^ & grad(i ,Jend+1)- &
598!^ & MIN(Cx,0.0_r8)* &
599!^ & grad(i+1,Jend+1))/ &
600!^ & (cff+Ce)
601!^
602 tl_t(i,jend+1,k,nout,itrc)=(cff* &
603 & tl_t(i,jend+1,k,nstp, &
604 & itrc)+ &
605 & ce* &
606 & tl_t(i,jend ,k,nout, &
607 & itrc)- &
608 & max(cx,0.0_r8)* &
609 & tl_grad(i ,jend+1)- &
610 & min(cx,0.0_r8)* &
611 & tl_grad(i+1,jend+1))/ &
612 & (cff+ce)
613
614 IF (tl_lbc(inorth,istvar(itrc),ng)%nudging) THEN
615!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend+1,k,nout,itrc)+ &
616!^ & tau* &
617!^ & (BOUNDARY(ng)% &
618!^ & t_north(i,k,itrc)- &
619!^ & t(i,Jend+1,k,nstp,itrc))
620!^
621 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend+1,k,nout, &
622 & itrc)- &
623 & tau* &
624 & tl_t(i,jend+1,k,nstp, &
625 & itrc)
626 END IF
627# ifdef MASKING
628!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend+1,k,nout,itrc)* &
629!^ & GRID(ng)%rmask(i,Jend+1)
630!^
631 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend+1,k,nout, &
632 & itrc)* &
633 & grid(ng)%rmask(i,jend+1)
634# endif
635 END IF
636 END DO
637 END DO
638 END IF
639!
640! Northern edge, clamped boundary condition.
641!
642 ELSE IF (tl_lbc(inorth,istvar(itrc),ng)%clamped) THEN
643 DO k=1,n(ng)
644 DO i=istr,iend
645 IF (lbc_apply(ng)%north(i)) THEN
646!^ t(i,Jend+1,k,nout,itrc)=BOUNDARY(ng)%t_north(i,k,itrc)
647!^
648 tl_t(i,jend+1,k,nout,itrc)=boundary(ng)% &
649 & tl_t_north(i,k,itrc)
650# ifdef MASKING
651!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend+1,k,nout,itrc)* &
652!^ & GRID(ng)%rmask(i,Jend+1)
653!^
654 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend+1,k,nout,itrc)* &
655 & grid(ng)%rmask(i,jend+1)
656# endif
657 END IF
658 END DO
659 END DO
660!
661! Northern edge, gradient boundary condition.
662!
663 ELSE IF (tl_lbc(inorth,istvar(itrc),ng)%gradient) THEN
664 DO k=1,n(ng)
665 DO i=istr,iend
666 IF (lbc_apply(ng)%north(i)) THEN
667!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend,k,nout,itrc)
668!^
669 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend,k,nout,itrc)
670# ifdef MASKING
671!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend+1,k,nout,itrc)* &
672!^ & GRID(ng)%rmask(i,Jend+1)
673!^
674 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend+1,k,nout,itrc)* &
675 & grid(ng)%rmask(i,jend+1)
676# endif
677 END IF
678 END DO
679 END DO
680!
681! Northern edge, closed boundary condition.
682!
683 ELSE IF (tl_lbc(inorth,istvar(itrc),ng)%closed) THEN
684 DO k=1,n(ng)
685 DO i=istr,iend
686 IF (lbc_apply(ng)%north(i)) THEN
687!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend,k,nout,itrc)
688!^
689 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend,k,nout,itrc)
690# ifdef MASKING
691!^ t(i,Jend+1,k,nout,itrc)=t(i,Jend+1,k,nout,itrc)* &
692!^ & GRID(ng)%rmask(i,Jend+1)
693!^
694 tl_t(i,jend+1,k,nout,itrc)=tl_t(i,jend+1,k,nout,itrc)* &
695 & grid(ng)%rmask(i,jend+1)
696# endif
697 END IF
698 END DO
699 END DO
700 END IF
701 END IF
702!
703!-----------------------------------------------------------------------
704! Boundary corners.
705!-----------------------------------------------------------------------
706!
707 IF (.not.(ewperiodic(ng).or.nsperiodic(ng))) THEN
708 IF (domain(ng)%SouthWest_Corner(tile)) THEN
709 IF (lbc_apply(ng)%south(istr-1).and. &
710 & lbc_apply(ng)%west (jstr-1)) THEN
711 DO k=1,n(ng)
712!^ t(Istr-1,Jstr-1,k,nout,itrc)=0.5_r8* &
713!^ & (t(Istr ,Jstr-1,k,nout, &
714!^ & itrc)+ &
715!^ & t(Istr-1,Jstr ,k,nout, &
716!^ & itrc))
717!^
718 tl_t(istr-1,jstr-1,k,nout,itrc)=0.5_r8* &
719 & (tl_t(istr ,jstr-1,k, &
720 & nout,itrc)+ &
721 & tl_t(istr-1,jstr ,k, &
722 & nout,itrc))
723 END DO
724 END IF
725 END IF
726 IF (domain(ng)%SouthEast_Corner(tile)) THEN
727 IF (lbc_apply(ng)%south(iend+1).and. &
728 & lbc_apply(ng)%east (jstr-1)) THEN
729 DO k=1,n(ng)
730!^ t(Iend+1,Jstr-1,k,nout,itrc)=0.5_r8* &
731!^ & (t(Iend ,Jstr-1,k,nout, &
732!^ & itrc)+ &
733!^ & t(Iend+1,Jstr ,k,nout, &
734!^ & itrc))
735!^
736 tl_t(iend+1,jstr-1,k,nout,itrc)=0.5_r8* &
737 & (tl_t(iend ,jstr-1,k, &
738 & nout,itrc)+ &
739 & tl_t(iend+1,jstr ,k, &
740 & nout,itrc))
741 END DO
742 END IF
743 END IF
744 IF (domain(ng)%NorthWest_Corner(tile)) THEN
745 IF (lbc_apply(ng)%north(istr-1).and. &
746 & lbc_apply(ng)%west (jend+1)) THEN
747 DO k=1,n(ng)
748!^ t(Istr-1,Jend+1,k,nout,itrc)=0.5_r8* &
749!^ & (t(Istr-1,Jend ,k,nout, &
750!^ & itrc)+ &
751!^ & t(Istr ,Jend+1,k,nout, &
752!^ & itrc))
753!^
754 tl_t(istr-1,jend+1,k,nout,itrc)=0.5_r8* &
755 & (tl_t(istr-1,jend ,k, &
756 & nout,itrc)+ &
757 & tl_t(istr ,jend+1,k, &
758 & nout,itrc))
759 END DO
760 END IF
761 END IF
762 IF (domain(ng)%NorthEast_Corner(tile)) THEN
763 IF (lbc_apply(ng)%north(iend+1).and. &
764 & lbc_apply(ng)%east (jend+1)) THEN
765 DO k=1,n(ng)
766!^ t(Iend+1,Jend+1,k,nout,itrc)=0.5_r8* &
767!^ & (t(Iend+1,Jend ,k,nout, &
768!^ & itrc)+ &
769!^ & t(Iend ,Jend+1,k,nout, &
770!^ & itrc))
771!^
772 tl_t(iend+1,jend+1,k,nout,itrc)=0.5_r8* &
773 & (tl_t(iend+1,jend ,k, &
774 & nout,itrc)+ &
775 & tl_t(iend ,jend+1,k, &
776 & nout,itrc))
777 END DO
778 END IF
779 END IF
780 END IF
781
782 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 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_param::n, mod_scalars::nsperiodic, mod_scalars::obcfac, mod_param::tl_lbc, mod_scalars::tobc_in, and mod_scalars::tobc_out.

Referenced by rp_ini_fields_mod::rp_ini_fields_tile(), rp_pre_step3d_mod::rp_pre_step3d_tile(), rp_step3d_t_mod::rp_step3d_t_tile(), and rp_t3dbc().

Here is the caller graph for this function: