ROMS
Loading...
Searching...
No Matches
conv_bry2d.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#if defined NONLINEAR && 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! A 2D boundary state variable to diffuse !
40! !
41! On Output: !
42! !
43! A Convolved 2D boundary state variable !
44! !
45! Routines: !
46! !
47! conv_r2d_bry_tile Nonlinear 2D boundary convolution at RHO-points!
48! conv_u2d_bry_tile Nonlinear 2D boundary convolution at U-points !
49! conv_v2d_bry_tile Nonlinear 2D boundary convolution at V-points !
50! !
51!=======================================================================
52!
53 implicit none
54!
55 PUBLIC
56!
57 CONTAINS
58!
59!***********************************************************************
60 SUBROUTINE conv_r2d_bry_tile (ng, tile, model, boundary, &
61 & edge, LBij, UBij, &
62 & LBi, UBi, LBj, UBj, &
63 & IminS, ImaxS, JminS, JmaxS, &
64 & Nghost, NHsteps, DTsizeH, &
65 & Kh, &
66 & pm, pn, pmon_u, pnom_v, &
67# ifdef MASKING
68 & rmask, umask, vmask, &
69# endif
70 & A)
71!***********************************************************************
72!
73 USE mod_param
74 USE mod_scalars
75!
77# ifdef DISTRIBUTE
79# endif
80!
81! Imported variable declarations.
82!
83 integer, intent(in) :: ng, tile, model, boundary
84 integer, intent(in) :: edge(4)
85 integer, intent(in) :: LBij, UBij
86 integer, intent(in) :: LBi, UBi, LBj, UBj
87 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
88 integer, intent(in) :: Nghost, NHsteps
89
90 real(r8), intent(in) :: DTsizeH
91!
92# ifdef ASSUMED_SHAPE
93 real(r8), intent(in) :: pm(LBi:,LBj:)
94 real(r8), intent(in) :: pn(LBi:,LBj:)
95 real(r8), intent(in) :: pmon_u(LBi:,LBj:)
96 real(r8), intent(in) :: pnom_v(LBi:,LBj:)
97# ifdef MASKING
98 real(r8), intent(in) :: rmask(LBi:,LBj:)
99 real(r8), intent(in) :: umask(LBi:,LBj:)
100 real(r8), intent(in) :: vmask(LBi:,LBj:)
101# endif
102 real(r8), intent(in) :: Kh(LBi:,LBj:)
103 real(r8), intent(inout) :: A(LBij:)
104# else
105 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
106 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
107 real(r8), intent(in) :: pmon_u(LBi:UBi,LBj:UBj)
108 real(r8), intent(in) :: pnom_v(LBi:UBi,LBj:UBj)
109# ifdef MASKING
110 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
111 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
112 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
113# endif
114 real(r8), intent(in) :: Kh(LBi:UBi,LBj:UBj)
115 real(r8), intent(inout) :: A(LBij:UBij)
116# endif
117!
118! Local variable declarations.
119!
120 logical, dimension(4) :: Lconvolve
121
122 integer :: Nnew, Nold, Nsav, i, j, step
123
124 real(r8), dimension(LBij:UBij,2) :: Awrk
125
126 real(r8), dimension(JminS:JmaxS) :: FE
127 real(r8), dimension(IminS:ImaxS) :: FX
128 real(r8), dimension(LBij:UBij) :: Hfac
129
130# include "set_bounds.h"
131!
132!-----------------------------------------------------------------------
133! Space convolution of the diffusion equation for a 1D state variable
134! at RHO-points.
135!-----------------------------------------------------------------------
136!
137 lconvolve(iwest )=domain(ng)%Western_Edge (tile)
138 lconvolve(ieast )=domain(ng)%Eastern_Edge (tile)
139 lconvolve(isouth)=domain(ng)%Southern_Edge(tile)
140 lconvolve(inorth)=domain(ng)%Northern_Edge(tile)
141!
142! Compute metrics factor.
143!
144 IF (lconvolve(boundary)) THEN
145 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
146 i=edge(boundary)
147 DO j=jstr,jend
148 hfac(j)=dtsizeh*pm(i,j)*pn(i,j)
149 END DO
150 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
151 j=edge(boundary)
152 DO i=istr,iend
153 hfac(i)=dtsizeh*pm(i,j)*pn(i,j)
154 END DO
155 END IF
156 END IF
157!
158! Set integration indices and initial conditions.
159!
160 nold=1
161 nnew=2
162
163 CALL bc_r2d_bry_tile (ng, tile, boundary, &
164 & lbij, ubij, &
165 & a)
166# ifdef DISTRIBUTE
167 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
168 & lbij, ubij, &
169 & nghost, &
170 & ewperiodic(ng), nsperiodic(ng), &
171 & a)
172# endif
173 IF (lconvolve(boundary)) THEN
174 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
175 DO j=jstr-1,jend+1
176 awrk(j,nold)=a(j)
177 END DO
178 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
179 DO i=istr-1,iend+1
180 awrk(i,nold)=a(i)
181 END DO
182 END IF
183 END IF
184!
185!-----------------------------------------------------------------------
186! Integrate horizontal diffusion terms.
187!-----------------------------------------------------------------------
188!
189 DO step=1,nhsteps
190!
191! Compute XI- and ETA-components of diffusive flux.
192!
193 IF (lconvolve(boundary)) THEN
194 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
195 i=edge(boundary)
196 DO j=jstr,jend+1
197 fe(j)=pnom_v(i,j)*0.5_r8*(kh(i,j-1)+kh(i,j))* &
198 & (awrk(j ,nold)- &
199 & awrk(j-1,nold))
200# ifdef MASKING
201 fe(j)=fe(j)*vmask(i,j)
202# endif
203 END DO
204 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
205 j=edge(boundary)
206 DO i=istr,iend+1
207 fx(i)=pmon_u(i,j)*0.5_r8*(kh(i-1,j)+kh(i,j))* &
208 & (awrk(i ,nold)- &
209 & awrk(i-1,nold))
210# ifdef MASKING
211 fx(i)=fx(i)*umask(i,j)
212# endif
213 END DO
214 END IF
215 END IF
216!
217! Time-step horizontal diffusion terms.
218!
219 IF (lconvolve(boundary)) THEN
220 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
221 DO j=jstr,jend
222 awrk(j,nnew)=awrk(j,nold)+ &
223 & hfac(j)* &
224 & (fe(j+1)-fe(j))
225 END DO
226 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
227 DO i=istr,iend
228 awrk(i,nnew)=awrk(i,nold)+ &
229 & hfac(i)* &
230 & (fx(i+1)-fx(i))
231 END DO
232 END IF
233 END IF
234!
235! Apply boundary conditions. If applicable, exchange boundary data.
236!
237 CALL bc_r2d_bry_tile (ng, tile, boundary, &
238 & lbij, ubij, &
239 & awrk(:,nnew))
240# ifdef DISTRIBUTE
241 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
242 & lbij, ubij, &
243 & nghost, &
244 & ewperiodic(ng), nsperiodic(ng), &
245 & awrk(:,nnew))
246# endif
247!
248! Update integration indices.
249!
250 nsav=nold
251 nold=nnew
252 nnew=nsav
253 END DO
254!
255!-----------------------------------------------------------------------
256! Load convolved solution.
257!-----------------------------------------------------------------------
258!
259 IF (lconvolve(boundary)) THEN
260 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
261 DO j=jstr,jend
262 a(j)=awrk(j,nold)
263 END DO
264 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
265 DO i=istr,iend
266 a(i)=awrk(i,nold)
267 END DO
268 END IF
269 END IF
270 CALL bc_r2d_bry_tile (ng, tile, boundary, &
271 & lbij, ubij, &
272 & a)
273# ifdef DISTRIBUTE
274 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
275 & lbij, ubij, &
276 & nghost, &
277 & ewperiodic(ng), nsperiodic(ng), &
278 & a)
279# endif
280
281 RETURN
282 END SUBROUTINE conv_r2d_bry_tile
283
284!
285!***********************************************************************
286 SUBROUTINE conv_u2d_bry_tile (ng, tile, model, boundary, &
287 & edge, LBij, UBij, &
288 & LBi, UBi, LBj, UBj, &
289 & IminS, ImaxS, JminS, JmaxS, &
290 & Nghost, NHsteps, DTsizeH, &
291 & Kh, &
292 & pm, pn, pmon_r, pnom_p, &
293# ifdef MASKING
294 & umask, pmask, &
295# endif
296 & A)
297!***********************************************************************
298!
299 USE mod_param
300 USE mod_scalars
301!
303# ifdef DISTRIBUTE
305# endif
306!
307! Imported variable declarations.
308!
309 integer, intent(in) :: ng, tile, model, boundary
310 integer, intent(in) :: edge(4)
311 integer, intent(in) :: LBij, UBij
312 integer, intent(in) :: LBi, UBi, LBj, UBj
313 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
314 integer, intent(in) :: Nghost, NHsteps
315
316 real(r8), intent(in) :: DTsizeH
317!
318# ifdef ASSUMED_SHAPE
319 real(r8), intent(in) :: pm(LBi:,LBj:)
320 real(r8), intent(in) :: pn(LBi:,LBj:)
321 real(r8), intent(in) :: pmon_r(LBi:,LBj:)
322 real(r8), intent(in) :: pnom_p(LBi:,LBj:)
323# ifdef MASKING
324 real(r8), intent(in) :: umask(LBi:,LBj:)
325 real(r8), intent(in) :: pmask(LBi:,LBj:)
326# endif
327 real(r8), intent(in) :: Kh(LBi:,LBj:)
328 real(r8), intent(inout) :: A(LBij:)
329# else
330 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
331 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
332 real(r8), intent(in) :: pmon_r(LBi:UBi,LBj:UBj)
333 real(r8), intent(in) :: pnom_p(LBi:UBi,LBj:UBj)
334# ifdef MASKING
335 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
336 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
337# endif
338 real(r8), intent(in) :: Kh(LBi:UBi,LBj:UBj)
339 real(r8), intent(inout) :: A(LBij:UBij)
340# endif
341!
342! Local variable declarations.
343!
344 logical, dimension(4) :: Lconvolve
345
346 integer :: Nnew, Nold, Nsav, i, j, step
347
348 real(r8) :: cff
349
350 real(r8), dimension(LBij:UBij,2) :: Awrk
351
352 real(r8), dimension(JminS:JmaxS) :: FE
353 real(r8), dimension(IminS:ImaxS) :: FX
354 real(r8), dimension(LBij:UBij) :: Hfac
355
356# include "set_bounds.h"
357!
358!-----------------------------------------------------------------------
359! Space convolution of the diffusion equation for a 1D state variable
360! at U-points.
361!-----------------------------------------------------------------------
362!
363 lconvolve(iwest )=domain(ng)%Western_Edge (tile)
364 lconvolve(ieast )=domain(ng)%Eastern_Edge (tile)
365 lconvolve(isouth)=domain(ng)%Southern_Edge(tile)
366 lconvolve(inorth)=domain(ng)%Northern_Edge(tile)
367!
368! Compute metrics factor.
369!
370 cff=dtsizeh*0.25_r8
371 IF (lconvolve(boundary)) THEN
372 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
373 i=edge(boundary)
374 DO j=jstr,jend
375 hfac(j)=cff*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))
376 END DO
377 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
378 j=edge(boundary)
379 DO i=istru,iend
380 hfac(i)=cff*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))
381 END DO
382 END IF
383 END IF
384!
385! Set integration indices and initial conditions.
386!
387 nold=1
388 nnew=2
389
390 CALL bc_u2d_bry_tile (ng, tile, boundary, &
391 & lbij, ubij, &
392 & a)
393# ifdef DISTRIBUTE
394 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
395 & lbij, ubij, &
396 & nghost, &
397 & ewperiodic(ng), nsperiodic(ng), &
398 & a)
399# endif
400 IF (lconvolve(boundary)) THEN
401 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
402 DO j=jstr-1,jend+1
403 awrk(j,nold)=a(j)
404 END DO
405 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
406 DO i=istru-1,iend+1
407 awrk(i,nold)=a(i)
408 END DO
409 END IF
410 END IF
411!
412!-----------------------------------------------------------------------
413! Integrate horizontal diffusion terms.
414!-----------------------------------------------------------------------
415!
416 DO step=1,nhsteps
417!
418! Compute XI- and ETA-components of diffusive flux.
419!
420 IF (lconvolve(boundary)) THEN
421 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
422 i=edge(boundary)
423 DO j=jstr,jend+1
424 fe(j)=pnom_p(i,j)*0.25_r8*(kh(i-1,j )+kh(i,j )+ &
425 & kh(i-1,j-1)+kh(i,j-1))* &
426 & (awrk(j ,nold)- &
427 & awrk(j-1,nold))
428# ifdef MASKING
429 fe(j)=fe(j)*pmask(i,j)
430# endif
431 END DO
432 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
433 j=edge(boundary)
434 DO i=istru-1,iend
435 fx(i)=pmon_r(i,j)*kh(i,j)* &
436 & (awrk(i+1,nold)- &
437 & awrk(i ,nold))
438 END DO
439 END IF
440 END IF
441!
442! Time-step horizontal diffusion terms.
443!
444 IF (lconvolve(boundary)) THEN
445 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
446 DO j=jstr,jend
447 awrk(j,nnew)=awrk(j,nold)+ &
448 & hfac(j)* &
449 & (fe(j+1)-fe(j))
450 END DO
451 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
452 DO i=istru,iend
453 awrk(i,nnew)=awrk(i,nold)+ &
454 & hfac(i)* &
455 & (fx(i)-fx(i-1))
456 END DO
457 END IF
458 END IF
459!
460! Apply boundary conditions. If applicable, exchange boundary data.
461!
462 CALL bc_u2d_bry_tile (ng, tile, boundary, &
463 & lbij, ubij, &
464 & awrk(:,nnew))
465# ifdef DISTRIBUTE
466 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
467 & lbij, ubij, &
468 & nghost, &
469 & ewperiodic(ng), nsperiodic(ng), &
470 & awrk(:,nnew))
471# endif
472!
473! Update integration indices.
474!
475 nsav=nold
476 nold=nnew
477 nnew=nsav
478 END DO
479!
480!-----------------------------------------------------------------------
481! Load convolved solution.
482!-----------------------------------------------------------------------
483!
484 IF (lconvolve(boundary)) THEN
485 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
486 DO j=jstr,jend
487 a(j)=awrk(j,nold)
488 END DO
489 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
490 DO i=istru,iend
491 a(i)=awrk(i,nold)
492 END DO
493 END IF
494 END IF
495 CALL bc_u2d_bry_tile (ng, tile, boundary, &
496 & lbij, ubij, &
497 & a)
498# ifdef DISTRIBUTE
499 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
500 & lbij, ubij, &
501 & nghost, &
502 & ewperiodic(ng), nsperiodic(ng), &
503 & a)
504# endif
505
506 RETURN
507 END SUBROUTINE conv_u2d_bry_tile
508
509!
510!***********************************************************************
511 SUBROUTINE conv_v2d_bry_tile (ng, tile, model, boundary, &
512 & edge, LBij, UBij, &
513 & LBi, UBi, LBj, UBj, &
514 & IminS, ImaxS, JminS, JmaxS, &
515 & Nghost, NHsteps, DTsizeH, &
516 & Kh, &
517 & pm, pn, pmon_p, pnom_r, &
518# ifdef MASKING
519 & vmask, pmask, &
520# endif
521 & A)
522!***********************************************************************
523!
524 USE mod_param
525 USE mod_scalars
526!
528# ifdef DISTRIBUTE
530# endif
531!
532! Imported variable declarations.
533!
534 integer, intent(in) :: ng, tile, model, boundary
535 integer, intent(in) :: edge(4)
536 integer, intent(in) :: LBij, UBij
537 integer, intent(in) :: LBi, UBi, LBj, UBj
538 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
539 integer, intent(in) :: Nghost, NHsteps
540
541 real(r8), intent(in) :: DTsizeH
542!
543# ifdef ASSUMED_SHAPE
544 real(r8), intent(in) :: pm(LBi:,LBj:)
545 real(r8), intent(in) :: pn(LBi:,LBj:)
546 real(r8), intent(in) :: pmon_p(LBi:,LBj:)
547 real(r8), intent(in) :: pnom_r(LBi:,LBj:)
548# ifdef MASKING
549 real(r8), intent(in) :: vmask(LBi:,LBj:)
550 real(r8), intent(in) :: pmask(LBi:,LBj:)
551# endif
552 real(r8), intent(in) :: Kh(LBi:,LBj:)
553 real(r8), intent(inout) :: A(LBij:)
554# else
555 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
556 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
557 real(r8), intent(in) :: pmon_p(LBi:UBi,LBj:UBj)
558 real(r8), intent(in) :: pnom_r(LBi:UBi,LBj:UBj)
559# ifdef MASKING
560 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
561 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
562# endif
563 real(r8), intent(in) :: Kh(LBi:UBi,LBj:UBj)
564 real(r8), intent(inout) :: A(LBij:UBij)
565# endif
566!
567! Local variable declarations.
568!
569 logical, dimension(4) :: Lconvolve
570
571 integer :: Nnew, Nold, Nsav, i, j, step
572
573 real(r8) :: cff
574
575 real(r8), dimension(LBij:UBij,2) :: Awrk
576
577 real(r8), dimension(JminS:JmaxS) :: FE
578 real(r8), dimension(IminS:ImaxS) :: FX
579 real(r8), dimension(LBij:UBij) :: Hfac
580
581# include "set_bounds.h"
582!
583!-----------------------------------------------------------------------
584! Space convolution of the diffusion equation for a 2D state variable
585! at RHO-points.
586!-----------------------------------------------------------------------
587!
588 lconvolve(iwest )=domain(ng)%Western_Edge (tile)
589 lconvolve(ieast )=domain(ng)%Eastern_Edge (tile)
590 lconvolve(isouth)=domain(ng)%Southern_Edge(tile)
591 lconvolve(inorth)=domain(ng)%Northern_Edge(tile)
592!
593! Compute metrics factor.
594!
595 cff=dtsizeh*0.25_r8
596 IF (lconvolve(boundary)) THEN
597 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
598 i=edge(boundary)
599 DO j=jstrv,jend
600 hfac(j)=cff*(pm(i,j-1)+pm(i,j))*(pn(i,j-1)+pn(i,j))
601 END DO
602 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
603 j=edge(boundary)
604 DO i=istr,iend
605 hfac(i)=cff*(pm(i,j-1)+pm(i,j))*(pn(i,j-1)+pn(i,j))
606 END DO
607 END IF
608 END IF
609!
610! Set integration indices and initial conditions.
611!
612 nold=1
613 nnew=2
614
615 CALL bc_v2d_bry_tile (ng, tile, boundary, &
616 & lbij, ubij, &
617 & a)
618# ifdef DISTRIBUTE
619 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
620 & lbij, ubij, &
621 & nghost, &
622 & ewperiodic(ng), nsperiodic(ng), &
623 & a)
624# endif
625 IF (lconvolve(boundary)) THEN
626 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
627 DO j=jstrv-1,jend+1
628 awrk(j,nold)=a(j)
629 END DO
630 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
631 DO i=istr-1,iend+1
632 awrk(i,nold)=a(i)
633 END DO
634 END IF
635 END IF
636!
637!-----------------------------------------------------------------------
638! Integrate horizontal diffusion terms.
639!-----------------------------------------------------------------------
640!
641 DO step=1,nhsteps
642!
643! Compute XI- and ETA-components of diffusive flux.
644!
645 IF (lconvolve(boundary)) THEN
646 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
647 i=edge(boundary)
648 DO j=jstrv-1,jend
649 fe(j)=pnom_r(i,j)*kh(i,j)* &
650 & (awrk(j+1,nold)- &
651 & awrk(j ,nold))
652 END DO
653 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
654 j=edge(boundary)
655 DO i=istr,iend+1
656 fx(i)=pmon_p(i,j)*0.25_r8*(kh(i-1,j )+kh(i,j )+ &
657 & kh(i-1,j-1)+kh(i,j-1))* &
658 & (awrk(i ,nold)- &
659 & awrk(i-1,nold))
660# ifdef MASKING
661 fx(i)=fx(i)*pmask(i,j)
662# endif
663 END DO
664 END IF
665 END IF
666!
667! Time-step horizontal diffusion terms.
668!
669 IF (lconvolve(boundary)) THEN
670 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
671 DO j=jstrv,jend
672 awrk(j,nnew)=awrk(j,nold)+ &
673 & hfac(j)* &
674 & (fe(j)-fe(j-1))
675 END DO
676 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
677 DO i=istr,iend
678 awrk(i,nnew)=awrk(i,nold)+ &
679 & hfac(i)* &
680 & (fx(i+1)-fx(i))
681 END DO
682 END IF
683 END IF
684!
685! Apply boundary conditions. If applicable, exchange boundary data.
686!
687 CALL bc_v2d_bry_tile (ng, tile, boundary, &
688 & lbij, ubij, &
689 & awrk(:,nnew))
690# ifdef DISTRIBUTE
691 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
692 & lbij, ubij, &
693 & nghost, &
694 & ewperiodic(ng), nsperiodic(ng), &
695 & awrk(:,nnew))
696# endif
697!
698! Update integration indices.
699!
700 nsav=nold
701 nold=nnew
702 nnew=nsav
703 END DO
704!
705!-----------------------------------------------------------------------
706! Load convolved solution.
707!-----------------------------------------------------------------------
708!
709 IF (lconvolve(boundary)) THEN
710 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
711 DO j=jstrv,jend
712 a(j)=awrk(j,nold)
713 END DO
714 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
715 DO i=istr,iend
716 a(i)=awrk(i,nold)
717 END DO
718 END IF
719 END IF
720 CALL bc_v2d_bry_tile (ng, tile, boundary, &
721 & lbij, ubij, &
722 & a)
723# ifdef DISTRIBUTE
724 CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
725 & lbij, ubij, &
726 & nghost, &
727 & ewperiodic(ng), nsperiodic(ng), &
728 & a)
729# endif
730
731 RETURN
732 END SUBROUTINE conv_v2d_bry_tile
733#endif
734 END MODULE 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
subroutine 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, a)
Definition conv_bry2d.F:71
subroutine 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, a)
Definition conv_bry2d.F:522
subroutine 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, a)
Definition conv_bry2d.F:297
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)