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

Functions/Subroutines

subroutine, public ad_t3drelax (ng, tile)
 
subroutine ad_t3drelax_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nnew, umask, vmask, hz, pmon_u, pnom_v, pm, pn, t, ad_t)
 

Function/Subroutine Documentation

◆ ad_t3drelax()

subroutine, public ad_t3drelax_mod::ad_t3drelax ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 28 of file ad_t3drelax.F.

29!***********************************************************************
30!
31 USE mod_param
32 USE mod_grid
33 USE mod_ocean
34 USE mod_stepping
35!
36! Imported variable declarations.
37!
38 integer, intent(in) :: ng, tile
39!
40! Local variable declarations.
41!
42 character (len=*), parameter :: MyFile = &
43 & __FILE__
44!
45# include "tile.h"
46!
47# ifdef PROFILE
48 CALL wclock_on (ng, iadm, 24, __line__, myfile)
49# endif
50 CALL ad_t3drelax_tile (ng, tile, &
51 & lbi, ubi, lbj, ubj, &
52 & imins, imaxs, jmins, jmaxs, &
53 & nrhs(ng), nnew(ng), &
54# ifdef MASKING
55 & grid(ng) % umask, &
56 & grid(ng) % vmask, &
57# endif
58 & grid(ng) % Hz, &
59 & grid(ng) % pmon_u, &
60 & grid(ng) % pnom_v, &
61 & grid(ng) % pm, &
62 & grid(ng) % pn, &
63 & ocean(ng) % t, &
64 & ocean(ng) % ad_t)
65# ifdef PROFILE
66 CALL wclock_off (ng, iadm, 24, __line__, myfile)
67# endif
68!
69 RETURN
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
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

References ad_t3drelax_tile(), mod_grid::grid, mod_param::iadm, mod_stepping::nnew, mod_stepping::nrhs, mod_ocean::ocean, wclock_off(), and wclock_on().

Referenced by ad_rhs3d_mod::ad_rhs3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ad_t3drelax_tile()

subroutine ad_t3drelax_mod::ad_t3drelax_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) nrhs,
integer, intent(in) nnew,
real(r8), dimension(lbi:,lbj:), intent(in) umask,
real(r8), dimension(lbi:,lbj:), intent(in) vmask,
real(r8), dimension(lbi:,lbj:,:), intent(in) hz,
real(r8), dimension(lbi:,lbj:), intent(in) pmon_u,
real(r8), dimension(lbi:,lbj:), intent(in) pnom_v,
real(r8), dimension(lbi:,lbj:), intent(in) pm,
real(r8), dimension(lbi:,lbj:), intent(in) pn,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(in) t,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) ad_t )
private

Definition at line 73 of file ad_t3drelax.F.

