ROMS
Loading...
Searching...
No Matches
tl_conv_bry2d.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#if defined TANGENT && defined FOUR_DVAR && defined ADJUST_BOUNDARY
5!
6!git $Id$
7!================================================== Hernan G. Arango ===
8! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
9! Licensed under a MIT/X style license !
10! See License_ROMS.md !
11!=======================================================================
12! !
13! These routines applies the background error covariance to data !
14! assimilation fields via the space convolution of the diffusion !
15! equation (filter) for 3D state variables. The diffusion filter !
16! is solved using an explicit (inefficient) algorithm. !
17! !
18! For Gaussian (bell-shaped) correlations, the space convolution !
19! of the diffusion operator is an efficient way to estimate the !
20! finite domain error covariances. !
21! !
22! On Input: !
23! !
24! ng Nested grid number !
25! tile Tile partition !
26! model Calling model identifier !
27! boundary Boundary edge to convolve !
28! edge Boundary edges index !
29! LBij Lower bound MIN(I,J)-dimension !
30! LBij Lower bound MAX(I,J)-dimension !
31! LBi I-dimension Lower bound !
32! UBi I-dimension Upper bound !
33! LBj J-dimension Lower bound !
34! UBj J-dimension Upper bound !
35! Nghost Number of ghost points !
36! NHsteps Number of horizontal diffusion integration steps !
37! DTsizeH Horizontal diffusion pseudo time-step size !
38! Kh Horizontal diffusion coefficients !
39! tl_A 2D boundary state variable to diffuse !
40! !
41! On Output: !
42! !
43! tl_A Convolved 2D boundary state variable !
44! !
45! Routines: !
46! !
47! tl_conv_r2d_bry_tile Tangent linear 2D boundary convolution at !
48! RHO-points !
49! tl_conv_u2d_bry_tile Tangent linear 2D boundary convolution at !
50! U-points !
51! tl_conv_v2d_bry_tile Tangent linear 2D boundary convolution at !
52! V-points !
53! !
54!=======================================================================
55!
56 implicit none
57!
58 PUBLIC
59!
60 CONTAINS
61!
62!***********************************************************************
63 SUBROUTINE tl_conv_r2d_bry_tile (ng, tile, model, boundary, &
64 & edge, LBij, UBij, &
65 & LBi, UBi, LBj, UBj, &
66 & IminS, ImaxS, JminS, JmaxS, &
67 & Nghost, NHsteps, DTsizeH, &
68 & Kh, &
69 & pm, pn, pmon_u, pnom_v, &
70# ifdef MASKING
71 & rmask, umask, vmask, &
72# endif
73 & tl_A)
74!***********************************************************************
75!
76 USE mod_param
77 USE mod_scalars
78!
80# ifdef DISTRIBUTE
82# endif
83!
84! Imported variable declarations.
85!
86 integer, intent(in) :: ng, tile, model, boundary
87 integer, intent(in) :: edge(4)
88 integer, intent(in) :: LBij, UBij
89 integer, intent(in) :: LBi, UBi, LBj, UBj
90 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
91 integer, intent(in) :: Nghost, NHsteps
92
93 real(r8), intent(in) :: DTsizeH
94!
95# ifdef ASSUMED_SHAPE
96 real(r8), intent(in) :: pm(LBi:,LBj:)
97 real(r8), intent(in) :: pn(LBi:,LBj:)
98 real(r8), intent(in) :: pmon_u(LBi:,LBj:)
99 real(r8), intent(in) :: pnom_v(LBi:,LBj:)
100# ifdef MASKING
101 real(r8), intent(in) :: rmask(LBi:,LBj:)
102 real(r8), intent(in) :: umask(LBi:,LBj:)
103 real(r8), intent(in) :: vmask(LBi:,LBj:)
104# endif
105 real(r8), intent(in) :: Kh(LBi:,LBj:)
106 real(r8), intent(inout) :: tl_A(LBij:)
107# else
108 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
109 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
110 real(r8), intent(in) :: pmon_u(LBi:UBi,LBj:UBj)
111 real(r8), intent(in) :: pnom_v(LBi:UBi,LBj:UBj)
112# ifdef MASKING
113 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
114 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
115 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
116# endif
117 real(r8), intent(in) :: Kh(LBi:UBi,LBj:UBj)
118 real(r8), intent(inout) :: tl_A(LBij:UBij)
119# endif
120!
121! Local variable declarations.
122!
123 logical, dimension(4) :: Lconvolve
124
125 integer :: Nnew, Nold, Nsav, i, j, step
126
127 real(r8), dimension(LBij:UBij,2) :: tl_Awrk
128
129 real(r8), dimension(JminS:JmaxS) :: tl_FE
130 real(r8), dimension(IminS:ImaxS) :: tl_FX
131 real(r8), dimension(LBij:UBij) :: Hfac
132
133# include "set_bounds.h"
134!
135!-----------------------------------------------------------------------
136! Space convolution of the diffusion equation for a 1D state variable
137! at RHO-points.
138!-----------------------------------------------------------------------
139!
140 lconvolve(iwest )=domain(ng)%Western_Edge (tile)
141 lconvolve(ieast )=domain(ng)%Eastern_Edge (tile)
142 lconvolve(isouth)=domain(ng)%Southern_Edge(tile)
143 lconvolve(inorth)=domain(ng)%Northern_Edge(tile)
144!
145! Compute metrics factor.
146!
147 IF (lconvolve(boundary)) THEN
148 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
149 i=edge(boundary)
150 DO j=jstr,jend
151 hfac(j)=dtsizeh*pm(i,j)*pn(i,j)
152 END DO
153 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
154 j=edge(boundary)
155 DO i=istr,iend
156 hfac(i)=dtsizeh*pm(i,j)*pn(i,j)
157 END DO
158 END IF
159 END IF
160!
161! Set integration indices and initial conditions.
162!
163 nold=1
164 nnew=2
165
166!^ CALL bc_r2d_bry_tile (ng, tile, boundary, &
167!^ & LBij, UBij, &
168!^ & A)
169!^
170 CALL bc_r2d_bry_tile (ng, tile, boundary, &
171 & lbij, ubij, &
172 & tl_a)
173# ifdef DISTRIBUTE
174!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
175!^ & LBij, UBij, &
176!^ & Nghost, &
177!^ & EWperiodic(ng), NSperiodic(ng), &
178!^ & A)
179!^
180 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
181 & lbij, ubij, &
182 & nghost, &
183 & ewperiodic(ng), nsperiodic(ng), &
184 & tl_a)
185# endif
186 IF (lconvolve(boundary)) THEN
187 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
188 DO j=jstr-1,jend+1
189!^ Awrk(j,Nold)=A(j)
190!^
191 tl_awrk(j,nold)=tl_a(j)
192 END DO
193 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
194 DO i=istr-1,iend+1
195!^ Awrk(i,Nold)=A(i)
196!^
197 tl_awrk(i,nold)=tl_a(i)
198 END DO
199 END IF
200 END IF
201!
202!-----------------------------------------------------------------------
203! Integrate horizontal diffusion terms.
204!-----------------------------------------------------------------------
205!
206 DO step=1,nhsteps
207!
208! Compute XI- and ETA-components of diffusive flux.
209!
210 IF (lconvolve(boundary)) THEN
211 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
212 i=edge(boundary)
213 DO j=jstr,jend+1
214!^ FE(j)=pnom_v(i,j)*0.5_r8*(Kh(i,j-1)+Kh(i,j))* &
215!^ & (Awrk(j ,Nold)- &
216!^ & Awrk(j-1,Nold))
217!^
218 tl_fe(j)=pnom_v(i,j)*0.5_r8*(kh(i,j-1)+kh(i,j))* &
219 & (tl_awrk(j ,nold)- &
220 & tl_awrk(j-1,nold))
221# ifdef MASKING
222!^ FE(j)=FE(j)*vmask(i,j)
223!^
224 tl_fe(j)=tl_fe(j)*vmask(i,j)
225# endif
226 END DO
227 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
228 j=edge(boundary)
229 DO i=istr,iend+1
230!^ FX(i)=pmon_u(i,j)*0.5_r8*(Kh(i-1,j)+Kh(i,j))* &
231!^ & (Awrk(i ,Nold)- &
232!^ & Awrk(i-1,Nold))
233!^
234 tl_fx(i)=pmon_u(i,j)*0.5_r8*(kh(i-1,j)+kh(i,j))* &
235 & (tl_awrk(i ,nold)- &
236 & tl_awrk(i-1,nold))
237# ifdef MASKING
238!^ FX(i)=FX(i)*umask(i,j)
239!^
240 tl_fx(i)=tl_fx(i)*umask(i,j)
241# endif
242 END DO
243 END IF
244 END IF
245!
246! Time-step horizontal diffusion terms.
247!
248 IF (lconvolve(boundary)) THEN
249 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
250 DO j=jstr,jend
251!^ Awrk(j,Nnew)=Awrk(j,Nold)+ &
252!^ & Hfac(j)* &
253!^ & (FE(j+1)-FE(j))
254!^
255 tl_awrk(j,nnew)=tl_awrk(j,nold)+ &
256 & hfac(j)* &
257 & (tl_fe(j+1)-tl_fe(j))
258 END DO
259 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
260 DO i=istr,iend
261!^ Awrk(i,Nnew)=Awrk(i,Nold)+ &
262!^ & Hfac(i)* &
263!^ & (FX(i+1)-FX(i))
264!^
265 tl_awrk(i,nnew)=tl_awrk(i,nold)+ &
266 & hfac(i)* &
267 & (tl_fx(i+1)-tl_fx(i))
268 END DO
269 END IF
270 END IF
271!
272! Apply boundary conditions. If applicable, exchange boundary data.
273!
274!^ CALL bc_r2d_bry_tile (ng, tile, boundary, &
275!^ & LBij, UBij, &
276!^ & Awrk(:,Nnew))
277!^
278 CALL bc_r2d_bry_tile (ng, tile, boundary, &
279 & lbij, ubij, &
280 & tl_awrk(:,nnew))
281# ifdef DISTRIBUTE
282!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
283!^ & LBij, UBij, &
284!^ & Nghost, &
285!^ & EWperiodic(ng), NSperiodic(ng), &
286!^ & Awrk(:,Nnew))
287!^
288 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
289 & lbij, ubij, &
290 & nghost, &
291 & ewperiodic(ng), nsperiodic(ng), &
292 & tl_awrk(:,nnew))
293# endif
294!
295! Update integration indices.
296!
297 nsav=nold
298 nold=nnew
299 nnew=nsav
300 END DO
301!
302!-----------------------------------------------------------------------
303! Load convolved solution.
304!-----------------------------------------------------------------------
305!
306 IF (lconvolve(boundary)) THEN
307 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
308 DO j=jstr,jend
309!^ A(j)=Awrk(j,Nold)
310!^
311 tl_a(j)=tl_awrk(j,nold)
312 END DO
313 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
314 DO i=istr,iend
315!^ A(i)=Awrk(i,Nold)
316!^
317 tl_a(i)=tl_awrk(i,nold)
318 END DO
319 END IF
320 END IF
321!^ CALL bc_r2d_bry_tile (ng, tile, boundary, &
322!^ & LBij, UBij, &
323!^ & A)
324!^
325 CALL bc_r2d_bry_tile (ng, tile, boundary, &
326 & lbij, ubij, &
327 & tl_a)
328# ifdef DISTRIBUTE
329!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
330!^ & LBij, UBij, &
331!^ & Nghost, &
332!^ & EWperiodic(ng), NSperiodic(ng), &
333!^ & A)
334!^
335 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
336 & lbij, ubij, &
337 & nghost, &
338 & ewperiodic(ng), nsperiodic(ng), &
339 & tl_a)
340# endif
341
342 RETURN
343 END SUBROUTINE tl_conv_r2d_bry_tile
344
345!
346!***********************************************************************
347 SUBROUTINE tl_conv_u2d_bry_tile (ng, tile, model, boundary, &
348 & edge, LBij, UBij, &
349 & LBi, UBi, LBj, UBj, &
350 & IminS, ImaxS, JminS, JmaxS, &
351 & Nghost, NHsteps, DTsizeH, &
352 & Kh, &
353 & pm, pn, pmon_r, pnom_p, &
354# ifdef MASKING
355 & umask, pmask, &
356# endif
357 & tl_A)
358!***********************************************************************
359!
360 USE mod_param
361 USE mod_scalars
362!
364# ifdef DISTRIBUTE
366# endif
367!
368! Imported variable declarations.
369!
370 integer, intent(in) :: ng, tile, model, boundary
371 integer, intent(in) :: edge(4)
372 integer, intent(in) :: LBij, UBij
373 integer, intent(in) :: LBi, UBi, LBj, UBj
374 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
375 integer, intent(in) :: Nghost, NHsteps
376
377 real(r8), intent(in) :: DTsizeH
378!
379# ifdef ASSUMED_SHAPE
380 real(r8), intent(in) :: pm(LBi:,LBj:)
381 real(r8), intent(in) :: pn(LBi:,LBj:)
382 real(r8), intent(in) :: pmon_r(LBi:,LBj:)
383 real(r8), intent(in) :: pnom_p(LBi:,LBj:)
384# ifdef MASKING
385 real(r8), intent(in) :: umask(LBi:,LBj:)
386 real(r8), intent(in) :: pmask(LBi:,LBj:)
387# endif
388 real(r8), intent(in) :: Kh(LBi:,LBj:)
389 real(r8), intent(inout) :: tl_A(LBij:)
390# else
391 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
392 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
393 real(r8), intent(in) :: pmon_r(LBi:UBi,LBj:UBj)
394 real(r8), intent(in) :: pnom_p(LBi:UBi,LBj:UBj)
395# ifdef MASKING
396 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
397 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
398# endif
399 real(r8), intent(in) :: Kh(LBi:UBi,LBj:UBj)
400 real(r8), intent(inout) :: tl_A(LBij:UBij)
401# endif
402!
403! Local variable declarations.
404!
405 logical, dimension(4) :: Lconvolve
406
407 integer :: Nnew, Nold, Nsav, i, j, step
408
409 real(r8) :: cff
410
411 real(r8), dimension(LBij:UBij,2) :: tl_Awrk
412
413 real(r8), dimension(JminS:JmaxS) :: tl_FE
414 real(r8), dimension(IminS:ImaxS) :: tl_FX
415 real(r8), dimension(LBij:UBij) :: Hfac
416
417# include "set_bounds.h"
418!
419!-----------------------------------------------------------------------
420! Space convolution of the diffusion equation for a 1D state variable
421! at U-points.
422!-----------------------------------------------------------------------
423!
424 lconvolve(iwest )=domain(ng)%Western_Edge (tile)
425 lconvolve(ieast )=domain(ng)%Eastern_Edge (tile)
426 lconvolve(isouth)=domain(ng)%Southern_Edge(tile)
427 lconvolve(inorth)=domain(ng)%Northern_Edge(tile)
428!
429! Compute metrics factor.
430!
431 cff=dtsizeh*0.25_r8
432 IF (lconvolve(boundary)) THEN
433 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
434 i=edge(boundary)
435 DO j=jstr,jend
436 hfac(j)=cff*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))
437 END DO
438 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
439 j=edge(boundary)
440 DO i=istru,iend
441 hfac(i)=cff*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))
442 END DO
443 END IF
444 END IF
445!
446! Set integration indices and initial conditions.
447!
448 nold=1
449 nnew=2
450
451!^ CALL bc_u2d_bry_tile (ng, tile, boundary, &
452!^ & LBij, UBij, &
453!^ & A)
454!^
455 CALL bc_u2d_bry_tile (ng, tile, boundary, &
456 & lbij, ubij, &
457 & tl_a)
458# ifdef DISTRIBUTE
459!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
460!^ & LBij, UBij, &
461!^ & Nghost, &
462!^ & EWperiodic(ng), NSperiodic(ng), &
463!^ & A)
464!^
465 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
466 & lbij, ubij, &
467 & nghost, &
468 & ewperiodic(ng), nsperiodic(ng), &
469 & tl_a)
470# endif
471 IF (lconvolve(boundary)) THEN
472 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
473 DO j=jstr-1,jend+1
474!^ Awrk(j,Nold)=A(j)
475!^
476 tl_awrk(j,nold)=tl_a(j)
477 END DO
478 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
479 DO i=istru-1,iend+1
480!^ Awrk(i,Nold)=A(i)
481!^
482 tl_awrk(i,nold)=tl_a(i)
483 END DO
484 END IF
485 END IF
486!
487!-----------------------------------------------------------------------
488! Integrate horizontal diffusion terms.
489!-----------------------------------------------------------------------
490!
491 DO step=1,nhsteps
492!
493! Compute XI- and ETA-components of diffusive flux.
494!
495 IF (lconvolve(boundary)) THEN
496 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
497 i=edge(boundary)
498 DO j=jstr,jend+1
499!^ FE(j)=pnom_p(i,j)*0.25_r8*(Kh(i-1,j )+Kh(i,j )+ &
500!^ & Kh(i-1,j-1)+Kh(i,j-1))* &
501!^ & (Awrk(j ,Nold)- &
502!^ & Awrk(j-1,Nold))
503!^
504 tl_fe(j)=pnom_p(i,j)*0.25_r8*(kh(i-1,j )+kh(i,j )+ &
505 & kh(i-1,j-1)+kh(i,j-1))* &
506 & (tl_awrk(j ,nold)- &
507 & tl_awrk(j-1,nold))
508# ifdef MASKING
509!^ FE(j)=FE(j)*pmask(i,j)
510!^
511 tl_fe(j)=tl_fe(j)*pmask(i,j)
512# endif
513 END DO
514 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
515 j=edge(boundary)
516 DO i=istru-1,iend
517!^ FX(i)=pmon_r(i,j)*Kh(i,j)* &
518!^ & (Awrk(i+1,Nold)- &
519!^ & Awrk(i ,Nold))
520!^
521 tl_fx(i)=pmon_r(i,j)*kh(i,j)* &
522 & (tl_awrk(i+1,nold)- &
523 & tl_awrk(i ,nold))
524 END DO
525 END IF
526 END IF
527!
528! Time-step horizontal diffusion terms.
529!
530 IF (lconvolve(boundary)) THEN
531 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
532 DO j=jstr,jend
533!^ Awrk(j,Nnew)=Awrk(j,Nold)+ &
534!^ & Hfac(j)* &
535!^ & (FE(j+1)-FE(j))
536!^
537 tl_awrk(j,nnew)=tl_awrk(j,nold)+ &
538 & hfac(j)* &
539 & (tl_fe(j+1)-tl_fe(j))
540 END DO
541 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
542 DO i=istru,iend
543!^ Awrk(i,Nnew)=Awrk(i,Nold)+ &
544!^ & Hfac(i)* &
545!^ & (FX(i)-FX(i-1))
546!^
547 tl_awrk(i,nnew)=tl_awrk(i,nold)+ &
548 & hfac(i)* &
549 & (tl_fx(i)-tl_fx(i-1))
550 END DO
551 END IF
552 END IF
553!
554! Apply boundary conditions. If applicable, exchange boundary data.
555!
556!^ CALL bc_u2d_bry_tile (ng, tile, boundary, &
557!^ & LBij, UBij, &
558!^ & Awrk(:,Nnew))
559!^
560 CALL bc_u2d_bry_tile (ng, tile, boundary, &
561 & lbij, ubij, &
562 & tl_awrk(:,nnew))
563# ifdef DISTRIBUTE
564!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
565!^ & LBij, UBij, &
566!^ & Nghost, &
567!^ & EWperiodic(ng), NSperiodic(ng), &
568!^ & Awrk(:,Nnew))
569!^
570 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
571 & lbij, ubij, &
572 & nghost, &
573 & ewperiodic(ng), nsperiodic(ng), &
574 & tl_awrk(:,nnew))
575# endif
576!
577! Update integration indices.
578!
579 nsav=nold
580 nold=nnew
581 nnew=nsav
582 END DO
583!
584!-----------------------------------------------------------------------
585! Load convolved solution.
586!-----------------------------------------------------------------------
587!
588 IF (lconvolve(boundary)) THEN
589 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
590 DO j=jstr,jend
591!^ A(j)=Awrk(j,Nold)
592!^
593 tl_a(j)=tl_awrk(j,nold)
594 END DO
595 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
596 DO i=istru,iend
597!^ A(i)=Awrk(i,Nold)
598!^
599 tl_a(i)=tl_awrk(i,nold)
600 END DO
601 END IF
602 END IF
603!^ CALL bc_u2d_bry_tile (ng, tile, boundary, &
604!^ & LBij, UBij, &
605!^ & A)
606!^
607 CALL bc_u2d_bry_tile (ng, tile, boundary, &
608 & lbij, ubij, &
609 & tl_a)
610# ifdef DISTRIBUTE
611!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
612!^ & LBij, UBij, &
613!^ & Nghost, &
614!^ & EWperiodic(ng), NSperiodic(ng), &
615!^ & A)
616!^
617 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
618 & lbij, ubij, &
619 & nghost, ewperiodic(ng), nsperiodic(ng), &
620 & tl_a)
621# endif
622
623 RETURN
624 END SUBROUTINE tl_conv_u2d_bry_tile
625
626!
627!***********************************************************************
628 SUBROUTINE tl_conv_v2d_bry_tile (ng, tile, model, boundary, &
629 & edge, LBij, UBij, &
630 & LBi, UBi, LBj, UBj, &
631 & IminS, ImaxS, JminS, JmaxS, &
632 & Nghost, NHsteps, DTsizeH, &
633 & Kh, &
634 & pm, pn, pmon_p, pnom_r, &
635# ifdef MASKING
636 & vmask, pmask, &
637# endif
638 & tl_A)
639!***********************************************************************
640!
641 USE mod_param
642 USE mod_scalars
643!
645# ifdef DISTRIBUTE
647# endif
648!
649! Imported variable declarations.
650!
651 integer, intent(in) :: ng, tile, model, boundary
652 integer, intent(in) :: edge(4)
653 integer, intent(in) :: LBij, UBij
654 integer, intent(in) :: LBi, UBi, LBj, UBj
655 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
656 integer, intent(in) :: Nghost, NHsteps
657
658 real(r8), intent(in) :: DTsizeH
659!
660# ifdef ASSUMED_SHAPE
661 real(r8), intent(in) :: pm(LBi:,LBj:)
662 real(r8), intent(in) :: pn(LBi:,LBj:)
663 real(r8), intent(in) :: pmon_p(LBi:,LBj:)
664 real(r8), intent(in) :: pnom_r(LBi:,LBj:)
665# ifdef MASKING
666 real(r8), intent(in) :: vmask(LBi:,LBj:)
667 real(r8), intent(in) :: pmask(LBi:,LBj:)
668# endif
669 real(r8), intent(in) :: Kh(LBi:,LBj:)
670 real(r8), intent(inout) :: tl_A(LBij:)
671# else
672 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
673 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
674 real(r8), intent(in) :: pmon_p(LBi:UBi,LBj:UBj)
675 real(r8), intent(in) :: pnom_r(LBi:UBi,LBj:UBj)
676# ifdef MASKING
677 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
678 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
679# endif
680 real(r8), intent(in) :: Kh(LBi:UBi,LBj:UBj)
681 real(r8), intent(inout) :: tl_A(LBij:UBij)
682# endif
683!
684! Local variable declarations.
685!
686 logical, dimension(4) :: Lconvolve
687
688 integer :: Nnew, Nold, Nsav, i, j, step
689
690 real(r8) :: cff
691
692 real(r8), dimension(LBij:UBij,2) :: tl_Awrk
693
694 real(r8), dimension(JminS:JmaxS) :: tl_FE
695 real(r8), dimension(IminS:ImaxS) :: tl_FX
696 real(r8), dimension(LBij:UBij) :: Hfac
697
698# include "set_bounds.h"
699!
700!-----------------------------------------------------------------------
701! Space convolution of the diffusion equation for a 2D state variable
702! at RHO-points.
703!-----------------------------------------------------------------------
704!
705 lconvolve(iwest )=domain(ng)%Western_Edge (tile)
706 lconvolve(ieast )=domain(ng)%Eastern_Edge (tile)
707 lconvolve(isouth)=domain(ng)%Southern_Edge(tile)
708 lconvolve(inorth)=domain(ng)%Northern_Edge(tile)
709!
710! Compute metrics factor.
711!
712 cff=dtsizeh*0.25_r8
713 IF (lconvolve(boundary)) THEN
714 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
715 i=edge(boundary)
716 DO j=jstrv,jend
717 hfac(j)=cff*(pm(i,j-1)+pm(i,j))*(pn(i,j-1)+pn(i,j))
718 END DO
719 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
720 j=edge(boundary)
721 DO i=istr,iend
722 hfac(i)=cff*(pm(i,j-1)+pm(i,j))*(pn(i,j-1)+pn(i,j))
723 END DO
724 END IF
725 END IF
726!
727! Set integration indices and initial conditions.
728!
729 nold=1
730 nnew=2
731
732!^ CALL bc_v2d_bry_tile (ng, tile, boundary, &
733!^ & LBij, UBij, &
734!^ & A)
735!^
736 CALL bc_v2d_bry_tile (ng, tile, boundary, &
737 & lbij, ubij, &
738 & tl_a)
739# ifdef DISTRIBUTE
740!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
741!^ & LBij, UBij, &
742!^ & Nghost, &
743!^ & EWperiodic(ng), NSperiodic(ng), &
744!^ & A)
745!^
746 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
747 & lbij, ubij, &
748 & nghost, &
749 & ewperiodic(ng), nsperiodic(ng), &
750 & tl_a)
751# endif
752 IF (lconvolve(boundary)) THEN
753 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
754 DO j=jstrv-1,jend+1
755!^ Awrk(j,Nold)=A(j)
756!^
757 tl_awrk(j,nold)=tl_a(j)
758 END DO
759 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
760 DO i=istr-1,iend+1
761!^ Awrk(i,Nold)=A(i)
762!^
763 tl_awrk(i,nold)=tl_a(i)
764 END DO
765 END IF
766 END IF
767!
768!-----------------------------------------------------------------------
769! Integrate horizontal diffusion terms.
770!-----------------------------------------------------------------------
771!
772 DO step=1,nhsteps
773!
774! Compute XI- and ETA-components of diffusive flux.
775!
776 IF (lconvolve(boundary)) THEN
777 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
778 i=edge(boundary)
779 DO j=jstrv-1,jend
780!^ FE(j)=pnom_r(i,j)*Kh(i,j)* &
781!^ & (Awrk(j+1,Nold)- &
782!^ & Awrk(j ,Nold))
783!^
784 tl_fe(j)=pnom_r(i,j)*kh(i,j)* &
785 & (tl_awrk(j+1,nold)- &
786 & tl_awrk(j ,nold))
787 END DO
788 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
789 j=edge(boundary)
790 DO i=istr,iend+1
791!^ FX(i)=pmon_p(i,j)*0.25_r8*(Kh(i-1,j )+Kh(i,j )+ &
792!^ & Kh(i-1,j-1)+Kh(i,j-1))* &
793!^ & (Awrk(i ,Nold)- &
794!^ & Awrk(i-1,Nold))
795!^
796 tl_fx(i)=pmon_p(i,j)*0.25_r8*(kh(i-1,j )+kh(i,j )+ &
797 & kh(i-1,j-1)+kh(i,j-1))* &
798 & (tl_awrk(i ,nold)- &
799 & tl_awrk(i-1,nold))
800# ifdef MASKING
801!^ FX(i)=FX(i)*pmask(i,j)
802!^
803 tl_fx(i)=tl_fx(i)*pmask(i,j)
804# endif
805 END DO
806 END IF
807 END IF
808!
809! Time-step horizontal diffusion terms.
810!
811 IF (lconvolve(boundary)) THEN
812 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
813 DO j=jstrv,jend
814!^ Awrk(j,Nnew)=Awrk(j,Nold)+ &
815!^ & Hfac(j)* &
816!^ & (FE(j)-FE(j-1))
817!^
818 tl_awrk(j,nnew)=tl_awrk(j,nold)+ &
819 & hfac(j)* &
820 & (tl_fe(j)-tl_fe(j-1))
821 END DO
822 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
823 DO i=istr,iend
824!^ Awrk(i,Nnew)=Awrk(i,Nold)+ &
825!^ & Hfac(i)* &
826!^ & (FX(i+1)-FX(i))
827!^
828 tl_awrk(i,nnew)=tl_awrk(i,nold)+ &
829 & hfac(i)* &
830 & (tl_fx(i+1)-tl_fx(i))
831 END DO
832 END IF
833 END IF
834!
835! Apply boundary conditions. If applicable, exchange boundary data.
836!
837!^ CALL bc_v2d_bry_tile (ng, tile, boundary, &
838!^ & LBij, UBij, &
839!^ & Awrk(:,Nnew))
840!^
841 CALL bc_v2d_bry_tile (ng, tile, boundary, &
842 & lbij, ubij, &
843 & tl_awrk(:,nnew))
844# ifdef DISTRIBUTE
845!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
846!^ & LBij, UBij, &
847!^ & Nghost, &
848!^ & EWperiodic(ng), NSperiodic(ng), &
849!^ & Awrk(:,Nnew))
850!^
851 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
852 & lbij, ubij, &
853 & nghost, &
854 & ewperiodic(ng), nsperiodic(ng), &
855 & tl_awrk(:,nnew))
856# endif
857!
858! Update integration indices.
859!
860 nsav=nold
861 nold=nnew
862 nnew=nsav
863 END DO
864!
865!-----------------------------------------------------------------------
866! Load convolved solution.
867!-----------------------------------------------------------------------
868!
869 IF (lconvolve(boundary)) THEN
870 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
871 DO j=jstrv,jend
872!^ A(j)=Awrk(j,Nold)
873!^
874 tl_a(j)=tl_awrk(j,nold)
875 END DO
876 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
877 DO i=istr,iend
878!^ A(i)=Awrk(i,Nold)
879!^
880 tl_a(i)=tl_awrk(i,nold)
881 END DO
882 END IF
883 END IF
884!^ CALL bc_v2d_bry_tile (ng, tile, boundary, &
885!^ & LBij, UBij, &
886!^ & A)
887!^
888 CALL bc_v2d_bry_tile (ng, tile, boundary, &
889 & lbij, ubij, &
890 & tl_a)
891# ifdef DISTRIBUTE
892!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
893!^ & LBij, UBij, &
894!^ & Nghost, &
895!^ & EWperiodic(ng), NSperiodic(ng), &
896!^ & A)
897!^
898 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
899 & lbij, ubij, &
900 & nghost, &
901 & ewperiodic(ng), nsperiodic(ng), &
902 & tl_a)
903# endif
904
905 RETURN
906 END SUBROUTINE tl_conv_v2d_bry_tile
907#endif
908 END MODULE tl_conv_bry2d_mod
subroutine bc_r2d_bry_tile(ng, tile, boundary, lbij, ubij, a)
Definition bc_bry2d.F:30
subroutine bc_u2d_bry_tile(ng, tile, boundary, lbij, ubij, a)
Definition bc_bry2d.F:100
subroutine bc_v2d_bry_tile(ng, tile, boundary, lbij, ubij, a)
Definition bc_bry2d.F:170
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth
subroutine mp_exchange2d_bry(ng, tile, model, nvar, boundary, lbij, ubij, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine tl_conv_v2d_bry_tile(ng, tile, model, boundary, edge, lbij, ubij, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nghost, nhsteps, dtsizeh, kh, pm, pn, pmon_p, pnom_r, vmask, pmask, tl_a)
subroutine tl_conv_u2d_bry_tile(ng, tile, model, boundary, edge, lbij, ubij, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nghost, nhsteps, dtsizeh, kh, pm, pn, pmon_r, pnom_p, umask, pmask, tl_a)
subroutine tl_conv_r2d_bry_tile(ng, tile, model, boundary, edge, lbij, ubij, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nghost, nhsteps, dtsizeh, kh, pm, pn, pmon_u, pnom_v, rmask, umask, vmask, tl_a)