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

Functions/Subroutines

subroutine, public adsen_initial (ng, tile)
 
subroutine adsen_initial_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, knew, nstp, rscope, uscope, vscope, u_adsg, v_adsg, wvel_adsg, t_adsg, ubar_adsg, vbar_adsg, zeta_adsg, ad_u, ad_v, ad_wvel, ad_t, ad_ubar, ad_vbar, ad_zeta)
 

Function/Subroutine Documentation

◆ adsen_initial()

subroutine, public adsen_initial_mod::adsen_initial ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 29 of file adsen_initial.F.

30!***********************************************************************
31!
32 USE mod_param
33 USE mod_clima
34 USE mod_grid
35 USE mod_ocean
36 USE mod_stepping
37!
38! Imported variable declarations.
39!
40 integer, intent(in) :: ng, tile
41!
42! Local variable declarations.
43!
44# include "tile.h"
45!
46 CALL adsen_initial_tile (ng, tile, &
47 & lbi, ubi, lbj, ubj, &
48 & imins, imaxs, jmins, jmaxs, &
49 & knew(ng), &
50# ifdef SOLVE3D
51 & nstp(ng), &
52# endif
53 & grid(ng) % Rscope, &
54 & grid(ng) % Uscope, &
55 & grid(ng) % Vscope, &
56# ifdef SOLVE3D
57 & clima(ng) % u_adsG, &
58 & clima(ng) % v_adsG, &
59 & clima(ng) % wvel_adsG, &
60 & clima(ng) % t_adsG, &
61# endif
62 & clima(ng) % ubar_adsG, &
63 & clima(ng) % vbar_adsG, &
64 & clima(ng) % zeta_adsG, &
65# ifdef SOLVE3D
66 & ocean(ng) % ad_u, &
67 & ocean(ng) % ad_v, &
68 & ocean(ng) % ad_wvel, &
69 & ocean(ng) % ad_t, &
70# endif
71 & ocean(ng) % ad_ubar, &
72 & ocean(ng) % ad_vbar, &
73 & ocean(ng) % ad_zeta)
74
75 RETURN
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, dimension(:), allocatable knew
integer, dimension(:), allocatable nstp

References adsen_initial_tile(), mod_clima::clima, mod_grid::grid, mod_stepping::knew, mod_stepping::nstp, and mod_ocean::ocean.

Referenced by ad_initial().

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

◆ adsen_initial_tile()

subroutine adsen_initial_mod::adsen_initial_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) knew,
integer, intent(in) nstp,
real(r8), dimension(lbi:,lbj:), intent(in) rscope,
real(r8), dimension(lbi:,lbj:), intent(in) uscope,
real(r8), dimension(lbi:,lbj:), intent(in) vscope,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) u_adsg,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) v_adsg,
real(r8), dimension(lbi:,lbj:,:,:), intent(in) wvel_adsg,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(in) t_adsg,
real(r8), dimension(lbi:,lbj:,:), intent(in) ubar_adsg,
real(r8), dimension(lbi:,lbj:,:), intent(in) vbar_adsg,
real(r8), dimension(lbi:,lbj:,:), intent(in) zeta_adsg,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) ad_u,
real(r8), dimension(lbi:,lbj:,:,:), intent(inout) ad_v,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_wvel,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) ad_t,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_ubar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_vbar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_zeta )
private

Definition at line 79 of file adsen_initial.F.