83!***********************************************************************
84!
85 USE mod_param
86 USE mod_scalars
87!
88! Imported variable declarations.
89!
90 integer, intent(in) :: ng, tile
91 integer, intent(in) :: LBi, UBi, LBj, UBj
92 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
93 integer, intent(in) :: nrhs, nnew
94
95# ifdef ASSUMED_SHAPE
96# ifdef MASKING
97 real(r8), intent(in) :: umask(LBi:,LBj:)
98 real(r8), intent(in) :: vmask(LBi:,LBj:)
99# endif
100 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
101 real(r8), intent(in) :: pmon_u(LBi:,LBj:)
102 real(r8), intent(in) :: pnom_v(LBi:,LBj:)
103 real(r8), intent(in) :: pm(LBi:,LBj:)
104 real(r8), intent(in) :: pn(LBi:,LBj:)
105 real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
106
107 real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:)
108# else
109# ifdef MASKING
110 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
111 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
112# endif
113 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
114 real(r8), intent(in) :: pmon_u(LBi:UBi,LBj:UBj)
115 real(r8), intent(in) :: pnom_v(LBi:UBi,LBj:UBj)
116 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
117 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
118 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
119
120 real(r8), intent(inout) :: ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
121# endif
122!
123! Local variable declarations.
124!
125 integer :: i, itrc, j, k
126 real(r8) :: adfac
127
128 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_FE
129 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_FX
130
131# include "set_bounds.h"
132!
133!----------------------------------------------------------------------
134! Compute horizontal diffusion relaxation of tracers between current
135! and previous representer tangent linear Picard iteration trajectory
136! (basic state).
137!----------------------------------------------------------------------
138!
139 ad_fx=0.0_r8
140 ad_fe=0.0_r8
141!
142 DO itrc=1,nt(ng)
143 IF (tl_tdiff(itrc,ng).gt.0.0_r8) THEN
144 DO k=1,n(ng)
145!
146! Time-step horizontal adjoint diffusion relaxation term (m Tunits).
147!
148 DO j=jstr,jend
149 DO i=istr,iend
150!^ tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)+ &
151!^ & dt(ng)*pm(i,j)*pn(i,j)* &
152!^ & (FX(i+1,j)-FX(i,j)+ &
153!^ & FE(i,j+1)-FE(i,j))
154!^
155 adfac=dt(ng)*pm(i,j)*pn(i,j)*ad_t(i,j,k,nnew,itrc)
156 ad_fx(i+1,j)=ad_fx(i+1,j)+adfac
157 ad_fx(i,j)=ad_fx(i,j)-adfac
158 ad_fe(i,j+1)=ad_fe(i,j+1)+adfac
159 ad_fe(i,j)=ad_fe(i,j)+adfac
160 END DO
161 END DO
162!
163! Compute XI- and ETA-components of adjoint diffusive tracer flux
164! (T m3/s).
165!
166 DO j=jstr,jend+1
167 DO i=istr,iend
168# ifdef MASKING
169!^ FE(i,j)=FE(i,j)*vmask(i,j)
170!^
171 ad_fe(i,j)=ad_fe(i,j)*vmask(i,j)
172# endif
173!^ FE(i,j)=0.5_r8*tl_Tdiff(itrc,ng)*pnom_v(i,j)* &
174!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
175!^ & (tl_t(i,j ,k,nrhs,itrc)-t(i,j ,k,nrhs,itrc)- &
176!^ & tl_t(i,j-1,k,nrhs,itrc)+t(i,j-1,k,nrhs,itrc))
177!^
178 adfac=0.5_r8*tl_tdiff(itrc,ng)*pnom_v(i,j)* &
179 & (hz(i,j,k)+hz(i,j-1,k))*ad_fe(i,j)
180 ad_t(i,j ,k,nrhs,itrc)=ad_t(i,j ,k,nrhs,itrc)+adfac
181 ad_t(i,j-1,k,nrhs,itrc)=ad_t(i,j-1,k,nrhs,itrc)-adfac
182 ad_fe(i,j)=0.0_r8
183 END DO
184 END DO
185!
186 DO j=jstr,jend
187 DO i=istr,iend+1
188# ifdef MASKING
189!^ FX(i,j)=FX(i,j)*umask(i,j)
190!^
191 ad_fx(i,j)=ad_fx(i,j)*umask(i,j)
192# endif
193!^ FX(i,j)=0.5_r8*tl_Tdiff(itrc,ng)*pmon_u(i,j)* &
194!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
195!^ & (tl_t(i ,j,k,nrhs,itrc)-t(i ,j,k,nrhs,itrc)- &
196!^ & tl_t(i-1,j,k,nrhs,itrc)+t(i-1,j,k,nrhs,itrc))
197!^
198 adfac=0.5_r8*tl_tdiff(itrc,ng)*pmon_u(i,j)* &
199 & (hz(i,j,k)+hz(i-1,j,k))*ad_fx(i,j)
200 ad_t(i ,j,k,nrhs,itrc)=ad_t(i ,j,k,nrhs,itrc)+adfac
201 ad_t(i-1,j,k,nrhs,itrc)=ad_t(i-1,j,k,nrhs,itrc)-adfac
202 ad_fx(i,j)=0.0_r8
203 END DO
204 END DO
205 END DO
206 END IF
207 END DO
208!
209 RETURN
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer, dimension(:), allocatable nt
Definition mod_param.F:489
real(r8), dimension(:,:), allocatable tl_tdiff
real(dp), dimension(:), allocatable dt

References mod_scalars::dt, and mod_scalars::tl_tdiff.

Referenced by ad_t3drelax().

Here is the caller graph for this function: