ROMS
Loading...
Searching...
No Matches
ad_conv_bry2d.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#if defined ADJOINT && 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! ad_A 2D boundary state variable to diffuse !
40! !
41! On Output: !
42! !
43! ad_A Convolved 2D boundary state variable !
44! !
45! Routines: !
46! !
47! ad_conv_r2d_bry_tile Tangent linear 2D boundary convolution at !
48! RHO-points !
49! ad_conv_u2d_bry_tile Tangent linear 2D boundary convolution at !
50! U-points !
51! ad_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 ad_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 & ad_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) :: ad_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) :: ad_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) :: adfac
128
129 real(r8), dimension(LBij:UBij,2) :: ad_Awrk
130
131 real(r8), dimension(JminS:JmaxS) :: ad_FE
132 real(r8), dimension(IminS:ImaxS) :: ad_FX
133 real(r8), dimension(LBij:UBij) :: Hfac
134
135# include "set_bounds.h"
136!
137!-----------------------------------------------------------------------
138! Initialize adjoint private variables.
139!-----------------------------------------------------------------------
140!
141 ad_awrk(lbij:ubij,1:2)=0.0_r8
142
143 ad_fe(jmins:jmaxs)=0.0_r8
144 ad_fx(imins:imaxs)=0.0_r8
145!
146!-----------------------------------------------------------------------
147! Adjoint space convolution of the diffusion equation for a 1D state
148! variable at RHO-points.
149!-----------------------------------------------------------------------
150!
151 lconvolve(iwest )=domain(ng)%Western_Edge (tile)
152 lconvolve(ieast )=domain(ng)%Eastern_Edge (tile)
153 lconvolve(isouth)=domain(ng)%Southern_Edge(tile)
154 lconvolve(inorth)=domain(ng)%Northern_Edge(tile)
155!
156! Set integration indices and initial conditions.
157!
158 nold=1
159 nnew=2
160!
161! Compute metrics factor.
162!
163 IF (lconvolve(boundary)) THEN
164 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
165 i=edge(boundary)
166 DO j=jstr,jend
167 hfac(j)=dtsizeh*pm(i,j)*pn(i,j)
168 END DO
169 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
170 j=edge(boundary)
171 DO i=istr,iend
172 hfac(i)=dtsizeh*pm(i,j)*pn(i,j)
173 END DO
174 END IF
175 END IF
176!
177!-----------------------------------------------------------------------
178! Adjoint of load convolved solution.
179!-----------------------------------------------------------------------
180!
181# ifdef DISTRIBUTE
182!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
183!^ & LBij, UBij, &
184!^ & Nghost, &
185!^ & EWperiodic(ng), NSperiodic(ng), &
186!^ & tl_A)
187!^
188 CALL ad_mp_exchange2d_bry (ng, tile, model, 1, boundary, &
189 & lbij, ubij, &
190 & nghost, &
191 & ewperiodic(ng), nsperiodic(ng), &
192 & ad_a)
193# endif
194!^ CALL bc_r2d_bry_tile (ng, tile, boundary, &
195!^ & LBij, UBij, &
196!^ & tl_A)
197!^
198 CALL ad_bc_r2d_bry_tile (ng, tile, boundary, &
199 & lbij, ubij, &
200 & ad_a)
201 IF (lconvolve(boundary)) THEN
202 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
203 DO j=jstr,jend
204!^ tl_A(j)=tl_Awrk(j,Nold)
205!^
206 ad_awrk(j,nold)=ad_awrk(j,nold)+ad_a(j)
207 ad_a(j)=0.0_r8
208 END DO
209 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
210 DO i=istr,iend
211!^ tl_A(i)=tl_Awrk(i,Nold)
212!^
213 ad_awrk(i,nold)=ad_awrk(i,nold)+ad_a(i)
214 ad_a(i)=0.0_r8
215 END DO
216 END IF
217 END IF
218!
219!-----------------------------------------------------------------------
220! Integrate adjoint horizontal diffusion terms.
221!-----------------------------------------------------------------------
222!
223 DO step=1,nhsteps
224!
225! Update integration indices.
226!
227 nsav=nnew
228 nnew=nold
229 nold=nsav
230!
231! Apply adjoint boundary conditions. If applicable, exchange boundary
232! data.
233!
234# ifdef DISTRIBUTE
235!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
236!^ & LBij, UBij, &
237!^ & Nghost, &
238!^ & EWperiodic(ng), NSperiodic(ng), &
239!^ & tl_Awrk(:,Nnew))
240!^
241 CALL ad_mp_exchange2d_bry (ng, tile, model, 1, boundary, &
242 & lbij, ubij, &
243 & nghost, &
244 & ewperiodic(ng), nsperiodic(ng), &
245 & ad_awrk(:,nnew))
246# endif
247!^ CALL bc_r2d_bry_tile (ng, tile, boundary, &
248!^ & LBij, UBij, &
249!^ & tl_Awrk(:,Nnew))
250!^
251 CALL ad_bc_r2d_bry_tile (ng, tile, boundary, &
252 & lbij, ubij, &
253 & ad_awrk(:,nnew))
254!
255! Time-step adjoint horizontal diffusion terms.
256!
257 IF (lconvolve(boundary)) THEN
258 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
259 DO j=jstr,jend
260!^ tl_Awrk(j,Nnew)=tl_Awrk(j,Nold)+ &
261!^ & Hfac(j)* &
262!^ & (tl_FE(j+1)-tl_FE(j))
263!^
264 adfac=hfac(j)*ad_awrk(j,nnew)
265 ad_fe(j )=ad_fe(j )-adfac
266 ad_fe(j+1)=ad_fe(j+1)+adfac
267 ad_awrk(j,nold)=ad_awrk(j,nold)+ &
268 & ad_awrk(j,nnew)
269 ad_awrk(j,nnew)=0.0_r8
270 END DO
271 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
272 DO i=istr,iend
273!^ tl_Awrk(i,Nnew)=tl_Awrk(i,Nold)+ &
274!^ & Hfac(i)* &
275!^ & (tl_FX(i+1)-tl_FX(i))
276!^
277 adfac=hfac(i)*ad_awrk(i,nnew)
278 ad_fx(i )=ad_fx(i )-adfac
279 ad_fx(i+1)=ad_fx(i+1)+adfac
280 ad_awrk(i,nold)=ad_awrk(i,nold)+ &
281 & ad_awrk(i,nnew)
282 ad_awrk(i,nnew)=0.0_r8
283 END DO
284 END IF
285 END IF
286!
287! Compute XI- and ETA-components of adjoint diffusive flux.
288!
289 IF (lconvolve(boundary)) THEN
290 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
291 i=edge(boundary)
292 DO j=jstr,jend+1
293# ifdef MASKING
294!^ tl_FE(j)=tl_FE(j)*vmask(i,j)
295!^
296 ad_fe(j)=ad_fe(j)*vmask(i,j)
297# endif
298!^ tl_FE(j)=pnom_v(i,j)*0.5_r8*(Kh(i,j-1)+Kh(i,j))* &
299!^ & (tl_Awrk(j ,Nold)- &
300!^ & tl_Awrk(j-1,Nold))
301!^
302 adfac=pnom_v(i,j)*0.5_r8*(kh(i,j-1)+kh(i,j))*ad_fe(j)
303 ad_awrk(j-1,nold)=ad_awrk(j-1,nold)-adfac
304 ad_awrk(j ,nold)=ad_awrk(j ,nold)+adfac
305 ad_fe(j)=0.0_r8
306 END DO
307 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
308 j=edge(boundary)
309 DO i=istr,iend+1
310# ifdef MASKING
311!^ tl_FX(i)=tl_FX(i)*umask(i,j)
312!^
313 ad_fx(i)=ad_fx(i)*umask(i,j)
314# endif
315!^ tl_FX(i)=pmon_u(i,j)*0.5_r8*(Kh(i-1,j)+Kh(i,j))* &
316!^ & (tl_Awrk(i ,Nold)- &
317!^ & tl_Awrk(i-1,Nold))
318!^
319 adfac=pmon_u(i,j)*0.5_r8*(kh(i-1,j)+kh(i,j))*ad_fx(i)
320 ad_awrk(i-1,nold)=ad_awrk(i-1,nold)-adfac
321 ad_awrk(i ,nold)=ad_awrk(i ,nold)+adfac
322 ad_fx(i)=0.0_r8
323 END DO
324 END IF
325 END IF
326 END DO
327!
328! Set adjoint initial conditions.
329!
330 IF (lconvolve(boundary)) THEN
331 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
332 DO j=jstr-1,jend+1
333!^ tl_Awrk(j,Nold)=tl_A(j)
334!^
335 ad_a(j)=ad_a(j)+ad_awrk(j,nold)
336 ad_awrk(j,nold)=0.0_r8
337 END DO
338 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
339 DO i=istr-1,iend+1
340!^ tl_Awrk(i,Nold)=tl_A(i)
341!^
342 ad_a(i)=ad_a(i)+ad_awrk(i,nold)
343 ad_awrk(i,nold)=0.0_r8
344 END DO
345 END IF
346 END IF
347# ifdef DISTRIBUTE
348!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
349!^ & LBij, UBij, &
350!^ & Nghost, &
351!^ & EWperiodic(ng), NSperiodic(ng), &
352!^ & tl_A)
353!^
354 CALL ad_mp_exchange2d_bry (ng, tile, model, 1, boundary, &
355 & lbij, ubij, &
356 & nghost, &
357 & ewperiodic(ng), nsperiodic(ng), &
358 & ad_a)
359# endif
360!^ CALL bc_r2d_bry_tile (ng, tile, boundary, &
361!^ & LBij, UBij, &
362!^ & tl_A)
363!^
364 CALL ad_bc_r2d_bry_tile (ng, tile, boundary, &
365 & lbij, ubij, &
366 & ad_a)
367
368 RETURN
369 END SUBROUTINE ad_conv_r2d_bry_tile
370
371!
372!***********************************************************************
373 SUBROUTINE ad_conv_u2d_bry_tile (ng, tile, model, boundary, &
374 & edge, LBij, UBij, &
375 & LBi, UBi, LBj, UBj, &
376 & IminS, ImaxS, JminS, JmaxS, &
377 & Nghost, NHsteps, DTsizeH, &
378 & Kh, &
379 & pm, pn, pmon_r, pnom_p, &
380# ifdef MASKING
381 & umask, pmask, &
382# endif
383 & ad_A)
384!***********************************************************************
385!
386 USE mod_param
387 USE mod_scalars
388!
390# ifdef DISTRIBUTE
392# endif
393!
394! Imported variable declarations.
395!
396 integer, intent(in) :: ng, tile, model, boundary
397 integer, intent(in) :: edge(4)
398 integer, intent(in) :: LBij, UBij
399 integer, intent(in) :: LBi, UBi, LBj, UBj
400 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
401 integer, intent(in) :: Nghost, NHsteps
402
403 real(r8), intent(in) :: DTsizeH
404!
405# ifdef ASSUMED_SHAPE
406 real(r8), intent(in) :: pm(LBi:,LBj:)
407 real(r8), intent(in) :: pn(LBi:,LBj:)
408 real(r8), intent(in) :: pmon_r(LBi:,LBj:)
409 real(r8), intent(in) :: pnom_p(LBi:,LBj:)
410# ifdef MASKING
411 real(r8), intent(in) :: umask(LBi:,LBj:)
412 real(r8), intent(in) :: pmask(LBi:,LBj:)
413# endif
414 real(r8), intent(in) :: Kh(LBi:,LBj:)
415 real(r8), intent(inout) :: ad_A(LBij:)
416# else
417 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
418 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
419 real(r8), intent(in) :: pmon_r(LBi:UBi,LBj:UBj)
420 real(r8), intent(in) :: pnom_p(LBi:UBi,LBj:UBj)
421# ifdef MASKING
422 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
423 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
424# endif
425 real(r8), intent(in) :: Kh(LBi:UBi,LBj:UBj)
426 real(r8), intent(inout) :: ad_A(LBij:UBij)
427# endif
428!
429! Local variable declarations.
430!
431 logical, dimension(4) :: Lconvolve
432
433 integer :: Nnew, Nold, Nsav, i, j, step
434
435 real(r8) :: adfac, cff
436
437 real(r8), dimension(LBij:UBij,2) :: ad_Awrk
438
439 real(r8), dimension(JminS:JmaxS) :: ad_FE
440 real(r8), dimension(IminS:ImaxS) :: ad_FX
441 real(r8), dimension(LBij:UBij) :: Hfac
442
443# include "set_bounds.h"
444!
445!-----------------------------------------------------------------------
446! Initialize adjoint private variables.
447!-----------------------------------------------------------------------
448!
449 ad_awrk(lbij:ubij,1:2)=0.0_r8
450
451 ad_fe(jmins:jmaxs)=0.0_r8
452 ad_fx(imins:imaxs)=0.0_r8
453!
454!-----------------------------------------------------------------------
455! Adjoint space convolution of the diffusion equation for a 1D state
456! variable at U-points.
457!-----------------------------------------------------------------------
458!
459 lconvolve(iwest )=domain(ng)%Western_Edge (tile)
460 lconvolve(ieast )=domain(ng)%Eastern_Edge (tile)
461 lconvolve(isouth)=domain(ng)%Southern_Edge(tile)
462 lconvolve(inorth)=domain(ng)%Northern_Edge(tile)
463!
464! Set integration indices and initial conditions.
465!
466 nold=1
467 nnew=2
468!
469! Compute metrics factor.
470!
471 cff=dtsizeh*0.25_r8
472 IF (lconvolve(boundary)) THEN
473 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
474 i=edge(boundary)
475 DO j=jstr,jend
476 hfac(j)=cff*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))
477 END DO
478 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
479 j=edge(boundary)
480 DO i=istru,iend
481 hfac(i)=cff*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))
482 END DO
483 END IF
484 END IF
485!
486!-----------------------------------------------------------------------
487! Adjoint of load convolved solution.
488!-----------------------------------------------------------------------
489!
490# ifdef DISTRIBUTE
491!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
492!^ & LBij, UBij, &
493!^ & Nghost, &
494!^ & EWperiodic(ng), NSperiodic(ng), &
495!^ & tl_A)
496!^
497 CALL ad_mp_exchange2d_bry (ng, tile, model, 1, boundary, &
498 & lbij, ubij, &
499 & nghost, &
500 & ewperiodic(ng), nsperiodic(ng), &
501 & ad_a)
502# endif
503!^ CALL bc_u2d_bry_tile (ng, tile, boundary, &
504!^ & LBij, UBij, &
505!^ & tl_A)
506!^
507 CALL ad_bc_u2d_bry_tile (ng, tile, boundary, &
508 & lbij, ubij, &
509 & ad_a)
510 IF (lconvolve(boundary)) THEN
511 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
512 DO j=jstr,jend
513!^ tl_A(j)=tl_Awrk(j,Nold)
514!^
515 ad_awrk(j,nold)=ad_awrk(j,nold)+ad_a(j)
516 ad_a(j)=0.0_r8
517 END DO
518 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
519 DO i=istru,iend
520!^ tl_A(i)=tl_Awrk(i,Nold)
521!^
522 ad_awrk(i,nold)=ad_awrk(i,nold)+ad_a(i)
523 ad_a(i)=0.0_r8
524 END DO
525 END IF
526 END IF
527!
528!-----------------------------------------------------------------------
529! Integrate adjoint horizontal diffusion terms.
530!-----------------------------------------------------------------------
531!
532 DO step=1,nhsteps
533!
534! Update integration indices.
535!
536 nsav=nnew
537 nnew=nold
538 nold=nsav
539!
540! Apply adjoint boundary conditions. If applicable, exchange boundary
541! data.
542!
543# ifdef DISTRIBUTE
544!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
545!^ & LBij, UBij, &
546!^ & Nghost, &
547!^ & EWperiodic(ng), NSperiodic(ng), &
548!^ & tl_Awrk(:,Nnew))
549!^
550 CALL ad_mp_exchange2d_bry (ng, tile, model, 1, boundary, &
551 & lbij, ubij, &
552 & nghost, &
553 & ewperiodic(ng), nsperiodic(ng), &
554 & ad_awrk(:,nnew))
555# endif
556!^ CALL bc_u2d_bry_tile (ng, tile, boundary, &
557!^ & LBij, UBij, &
558!^ & tl_Awrk(:,Nnew))
559!^
560 CALL ad_bc_u2d_bry_tile (ng, tile, boundary, &
561 & lbij, ubij, &
562 & ad_awrk(:,nnew))
563!
564! Time-step adjoint horizontal diffusion terms.
565!
566 IF (lconvolve(boundary)) THEN
567 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
568 DO j=jstr,jend
569!^ tl_Awrk(j,Nnew)=tl_Awrk(j,Nold)+ &
570!^ & Hfac(j)* &
571!^ & (tl_FE(j+1)-tl_FE(j))
572!^
573 adfac=hfac(j)*ad_awrk(j,nnew)
574 ad_fe(j )=ad_fe(j )-adfac
575 ad_fe(j+1)=ad_fe(j+1)+adfac
576 ad_awrk(j,nold)=ad_awrk(j,nold)+ &
577 & ad_awrk(j,nnew)
578 ad_awrk(j,nnew)=0.0_r8
579 END DO
580 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
581 DO i=istru,iend
582!^ tl_Awrk(i,Nnew)=tl_Awrk(i,Nold)+ &
583!^ & Hfac(i)* &
584!^ & (tl_FX(i)-tl_FX(i-1))
585!^
586 adfac=hfac(i)*ad_awrk(i,nnew)
587 ad_fx(i-1)=ad_fx(i-1)-adfac
588 ad_fx(i )=ad_fx(i )+adfac
589 ad_awrk(i,nold)=ad_awrk(i,nold)+ &
590 & ad_awrk(i,nnew)
591 ad_awrk(i,nnew)=0.0_r8
592 END DO
593 END IF
594 END IF
595!
596! Compute XI- and ETA-components of diffusive flux.
597!
598 IF (lconvolve(boundary)) THEN
599 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
600 i=edge(boundary)
601 DO j=jstr,jend+1
602# ifdef MASKING
603!^ tl_FE(j)=tl_FE(j)*pmask(i,j)
604!^
605 ad_fe(j)=ad_fe(j)*pmask(i,j)
606# endif
607!^ tl_FE(j)=pnom_p(i,j)*0.25_r8*(Kh(i-1,j )+Kh(i,j )+ &
608!^ & Kh(i-1,j-1)+Kh(i,j-1))* &
609!^ & (tl_Awrk(j ,Nold)- &
610!^ & tl_Awrk(j-1,Nold))
611!^
612 adfac=pnom_p(i,j)*0.25_r8*(kh(i-1,j )+kh(i,j )+ &
613 & kh(i-1,j-1)+kh(i,j-1))* &
614 & ad_fe(j)
615 ad_awrk(j-1,nold)=ad_awrk(j-1,nold)-adfac
616 ad_awrk(j ,nold)=ad_awrk(j ,nold)+adfac
617 ad_fe(j)=0.0_r8
618 END DO
619 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
620 j=edge(boundary)
621 DO i=istru-1,iend
622!^ tl_FX(i)=pmon_r(i,j)*Kh(i,j)* &
623!^ & (tl_Awrk(i+1,Nold)- &
624!^ & tl_Awrk(i ,Nold))
625!^
626 adfac=pmon_r(i,j)*kh(i,j)*ad_fx(i)
627 ad_awrk(i ,nold)=ad_awrk(i ,nold)-adfac
628 ad_awrk(i+1,nold)=ad_awrk(i+1,nold)+adfac
629 ad_fx(i)=0.0_r8
630 END DO
631 END IF
632 END IF
633 END DO
634!
635! Set adjoint initial conditions.
636!
637 IF (lconvolve(boundary)) THEN
638 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
639 DO j=jstr-1,jend+1
640!^ tl_Awrk(j,Nold)=tl_A(j)
641!^
642 ad_a(j)=ad_a(j)+ad_awrk(j,nold)
643 ad_awrk(j,nold)=0.0_r8
644 END DO
645 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
646 DO i=istru-1,iend+1
647!^ tl_Awrk(i,Nold)=tl_A(i)
648!^
649 ad_a(i)=ad_a(i)+ad_awrk(i,nold)
650 ad_awrk(i,nold)=0.0_r8
651 END DO
652 END IF
653 END IF
654# ifdef DISTRIBUTE
655!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
656!^ & LBij, UBij, &
657!^ & Nghost, &
658!^ & EWperiodic(ng), NSperiodic(ng), &
659!^ & tl_A)
660!^
661 CALL ad_mp_exchange2d_bry (ng, tile, model, 1, boundary, &
662 & lbij, ubij, &
663 & nghost, &
664 & ewperiodic(ng), nsperiodic(ng), &
665 & ad_a)
666# endif
667!^ CALL bc_u2d_bry_tile (ng, tile, boundary, &
668!^ & LBij, UBij, &
669!^ & tl_A)
670!^
671 CALL ad_bc_u2d_bry_tile (ng, tile, boundary, &
672 & lbij, ubij, &
673 & ad_a)
674
675 RETURN
676 END SUBROUTINE ad_conv_u2d_bry_tile
677
678!
679!***********************************************************************
680 SUBROUTINE ad_conv_v2d_bry_tile (ng, tile, model, boundary, &
681 & edge, LBij, UBij, &
682 & LBi, UBi, LBj, UBj, &
683 & IminS, ImaxS, JminS, JmaxS, &
684 & Nghost, NHsteps, DTsizeH, &
685 & Kh, &
686 & pm, pn, pmon_p, pnom_r, &
687# ifdef MASKING
688 & vmask, pmask, &
689# endif
690 & ad_A)
691!***********************************************************************
692!
693 USE mod_param
694 USE mod_scalars
695!
697# ifdef DISTRIBUTE
699# endif
700!
701! Imported variable declarations.
702!
703 integer, intent(in) :: ng, tile, model, boundary
704 integer, intent(in) :: edge(4)
705 integer, intent(in) :: LBij, UBij
706 integer, intent(in) :: LBi, UBi, LBj, UBj
707 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
708 integer, intent(in) :: Nghost, NHsteps
709
710 real(r8), intent(in) :: DTsizeH
711!
712# ifdef ASSUMED_SHAPE
713 real(r8), intent(in) :: pm(LBi:,LBj:)
714 real(r8), intent(in) :: pn(LBi:,LBj:)
715 real(r8), intent(in) :: pmon_p(LBi:,LBj:)
716 real(r8), intent(in) :: pnom_r(LBi:,LBj:)
717# ifdef MASKING
718 real(r8), intent(in) :: vmask(LBi:,LBj:)
719 real(r8), intent(in) :: pmask(LBi:,LBj:)
720# endif
721 real(r8), intent(in) :: Kh(LBi:,LBj:)
722 real(r8), intent(inout) :: ad_A(LBij:)
723# else
724 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
725 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
726 real(r8), intent(in) :: pmon_p(LBi:UBi,LBj:UBj)
727 real(r8), intent(in) :: pnom_r(LBi:UBi,LBj:UBj)
728# ifdef MASKING
729 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
730 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
731# endif
732 real(r8), intent(in) :: Kh(LBi:UBi,LBj:UBj)
733 real(r8), intent(inout) :: ad_A(LBij:UBij)
734# endif
735!
736! Local variable declarations.
737!
738 logical, dimension(4) :: Lconvolve
739
740 integer :: Nnew, Nold, Nsav, i, j, step
741
742 real(r8) :: adfac, cff
743
744 real(r8), dimension(LBij:UBij,2) :: ad_Awrk
745
746 real(r8), dimension(JminS:JmaxS) :: ad_FE
747 real(r8), dimension(IminS:ImaxS) :: ad_FX
748 real(r8), dimension(LBij:UBij) :: Hfac
749
750# include "set_bounds.h"
751!
752!-----------------------------------------------------------------------
753! Initialize adjoint private variables.
754!-----------------------------------------------------------------------
755!
756 ad_awrk(lbij:ubij,1:2)=0.0_r8
757
758 ad_fe(jmins:jmaxs)=0.0_r8
759 ad_fx(imins:imaxs)=0.0_r8
760
761!-----------------------------------------------------------------------
762! Adjoint space convolution of the diffusion equation for a 2D state
763! variable at RHO-points.
764!-----------------------------------------------------------------------
765!
766 lconvolve(iwest )=domain(ng)%Western_Edge (tile)
767 lconvolve(ieast )=domain(ng)%Eastern_Edge (tile)
768 lconvolve(isouth)=domain(ng)%Southern_Edge(tile)
769 lconvolve(inorth)=domain(ng)%Northern_Edge(tile)
770!
771! Set integration indices and initial conditions.
772!
773 nold=1
774 nnew=2
775!
776! Compute metrics factor.
777!
778 cff=dtsizeh*0.25_r8
779 IF (lconvolve(boundary)) THEN
780 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
781 i=edge(boundary)
782 DO j=jstrv,jend
783 hfac(j)=cff*(pm(i,j-1)+pm(i,j))*(pn(i,j-1)+pn(i,j))
784 END DO
785 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
786 j=edge(boundary)
787 DO i=istr,iend
788 hfac(i)=cff*(pm(i,j-1)+pm(i,j))*(pn(i,j-1)+pn(i,j))
789 END DO
790 END IF
791 END IF
792!
793!-----------------------------------------------------------------------
794! Adjoint of load convolved solution.
795!-----------------------------------------------------------------------
796!
797# ifdef DISTRIBUTE
798!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
799!^ & LBij, UBij, &
800!^ & Nghost, &
801!^ & EWperiodic(ng), NSperiodic(ng), &
802!^ & tl_A)
803!^
804 CALL ad_mp_exchange2d_bry (ng, tile, model, 1, boundary, &
805 & lbij, ubij, &
806 & nghost, &
807 & ewperiodic(ng), nsperiodic(ng), &
808 & ad_a)
809# endif
810!^ CALL bc_v2d_bry_tile (ng, tile, boundary, &
811!^ & LBij, UBij, &
812!^ & tl_A)
813!^
814 CALL ad_bc_v2d_bry_tile (ng, tile, boundary, &
815 & lbij, ubij, &
816 & ad_a)
817 IF (lconvolve(boundary)) THEN
818 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
819 DO j=jstrv,jend
820!^ tl_A(j)=tl_Awrk(j,Nold)
821!^
822 ad_awrk(j,nold)=ad_awrk(j,nold)+ad_a(j)
823 ad_a(j)=0.0_r8
824 END DO
825 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
826 DO i=istr,iend
827!^ tl_A(i)=tl_Awrk(i,Nold)
828!^
829 ad_awrk(i,nold)=ad_awrk(i,nold)+ad_a(i)
830 ad_a(i)=0.0_r8
831 END DO
832 END IF
833 END IF
834!
835!-----------------------------------------------------------------------
836! Integrate adjoint horizontal diffusion terms.
837!-----------------------------------------------------------------------
838!
839 DO step=1,nhsteps
840!
841! Update integration indices.
842!
843 nsav=nnew
844 nnew=nold
845 nold=nsav
846!
847! Apply boundary conditions. If applicable, exchange boundary data.
848!
849# ifdef DISTRIBUTE
850!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
851!^ & LBij, UBij, &
852!^ & Nghost, &
853!^ & EWperiodic(ng), NSperiodic(ng), &
854!^ & tl_Awrk(:,Nnew))
855!^
856 CALL ad_mp_exchange2d_bry (ng, tile, model, 1, boundary, &
857 & lbij, ubij, &
858 & nghost, &
859 & ewperiodic(ng), nsperiodic(ng), &
860 & ad_awrk(:,nnew))
861# endif
862!^ CALL bc_v2d_bry_tile (ng, tile, boundary, &
863!^ & LBij, UBij, &
864!^ & tl_Awrk(:,Nnew))
865!^
866 CALL ad_bc_v2d_bry_tile (ng, tile, boundary, &
867 & lbij, ubij, &
868 & ad_awrk(:,nnew))
869!
870! Time-step adjoint horizontal diffusion terms.
871!
872 IF (lconvolve(boundary)) THEN
873 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
874 DO j=jstrv,jend
875!^ tl_Awrk(j,Nnew)=tl_Awrk(j,Nold)+ &
876!^ & Hfac(j)* &
877!^ & (tl_FE(j)-tl_FE(j-1))
878!^
879 adfac=hfac(j)*ad_awrk(j,nnew)
880 ad_fe(j-1)=ad_fe(j-1)-adfac
881 ad_fe(j )=ad_fe(j )+adfac
882 ad_awrk(j,nold)=ad_awrk(j,nold)+ &
883 & ad_awrk(j,nnew)
884 ad_awrk(j,nnew)=0.0_r8
885 END DO
886 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
887 DO i=istr,iend
888!^ tl_Awrk(i,Nnew)=tl_Awrk(i,Nold)+ &
889!^ & Hfac(i)* &
890!^ & (tl_FX(i+1)-tl_FX(i))
891!^
892 adfac=hfac(i)*ad_awrk(i,nnew)
893 ad_fx(i )=ad_fx(i )-adfac
894 ad_fx(i+1)=ad_fx(i+1)+adfac
895 ad_awrk(i,nold)=ad_awrk(i,nold)+ &
896 & ad_awrk(i,nnew)
897 ad_awrk(i,nnew)=0.0_r8
898 END DO
899 END IF
900 END IF
901!
902! Compute XI- and ETA-components of adjoint diffusive flux.
903!
904 IF (lconvolve(boundary)) THEN
905 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
906 i=edge(boundary)
907 DO j=jstrv-1,jend
908!^ tl_FE(j)=pnom_r(i,j)*Kh(i,j)* &
909!^ & (tl_Awrk(j+1,Nold)- &
910!^ & tl_Awrk(j ,Nold))
911!^
912 adfac=pnom_r(i,j)*kh(i,j)*ad_fe(j)
913 ad_awrk(j ,nold)=ad_awrk(j ,nold)-adfac
914 ad_awrk(j+1,nold)=ad_awrk(j+1,nold)+adfac
915 ad_fe(j)=0.0_r8
916 END DO
917 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
918 j=edge(boundary)
919 DO i=istr,iend+1
920# ifdef MASKING
921!^ tl_FX(i)=tl_FX(i)*pmask(i,j)
922!^
923 ad_fx(i)=ad_fx(i)*pmask(i,j)
924# endif
925!^ tl_FX(i)=pmon_p(i,j)*0.25_r8*(Kh(i-1,j )+Kh(i,j )+ &
926!^ & Kh(i-1,j-1)+Kh(i,j-1))* &
927!^ & (tl_Awrk(i ,Nold)- &
928!^ & tl_Awrk(i-1,Nold))
929!^
930 adfac=pmon_p(i,j)*0.25_r8*(kh(i-1,j )+kh(i,j )+ &
931 & kh(i-1,j-1)+kh(i,j-1))* &
932 & ad_fx(i)
933 ad_awrk(i-1,nold)=ad_awrk(i-1,nold)-adfac
934 ad_awrk(i ,nold)=ad_awrk(i ,nold)+adfac
935 ad_fx(i)=0.0_r8
936 END DO
937 END IF
938 END IF
939 END DO
940!
941! Set adjoint initial conditions.
942!
943 IF (lconvolve(boundary)) THEN
944 IF ((boundary.eq.iwest).or.(boundary.eq.ieast)) THEN
945 DO j=jstrv-1,jend+1
946!^ tl_Awrk(j,Nold)=tl_A(j)
947!^
948 ad_a(j)=ad_a(j)+ad_awrk(j,nold)
949 ad_awrk(j,nold)=0.0_r8
950 END DO
951 ELSE IF ((boundary.eq.isouth).or.(boundary.eq.inorth)) THEN
952 DO i=istr-1,iend+1
953!^ tl_Awrk(i,Nold)=tl_A(i)
954!^
955 ad_a(i)=ad_a(i)+ad_awrk(i,nold)
956 ad_awrk(i,nold)=0.0_r8
957 END DO
958 END IF
959 END IF
960# ifdef DISTRIBUTE
961!^ CALL mp_exchange2d_bry (ng, tile, model, 1, boundary, &
962!^ & LBij, UBij, &
963!^ & Nghost, &
964!^ & EWperiodic(ng), NSperiodic(ng), &
965!^ & tl_A)
966!^
967 CALL ad_mp_exchange2d_bry (ng, tile, model, 1, boundary, &
968 & lbij, ubij, &
969 & nghost, &
970 & ewperiodic(ng), nsperiodic(ng), &
971 & ad_a)
972# endif
973!^ CALL bc_v2d_bry_tile (ng, tile, boundary, &
974!^ & LBij, UBij, &
975!^ & tl_A)
976!^
977 CALL ad_bc_v2d_bry_tile (ng, tile, boundary, &
978 & lbij, ubij, &
979 & ad_a)
980
981 RETURN
982 END SUBROUTINE ad_conv_v2d_bry_tile
983#endif
984 END MODULE ad_conv_bry2d_mod
subroutine ad_bc_u2d_bry_tile(ng, tile, boundary, lbij, ubij, ad_a)
subroutine ad_bc_v2d_bry_tile(ng, tile, boundary, lbij, ubij, ad_a)
subroutine ad_bc_r2d_bry_tile(ng, tile, boundary, lbij, ubij, ad_a)
Definition ad_bc_bry2d.F:31
subroutine ad_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, ad_a)
subroutine ad_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, ad_a)
subroutine ad_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, ad_a)
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 ad_mp_exchange2d_bry(ng, tile, model, nvar, boundary, lbij, ubij, nghost, ew_periodic, ns_periodic, ad_a, ad_b, ad_c, ad_d)