ROMS
Loading...
Searching...
No Matches
ad_uv3drelax.F
Go to the documentation of this file.
1#include "cppdefs.h"
3
4#if defined TL_IOMS && defined RPM_RELAXATION && defined SOLVE3D
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! This routine is the adjoint of relaxing current representer tangent !
14! linear 3D momentum to previous Picard iteration solution (basic !
15! state) to improve stability and convergence. !
16! !
17!=======================================================================
18!
19 implicit none
20!
21 PRIVATE
22!
23 PUBLIC ad_uv3drelax
24!
25 CONTAINS
26!
27!***********************************************************************
28 SUBROUTINE ad_uv3drelax (ng, tile)
29!***********************************************************************
30!
31 USE mod_param
32 USE mod_coupling
33 USE mod_grid
34 USE mod_ocean
35 USE mod_stepping
36!
37! Imported variable declarations.
38!
39 integer, intent(in) :: ng, tile
40!
41! Local variable declarations.
42!
43 character (len=*), parameter :: myfile = &
44 & __FILE__
45!
46# include "tile.h"
47!
48# ifdef PROFILE
49 CALL wclock_on (ng, iadm, 30, __line__, myfile)
50# endif
51 CALL ad_uv3drelax_tile (ng, tile, &
52 & lbi, ubi, lbj, ubj, &
53 & imins, imaxs, jmins, jmaxs, &
54 & nrhs(ng), nnew(ng), &
55# ifdef MASKING
56 & grid(ng) % pmask, &
57# endif
58 & grid(ng) % Hz, &
59 & grid(ng) % pm, &
60 & grid(ng) % pmon_p, &
61 & grid(ng) % pmon_r, &
62 & grid(ng) % pn, &
63 & grid(ng) % pnom_p, &
64 & grid(ng) % pnom_r, &
65 & ocean(ng) % u, &
66 & ocean(ng) % v, &
67 & ocean(ng) % ad_u, &
68 & ocean(ng) % ad_v)
69# ifdef PROFILE
70 CALL wclock_off (ng, iadm, 30, __line__, myfile)
71# endif
72!
73 RETURN
74 END SUBROUTINE ad_uv3drelax
75
76!
77!***********************************************************************
78 SUBROUTINE ad_uv3drelax_tile (ng, tile, &
79 & LBi, UBi, LBj, UBj, &
80 & IminS, ImaxS, JminS, JmaxS, &
81 & nrhs, nnew, &
82# ifdef MASKING
83 & pmask, &
84# endif
85 & Hz, &
86 & pm, pmon_p, pmon_r, &
87 & pn, pnom_p, pnom_r, &
88 & u, v, &
89 & ad_u, ad_v)
90!***********************************************************************
91!
92 USE mod_param
93 USE mod_scalars
94!
95! Imported variable declarations.
96!
97 integer, intent(in) :: ng, tile
98 integer, intent(in) :: LBi, UBi, LBj, UBj
99 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
100 integer, intent(in) :: nrhs, nnew
101
102# ifdef ASSUMED_SHAPE
103# ifdef MASKING
104 real(r8), intent(in) :: pmask(LBi:,LBj:)
105# endif
106 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
107 real(r8), intent(in) :: pm(LBi:,LBj:)
108 real(r8), intent(in) :: pmon_p(LBi:,LBj:)
109 real(r8), intent(in) :: pmon_r(LBi:,LBj:)
110 real(r8), intent(in) :: pn(LBi:,LBj:)
111 real(r8), intent(in) :: pnom_p(LBi:,LBj:)
112 real(r8), intent(in) :: pnom_r(LBi:,LBj:)
113
114 real(r8), intent(in) :: u(LBi:,LBj:,:,:)
115 real(r8), intent(in) :: v(LBi:,LBj:,:,:)
116
117 real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:)
118 real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:)
119# else
120# ifdef MASKING
121 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
122# endif
123 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
124 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
125 real(r8), intent(in) :: pmon_p(LBi:UBi,LBj:UBj)
126 real(r8), intent(in) :: pmon_r(LBi:UBi,LBj:UBj)
127 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
128 real(r8), intent(in) :: pnom_p(LBi:UBi,LBj:UBj)
129 real(r8), intent(in) :: pnom_r(LBi:UBi,LBj:UBj)
130
131 real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
132 real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
133
134 real(r8), intent(inout) :: ad_u(LBi:UBi,LBj:UBj,N(ng),2)
135 real(r8), intent(inout) :: ad_v(LBi:UBi,LBj:UBj,N(ng),2)
136# endif
137!
138! Local variable declarations.
139!
140 integer :: i, j, k
141
142 real(r8) :: cff, adfac
143
144 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_UFe
145 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_UFx
146 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_VFe
147 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_VFx
148
149# include "set_bounds.h"
150!
151!-----------------------------------------------------------------------
152! Compute horizontal diffusion relaxation of 3D momentum between
153! current and previous representer tangent linear Picard iteration
154! trajectory (basic state).
155!-----------------------------------------------------------------------
156!
157 IF (tl_m3diff(ng).gt.0.0_r8) THEN
158!
159 ad_ufx=0.0_r8
160 ad_vfx=0.0_r8
161 ad_ufe=0.0_r8
162 ad_vfe=0.0_r8
163!
164 k_loop : DO k=1,n(ng)
165!
166! Time-step adjoint diffusive relaxation term. Notice that momentum
167! at this stage is HzU and HzV and has m2/s units.
168!
169 cff=dt(ng)*0.25_r8
170 DO j=jstr,jend
171 DO i=istru,iend
172!^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)+ &
173!^ & cff*(pm(i-1,j)+pm(i,j))* &
174!^ & (pn(i-1,j)+pn(i,j))* &
175!^ & (UFx(i,j)-UFx(i-1,j)+ &
176!^ & UFe(i,j+1)-UFe(i,j))
177!^
178 adfac=cff*(pm(i-1,j)+pm(i,j))* &
179 & (pn(i-1,j)+pn(i,j))*ad_u(i,j,k,nnew)
180 ad_ufx(i,j)=ad_ufx(i,j)+adfac
181 ad_ufx(i-1,j)=ad_ufx(i-1,j)-adfac
182 ad_ufe(i,j+1)=ad_ufe(i,j+1)+adfac
183 ad_ufe(i,j)=ad_ufe(i,j)-adfac
184 END DO
185 END DO
186 DO j=jstrv,jend
187 DO i=istr,iend
188!^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)+ &
189!^ & cff*(pm(i,j)+pm(i,j-1))* &
190!^ & (pn(i,j)+pn(i,j-1))* &
191!^ & (VFx(i+1,j)-VFx(i,j)+ &
192!^ & VFe(i,j)-VFe(i,j-1))
193!^
194 adfac=cff*(pm(i,j)+pm(i,j-1))* &
195 & (pn(i,j)+pn(i,j-1))*ad_v(i,j,k,nnew)
196 ad_vfx(i+1,j)=ad_vfx(i+1,j)+adfac
197 ad_vfx(i,j)=ad_vfx(i,j)-adfac
198 ad_vfe(i,j)=ad_vfe(i,j)+adfac
199 ad_vfe(i,j-1)=ad_vfe(i,j-1)-adfac
200 END DO
201 END DO
202!
203! Compute flux-components of the adjoint diffusive relaxation (m4/s2)
204! in XI- and ETA-directions.
205!
206 DO j=jstrv-1,jend
207 DO i=istr,iend
208!^ VFe(i,j)=tl_M3diff(ng)*pnom_r(i,j)* &
209!^ & Hz(i,j,k)* &
210!^ & (tl_v(i,j+1,k,nrhs)-v(i,j+1,k,nrhs)- &
211!^ & tl_v(i,j ,k,nrhs)+v(i,j ,k,nrhs))
212!^
213 adfac=tl_m3diff(ng)*pnom_r(i,j)*hz(i,j,k)*ad_vfe(i,j)
214 ad_v(i,j+1,k,nrhs)=ad_v(i,j+1,k,nrhs)+adfac
215 ad_v(i,j ,k,nrhs)=ad_v(i,j ,k,nrhs)-adfac
216 ad_vfe(i,j)=0.0_r8
217 END DO
218 END DO
219 DO j=jstrv,jend
220 DO i=istr,iend+1
221# ifdef MASKING
222!^ VFx(i,j)=VFx(i,j)*pmask(i,j)
223!^
224 ad_vfx(i,j)=ad_vfx(i,j)*pmask(i,j)
225# endif
226!^ VFx(i,j)=tl_M3diff(ng)*pmon_p(i,j)* &
227!^ & 0.25_r8*(Hz(i,j ,k)+Hz(i-1,j ,k)+ &
228!^ & Hz(i,j-1,k)+Hz(i-1,j-1,k))* &
229!^ & (tl_v(i ,j,k,nrhs)-v(i ,j,k,nrhs)- &
230!^ & tl_v(i-1,j,k,nrhs)+v(i-1,j,k,nrhs))
231!^
232 adfac=tl_m3diff(ng)*pmon_p(i,j)* &
233 & 0.25_r8*(hz(i,j ,k)+hz(i-1,j ,k)+ &
234 & hz(i,j-1,k)+hz(i-1,j-1,k))*ad_vfx(i,j)
235 ad_v(i ,j,k,nrhs)=ad_v(i ,j,k,nrhs)+adfac
236 ad_v(i-1,j,k,nrhs)=ad_v(i-1,j,k,nrhs)-adfac
237 ad_vfx(i,j)=0.0_r8
238 END DO
239 END DO
240 DO j=jstr,jend+1
241 DO i=istru,iend
242# ifdef MASKING
243!^ UFe(i,j)=UFe(i,j)*pmask(i,j)
244!^
245 ad_ufe(i,j)=ad_ufe(i,j)*pmask(i,j)
246# endif
247!^ UFe(i,j)=tl_M3diff(ng)*pnom_p(i,j)* &
248!^ & 0.25_r8*(Hz(i,j ,k)+Hz(i-1,j ,k)+ &
249!^ & Hz(i,j-1,k)+Hz(i-1,j-1,k))* &
250!^ & (tl_u(i,j ,k,nrhs)-u(i,j ,k,nrhs)- &
251!^ & tl_u(i,j-1,k,nrhs)+u(i,j-1,k,nrhs))
252!^
253 adfac=tl_m3diff(ng)*pnom_p(i,j)* &
254 & 0.25_r8*(hz(i,j ,k)+hz(i-1,j ,k)+ &
255 & hz(i,j-1,k)+hz(i-1,j-1,k))*ad_ufe(i,j)
256 ad_u(i,j ,k,nrhs)=ad_u(i,j ,k,nrhs)+adfac
257 ad_u(i,j-1,k,nrhs)=ad_u(i,j-1,k,nrhs)-adfac
258 ad_ufe(i,j)=0.0_r8
259 END DO
260 END DO
261 DO j=jstr,jend
262 DO i=istru-1,iend
263!^ UFx(i,j)=tl_M3diff(ng)*pmon_r(i,j)* &
264!^ & Hz(i,j,k)* &
265!^ & (tl_u(i+1,j,k,nrhs)-u(i+1,j,k,nrhs)- &
266!^ & tl_u(i ,j,k,nrhs)+u(i ,j,k,nrhs))
267!^
268 adfac=tl_m3diff(ng)*pmon_r(i,j)*hz(i,j,k)*ad_ufx(i,j)
269 ad_u(i+1,j,k,nrhs)=ad_u(i+1,j,k,nrhs)+adfac
270 ad_u(i ,j,k,nrhs)=ad_u(i ,j,k,nrhs)-adfac
271 ad_ufx(i,j)=0.0_r8
272 END DO
273 END DO
274 END DO k_loop
275 END IF
276!
277 RETURN
278 END SUBROUTINE ad_uv3drelax_tile
279#endif
280 END MODULE ad_uv3drelax_mod
subroutine ad_uv3drelax_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nnew, pmask, hz, pm, pmon_p, pmon_r, pn, pnom_p, pnom_r, u, v, ad_u, ad_v)
subroutine, public ad_uv3drelax(ng, tile)
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, parameter iadm
Definition mod_param.F:665
real(dp), dimension(:), allocatable dt
real(r8), dimension(:), allocatable tl_m3diff
integer, dimension(:), allocatable nrhs
integer, dimension(:), allocatable nnew
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3