ROMS
Loading...
Searching...
No Matches
propagator_hop.h
Go to the documentation of this file.
1 MODULE propagator_mod
2!
3!git $Id$
4!================================================== Hernan G. Arango ===
5! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
6! Licensed under a MIT/X style license !
7! See License_ROMS.md !
8!=======================================================================
9! !
10! Optimal Perturbation (Hessian singular vectors) Propagator: !
11! !
12! This routine solves the singular vectors of the propagator R(0,t) !
13! which measure the fastest growing of all possible perturbations !
14! over a given time interval. It involves a single integration of a !
15! perturbation forward in time with the tangent linear model over !
16! [0,t], followed by an integration of the result backward in time !
17! with the adjoint model over [t,0]. !
18! !
19! Reference: !
20! !
21! Moore, A.M. et al., 2004: A comprehensive ocean prediction and !
22! analysis system based on the tangent linear and adjoint of a !
23! regional ocean model, Ocean Modelling, 7, 227-258. !
24! !
25! Smith, K. D., A. M. Moore, and H. G. Arango, 2015: Estimates of !
26! ocean forecast error covariance derived from Hessian Singular !
27! Vectors, Ocean Modelling, 89, 104-121, !
28! doi:10.1016/j.ocemod.2015.03.003. !
29! !
30!=======================================================================
31!
32 USE mod_kinds
33!
34 implicit none
35!
36 PRIVATE
37 PUBLIC :: propagator_hop
38!
39 CONTAINS
40!
41!***********************************************************************
42 SUBROUTINE propagator_hop (RunInterval, state, ad_state)
43!***********************************************************************
44!
45 USE mod_param
46 USE mod_parallel
47#ifdef SOLVE3D
48 USE mod_coupling
49#endif
50 USE mod_iounits
51 USE mod_ocean
52 USE mod_scalars
53 USE mod_stepping
54!
55 USE close_io_mod, ONLY : close_inp
59 USE inner2state_mod, ONLY : ini_c_norm
60#ifdef SOLVE3D
61 USE set_depth_mod, ONLY : set_depth
62#endif
63 USE strings_mod, ONLY : founderror
64!
65! Imported variable declarations.
66!
67 real(dp), intent(in) :: runinterval
68!
69 TYPE (t_gst), intent(in) :: state(ngrids)
70 TYPE (t_gst), intent(inout) :: ad_state(ngrids)
71!
72! Local variable declarations.
73!
74 integer :: ng, tile
75 integer :: ktmp, ntmp, lini
76!
77 real(r8) :: statenorm(ngrids)
78!
79 character (len=*), parameter :: myfile = &
80 & __FILE__
81!
82!=======================================================================
83! Forward integration of the tangent linear model.
84!=======================================================================
85!
86 nrun=nrun+1
87 IF (master) THEN
88 DO ng=1,ngrids
89 WRITE (stdout,10) ' PROPAGATOR - Grid: ', ng, &
90 & ', Iteration: ', nrun, &
91 & ', number converged RITZ values: ', &
92 & nconv(ng)
93 END DO
94 END IF
95!
96! Initialize time stepping indices and counters.
97!
98 DO ng=1,ngrids
99 iif(ng)=1
100 indx1(ng)=1
101 kstp(ng)=1
102 krhs(ng)=1
103 knew(ng)=1
104 predictor_2d_step(ng)=.false.
105!
106 iic(ng)=0
107 nstp(ng)=1
108 nrhs(ng)=1
109 nnew(ng)=1
110!
111 synchro_flag(ng)=.true.
112 tdays(ng)=dstart
113 time(ng)=tdays(ng)*day2sec
114 ntstart(ng)=int((time(ng)-dstart*day2sec)/dt(ng))+1
115 ntend(ng)=ntimes(ng)
116 ntfirst(ng)=ntstart(ng)
117 END DO
118!
119 lini=1
120!
121!-----------------------------------------------------------------------
122! Clear tangent linear state variables. There is not need to clean
123! the basic state arrays since they were zeroth out at initialization
124! and bottom of previous iteration.
125!-----------------------------------------------------------------------
126!
127 DO ng=1,ngrids
128 DO tile=first_tile(ng),last_tile(ng),+1
129 CALL initialize_ocean (ng, tile, itlm)
130 END DO
131 END DO
132
133#ifdef SOLVE3D
134!
135!-----------------------------------------------------------------------
136! Compute basic state initial level thicknesses used for state norm
137! scaling. It uses zero time-averaged free-surface (rest state).
138! Therefore, the norm scaling is time invariant.
139!-----------------------------------------------------------------------
140!
141 DO ng=1,ngrids
142 DO tile=last_tile(ng),first_tile(ng),-1
143 CALL set_depth (ng, tile, itlm)
144 END DO
145 END DO
146#endif
147!
148!-----------------------------------------------------------------------
149! Compute tangent linear initial conditions from state vector.
150!-----------------------------------------------------------------------
151!
152 DO ng=1,ngrids
153 DO tile=first_tile(ng),last_tile(ng),+1
154 CALL tl_inner2state (ng, tile, lini, state(ng)%vector)
155 END DO
156 END DO
157!
158!-----------------------------------------------------------------------
159! Compute initial tangent linear state analysis error norm.
160!-----------------------------------------------------------------------
161!
162 DO ng=1,ngrids
163 DO tile=last_tile(ng),first_tile(ng),-1
164 CALL ini_c_norm (ng, tile, kstp(ng), nstp(ng), &
165 & statenorm(ng))
166 END DO
167 IF (master) THEN
168 WRITE (stdout,20) ' PROPAGATOR - Grid: ', ng, &
169 & ', Tangent Initial Norm: ', statenorm(ng)
170 END IF
171 END DO
172!
173!-----------------------------------------------------------------------
174! Read in initial forcing, climatology and assimilation data from
175! input NetCDF files. It loads the first relevant data record for
176! the time-interpolation between snapshots.
177!-----------------------------------------------------------------------
178!
179 DO ng=1,ngrids
180 CALL close_inp (ng, itlm)
181 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
182
183 CALL tl_get_idata (ng)
184 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
185
186 CALL tl_get_data (ng)
187 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
188 END DO
189!
190!-----------------------------------------------------------------------
191! Time-step the tangent linear model.
192!-----------------------------------------------------------------------
193!
194 DO ng=1,ngrids
195 IF (master) THEN
196 WRITE (stdout,30) 'TL', ng, ntstart(ng), ntend(ng)
197 END IF
198 time(ng)=time(ng)-dt(ng)
199 iic(ng)=ntstart(ng)-1
200 END DO
201
202#ifdef SOLVE3D
203 CALL tl_main3d (runinterval)
204#else
205 CALL tl_main2d (runinterval)
206#endif
207 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
208!
209!-----------------------------------------------------------------------
210! Clear nonlinear (basic state) and adjoint state variables.
211!-----------------------------------------------------------------------
212!
213 DO ng=1,ngrids
214 DO tile=first_tile(ng),last_tile(ng),+1
215 CALL initialize_ocean (ng, tile, inlm)
216 CALL initialize_ocean (ng, tile, iadm)
217#ifdef SOLVE3D
218 CALL initialize_coupling (ng, tile, 0)
219#endif
220 END DO
221 END DO
222
223#ifdef SOLVE3D
224!
225!-----------------------------------------------------------------------
226! Compute basic state final level thicknesses used for state norm
227! scaling. It uses zero time-averaged free-surface (rest state).
228! Therefore, the norm scaling is time invariant.
229!-----------------------------------------------------------------------
230!
231 DO ng=1,ngrids
232 DO tile=last_tile(ng),first_tile(ng),-1
233 CALL set_depth (ng, tile, itlm)
234 END DO
235 END DO
236#endif
237!
238!-----------------------------------------------------------------------
239! Compute final tangent linear energy norm.
240!-----------------------------------------------------------------------
241!
242 DO ng=1,ngrids
243 DO tile=first_tile(ng),last_tile(ng),+1
244 CALL tl_statenorm (ng, tile, knew(ng), nstp(ng), &
245 & statenorm(ng))
246 END DO
247 IF (master) THEN
248 WRITE (stdout,20) ' PROPAGATOR - Grid: ', ng, &
249 & ', Tangent Final Norm: ', statenorm(ng)
250 END IF
251 END DO
252!
253!=======================================================================
254! Backward integration with the adjoint model.
255!=======================================================================
256!
257! Initialize time stepping indices and counters.
258!
259 DO ng=1,ngrids
260 iif(ng)=1
261 indx1(ng)=1
262 ktmp=knew(ng)
263 kstp(ng)=1
264 krhs(ng)=3
265 knew(ng)=2
266 predictor_2d_step(ng)=.false.
267!
268 iic(ng)=0
269 ntmp=nstp(ng)
270 nstp(ng)=1
271 nrhs(ng)=1
272 nnew(ng)=2
273!
274 synchro_flag(ng)=.true.
275 tdays(ng)=dstart+dt(ng)*real(ntimes(ng),r8)*sec2day
276 time(ng)=tdays(ng)*day2sec
277 ntstart(ng)=ntimes(ng)+1
278 ntend(ng)=1
279 ntfirst(ng)=ntend(ng)
280 END DO
281!
282!-----------------------------------------------------------------------
283! Initialize adjoint model with the final tangent linear solution
284! scaled by the energy norm.
285!-----------------------------------------------------------------------
286!
287 DO ng=1,ngrids
288 DO tile=last_tile(ng),first_tile(ng),-1
289 CALL ad_ini_perturb (ng, tile, &
290 & ktmp, knew(ng), ntmp, nstp(ng))
291 END DO
292 END DO
293!
294!-----------------------------------------------------------------------
295! Read in initial forcing, climatology and assimilation data from
296! input NetCDF files. It loads the first relevant data record for
297! the time-interpolation between snapshots.
298!-----------------------------------------------------------------------
299!
300 DO ng=1,ngrids
301 CALL close_inp (ng, iadm)
302 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
303
304 CALL ad_get_idata (ng)
305 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
306
307 CALL ad_get_data (ng)
308 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
309 END DO
310!
311!-----------------------------------------------------------------------
312! Time-step the adjoint model backwards.
313!-----------------------------------------------------------------------
314!
315 DO ng=1,ngrids
316 IF (master) THEN
317 WRITE (stdout,30) 'AD', ng, ntstart(ng), ntend(ng)
318 END IF
319 time(ng)=time(ng)+dt(ng)
320 iic(ng)=ntstart(ng)+1
321 END DO
322
323#ifdef SOLVE3D
324 CALL ad_main3d (runinterval)
325#else
326 CALL ad_main2d (runinterval)
327#endif
328 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
329!
330!-----------------------------------------------------------------------
331! Clear nonlinear state (basic state) variables for next iteration
332! and to insure a rest state time averaged free-surface before adjoint
333! state norm scaling.
334!-----------------------------------------------------------------------
335!
336 DO ng=1,ngrids
337 DO tile=first_tile(ng),last_tile(ng),+1
338 CALL initialize_ocean (ng, tile, inlm)
339#ifdef SOLVE3D
340 CALL initialize_coupling (ng, tile, 0)
341#endif
342 END DO
343 END DO
344
345#ifdef SOLVE3D
346!
347!-----------------------------------------------------------------------
348! Compute basic state initial level thicknesses used for state norm
349! scaling. It uses zero free-surface (rest state). Therefore, the
350! norm scaling is time invariant.
351!-----------------------------------------------------------------------
352!
353 DO ng=1,ngrids
354 DO tile=last_tile(ng),first_tile(ng),-1
355 CALL set_depth (ng, tile, iadm)
356 END DO
357 END DO
358#endif
359!
360!-----------------------------------------------------------------------
361! Compute adjoint state vector.
362!-----------------------------------------------------------------------
363!
364 DO ng=1,ngrids
365 DO tile=first_tile(ng),last_tile(ng),+1
366 CALL ad_inner2state (ng, tile, lini, ad_state(ng)%vector)
367 END DO
368 END DO
369!
370 10 FORMAT (/,a,i2.2,a,i3.3,a,i3.3/)
371 20 FORMAT (/,a,i2.2,a,1p,e15.6,/)
372 30 FORMAT (/,1x,a,1x,'ROMS: started time-stepping:', &
373 & ' (Grid: ',i2.2,' TimeSteps: ',i8.8,' - ',i8.8,')')
374!
375 RETURN
376 END SUBROUTINE propagator_hop
377
378 END MODULE propagator_mod
subroutine ad_get_data(ng)
Definition ad_get_data.F:4
subroutine ad_get_idata(ng)
Definition ad_get_idata.F:4
subroutine ad_main2d
Definition ad_main2d.F:586
subroutine ad_main3d(runinterval)
Definition ad_main3d.F:4
subroutine, public close_inp(ng, model)
Definition close_io.F:92
subroutine, public tl_statenorm(ng, tile, kinp, ninp, statenorm)
Definition dotproduct.F:699
subroutine, public ad_ini_perturb(ng, tile, kinp, kout, ninp, nout)
subroutine, public tl_inner2state(ng, tile, lini, state)
Definition inner2state.F:63
subroutine, public ad_inner2state(ng, tile, lini, ad_state)
subroutine, public ini_c_norm(ng, tile, kinp, ninp, statenorm)
subroutine, public initialize_coupling(ng, tile, model)
integer stdout
integer, parameter r8
Definition mod_kinds.F:28
integer, parameter dp
Definition mod_kinds.F:25
subroutine, public initialize_ocean(ng, tile, model)
Definition mod_ocean.F:1526
integer, dimension(:), allocatable first_tile
logical master
integer, dimension(:), allocatable last_tile
integer, parameter inlm
Definition mod_param.F:662
integer, parameter iadm
Definition mod_param.F:665
integer ngrids
Definition mod_param.F:113
integer, parameter itlm
Definition mod_param.F:663
real(dp), parameter day2sec
integer, dimension(:), allocatable ntimes
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
integer, dimension(:), allocatable nconv
logical, dimension(:), allocatable synchro_flag
logical, dimension(:), allocatable predictor_2d_step
real(dp), dimension(:), allocatable tdays
real(dp) dstart
real(dp), parameter sec2day
integer, dimension(:), allocatable ntend
integer exit_flag
integer, dimension(:), allocatable indx1
integer, dimension(:), allocatable ntfirst
real(dp), dimension(:), allocatable time
integer, dimension(:), allocatable ntstart
integer nrun
integer, dimension(:), allocatable iif
integer noerror
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable knew
integer, dimension(:), allocatable nrhs
integer, dimension(:), allocatable nnew
integer, dimension(:), allocatable krhs
integer, dimension(:), allocatable nstp
subroutine, public propagator_hop(runinterval, state, ad_state)
subroutine, public set_depth(ng, tile, model)
Definition set_depth.F:34
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52
subroutine tl_get_data(ng)
Definition tl_get_data.F:4
subroutine tl_get_idata(ng)
Definition tl_get_idata.F:4
subroutine tl_main2d
Definition tl_main2d.F:429
subroutine tl_main3d(runinterval)
Definition tl_main3d.F:4