97!***********************************************************************
98!
99 USE mod_param
100 USE mod_ncparam
101 USE mod_scalars
102!
103! Imported variable declarations.
104!
105 integer, intent(in) :: ng, tile
106 integer, intent(in) :: LBi, UBi, LBj, UBj
107 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
108 integer, intent(in) :: knew
109# ifdef SOLVE3D
110 integer, intent(in) :: nstp
111# endif
112!
113# ifdef ASSUMED_SHAPE
114 real(r8), intent(in) :: Rscope(LBi:,LBj:)
115 real(r8), intent(in) :: Uscope(LBi:,LBj:)
116 real(r8), intent(in) :: Vscope(LBi:,LBj:)
117# ifdef SOLVE3D
118 real(r8), intent(in) :: u_adsG(LBi:,LBj:,:,:)
119 real(r8), intent(in) :: v_adsG(LBi:,LBj:,:,:)
120 real(r8), intent(in) :: wvel_adsG(LBi:,LBj:,:,:)
121 real(r8), intent(in) :: t_adsG(LBi:,LBj:,:,:,:)
122# endif
123 real(r8), intent(in) :: ubar_adsG(LBi:,LBj:,:)
124 real(r8), intent(in) :: vbar_adsG(LBi:,LBj:,:)
125 real(r8), intent(in) :: zeta_adsG(LBi:,LBj:,:)
126# ifdef SOLVE3D
127 real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:)
128 real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:)
129 real(r8), intent(inout) :: ad_wvel(LBi:,LBj:,:)
130 real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:)
131# endif
132 real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:)
133 real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:)
134 real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:)
135# else
136 real(r8), intent(in) :: Rscope(LBi:UBi,LBj:UBj)
137 real(r8), intent(in) :: Uscope(LBi:UBi,LBj:UBj)
138 real(r8), intent(in) :: Vscope(LBi:UBi,LBj:UBj)
139# ifdef SOLVE3D
140 real(r8), intent(in) :: u_adsG(LBi:UBi,LBj:UBj,N(ng),2)
141 real(r8), intent(in) :: v_adsG(LBi:UBi,LBj:UBj,N(ng),2)
142 real(r8), intent(in) :: wvel_adsG(LBi:UBi,LBj:UBj,N(ng),2)
143 real(r8), intent(in) :: t_adsG(LBi:UBi,LBj:UBj,N(ng),2,NT(ng))
144# endif
145 real(r8), intent(in) :: ubar_adsG(LBi:UBi,LBj:UBj,2)
146 real(r8), intent(in) :: vbar_adsG(LBi:UBi,LBj:UBj,2)
147 real(r8), intent(in) :: zeta_adsG(LBi:UBi,LBj:UBj,2)
148# ifdef SOLVE3D
149 real(r8), intent(inout) :: ad_u(LBi:UBi,LBj:UBj,N(ng),2)
150 real(r8), intent(inout) :: ad_v(LBi:UBi,LBj:UBj,N(ng),2)
151 real(r8), intent(inout) :: ad_wvel(LBi:UBi,LBj:UBj,0:N(ng))
152 real(r8), intent(inout) :: ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
153# endif
154 real(r8), intent(inout) :: ad_ubar(LBi:UBi,LBj:UBj,:)
155 real(r8), intent(inout) :: ad_vbar(LBi:UBi,LBj:UBj,:)
156 real(r8), intent(inout) :: ad_zeta(LBi:UBi,LBj:UBj,:)
157# endif
158!
159! Local variable declarations.
160!
161 integer :: i, itrc, j, k
162
163# include "set_bounds.h"
164!
165!------------------------------------------------------------------------
166! Initialize adjoint staye with the functional whose sensitivity is
167! required. Use functional loaded into first record of climatological
168! arrays.
169!------------------------------------------------------------------------
170!
171! Free-surface.
172!
173 IF (scalars(ng)%Lstate(isfsur)) THEN
174 DO j=jstrr,jendr
175 DO i=istrr,iendr
176 ad_zeta(i,j,knew)=zeta_adsg(i,j,1)*rscope(i,j)
177 END DO
178 END DO
179 END IF
180!
181! 2D Momentum.
182!
183 IF (scalars(ng)%Lstate(isubar)) THEN
184 DO j=jstrr,jendr
185 DO i=istr,iendr
186 ad_ubar(i,j,knew)=ubar_adsg(i,j,1)*uscope(i,j)
187 END DO
188 END DO
189 END IF
190!
191 IF (scalars(ng)%Lstate(isvbar)) THEN
192 DO j=jstr,jendr
193 DO i=istrr,iendr
194 ad_vbar(i,j,knew)=vbar_adsg(i,j,1)*vscope(i,j)
195 END DO
196 END DO
197 END IF
198# ifdef SOLVE3D
199!
200! 3D Momentum.
201!
202 IF (scalars(ng)%Lstate(isuvel)) THEN
203 DO k=kstrs(ng),kends(ng)
204 DO j=jstrr,jendr
205 DO i=istr,iendr
206 ad_u(i,j,k,nstp)=u_adsg(i,j,k,1)*uscope(i,j)
207 END DO
208 END DO
209 END DO
210 END IF
211!
212 IF (scalars(ng)%Lstate(isvvel)) THEN
213 DO k=kstrs(ng),kends(ng)
214 DO j=jstr,jendr
215 DO i=istrr,iendr
216 ad_v(i,j,k,nstp)=v_adsg(i,j,k,1)*vscope(i,j)
217 END DO
218 END DO
219 END DO
220 END IF
221!
222! Vertical velocity.
223!
224 IF (scalars(ng)%Lstate(isvvel)) THEN
225 DO k=kstrs(ng),kends(ng)
226 DO j=jstrr,jendr
227 DO i=istrr,iendr
228 ad_wvel(i,j,k)=wvel_adsg(i,j,k,1)*rscope(i,j)
229 END DO
230 END DO
231 END DO
232 END IF
233!
234! Tracers.
235!
236 DO itrc=1,nt(ng)
237 IF (scalars(ng)%Lstate(istvar(itrc))) THEN
238 DO k=kstrs(ng),kends(ng)
239 DO j=jstrr,jendr
240 DO i=istrr,iendr
241 ad_t(i,j,k,nstp,itrc)=t_adsg(i,j,k,1,itrc)*rscope(i,j)
242 END DO
243 END DO
244 END DO
245 END IF
246 END DO
247# endif
248
249 RETURN
integer isvvel
integer isvbar
integer, dimension(:), allocatable istvar
integer isuvel
integer isfsur
integer isubar
integer, dimension(:), allocatable nt
Definition mod_param.F:489
integer, dimension(:), allocatable kends
integer, dimension(:), allocatable kstrs
type(t_scalars), dimension(:), allocatable scalars
Definition mod_scalars.F:65

References mod_ncparam::isfsur, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::kends, mod_scalars::kstrs, and mod_scalars::scalars.

Referenced by adsen_initial().

Here is the caller graph for this function: