ROMS
Loading...
Searching...
No Matches
ntimestep.F
Go to the documentation of this file.
1#include "cppdefs.h"
2 SUBROUTINE ntimesteps (model, RunInterval, nl, Nsteps, Rsteps)
3!
4!git $Id$
5!=======================================================================
6! Copyright (c) 2002-2025 The ROMS Group !
7! Licensed under a MIT/X style license !
8! See License_ROMS.md !
9!================================================== Hernan G. Arango ===
10! !
11! This routine set the number of time-steps to compute. In nesting !
12! applications, the number of time-steps depends on nesting layer !
13! configuration. !
14! !
15! On Input: !
16! !
17! model Calling model identifier (integer) !
18! RunInterval Time interval (seconds) to advance forward or !
19! backwards the primitive equations (scalar) !
20! nl Nesting layer number (integer) !
21! !
22! On Input: !
23! !
24! nl Updated nesting layer number (integer) !
25! Nsteps Number of time-steps to solve for all grids in !
26! current nesting layer "nl" (integer) !
27! Rsteps Number of time-steps to complete RunInterval !
28! time window (integer) !
29! !
30!=======================================================================
31!
32 USE mod_param
33 USE mod_parallel
34 USE mod_iounits
35#ifdef NESTING
36 USE mod_nesting
37#endif
38 USE mod_scalars
39!
40! Imported variable declarations.
41!
42 integer, intent(in) :: model
43 integer, intent(inout) :: nl
44 integer, intent(out) :: Nsteps, Rsteps
45!
46 real(dp), intent(in) :: RunInterval
47!
48! Local variable declarations.
49!
50 integer :: extra, gn, ig, il, ng, ngm1
51 integer, dimension(Ngrids) :: WindowSteps, my_Nsteps
52
53#if defined MODEL_COUPLING && !defined MCT_LIB
54!
55 real(dp) :: ENDtime, NEXTtime
56#endif
57!
58!=======================================================================
59! Set number of time steps to execute for grid "ng".
60!=======================================================================
61!
62! Initialize.
63!
64 windowsteps=0
65 my_nsteps=0
66 nsteps=0
67 rsteps=0
68!
69! If appropriate, advance the nested layer counter.
70!
71 nl=nl+1
72
73#ifdef NESTING
74!
75! If refinement and telescoping refined grids, reset the nested layer
76! counter to advance all the grid to the current time of the coarser
77! of all grids (ng=1).
78!
79 IF (any(telescoping).and.(nl.gt.nestlayers)) THEN
80 DO il=1,nestlayers
81 DO ig=1,gridsinlayer(il)
82 ng=gridnumber(ig,il)
83 IF (telescoping(ng).and. &
84 & (refinestepscounter(ng).lt.refinesteps(ng))) THEN
85 nl=il
86 END IF
87 END DO
88 END DO
89 END IF
90#endif
91!
92 IF ((0.lt.nl).and.(nl.le.nestlayers)) THEN
93!
94! In composite and mosaic grids or all grids in the same nesting
95! layer, it is assumed that the donor and receiver grids have the
96! same time-step size. This is done to avoid the time interpolation
97! between donor and receiver grids. Only spatial interpolations
98! are possible in the current nesting design.
99!
100! In grid refinement, it is assumed that the donor and receiver grids
101! are an interger factor of the grid size and time-step size.
102!
103 windowsteps=0
104!
105! Loop over all grids in current layer nesting layer.
106!
107 DO ig=1,gridsinlayer(nl)
108 ng=gridnumber(ig,nl)
109
110#if defined MODEL_COUPLING && !defined MCT_LIB
111!
112! Set extra step parameter needed to finish the simulation due to ROMS
113! delayed output until the next half step. If RunInterval (seconds) is
114! less than full simulation interval because of model coupling, extra=1
115! for last coupling interval and extra=0 otherwise.
116!
117 nexttime=time(ng)+runinterval
118 endtime=initime(ng)+(ntimes(ng)-1)*dt(ng)
119 IF (nexttime.eq.endtime) THEN
120 extra=1
121 ELSE
122 extra=0
123 END IF
124
125#elif defined JEDI
126!
127! OOPS will advance ROMS kernels by smaller intervals, usually a single
128! timestep. The strategy here is different from that used for coupling
129! since OOPS controls every NLM, TLM, and ADM timestep. The exchange of
130! ROMS-to-JEDI state fields is also done one at a time and, in a delay
131! mode, to account for the modification to initial conditions and an
132! extra half step to complete the solution.
133!
134 IF (model.eq.inlm) THEN
135 nexttime=time(ng)+runinterval
136 endtime=initime(ng)+ntimes(ng)*dt(ng)
137 extra=0
138 ELSE IF (model.eq.itlm) THEN
139 nexttime=time(ng)+runinterval
140 endtime=initime(ng)+ntimes(ng)*dt(ng)
141 extra=0
142 ELSE IF (model.eq.iadm) THEN
143 nexttime=time(ng)-runinterval
144 endtime=initime(ng)+ntimes(ng)*dt(ng)
145 extra=0
146 END IF
147#else
148!
149! Here, extra=1, indicates that the RunInterval is the same as
150! simulation interval.
151!
152 extra=1
153#endif
154!
155! Determine number of steps in time-interval window.
156!
157 windowsteps(ig)=int((runinterval+0.5_r8*dt(ng))/dt(ng))
158!
159! Advancing model forward: Nonlinear, tangent linear, and representer
160! models.
161!
162 IF ((model.eq.inlm).or. &
163 & (model.eq.itlm).or. &
164 & (model.eq.irpm)) THEN
165
166#ifdef NESTING
167 IF (any(compositegrid(:,ng))) THEN
168 IF (step_counter(ng).le.(windowsteps(ig)+extra)) THEN
169 my_nsteps(ig)=1
170 step_counter(ng)=step_counter(ng)+1
171 ELSE
172 my_nsteps(ig)=0
173 END IF
174 ELSE IF (refinedgrid(ng).and.(refinescale(ng).eq.0)) THEN
175 refinestepscounter=0 ! The coarser grid, reset counters
176 IF (step_counter(ng).le.(windowsteps(ig)+extra)) THEN
177 my_nsteps(ig)=1
178 step_counter(ng)=step_counter(ng)+1
180 ELSE
181 my_nsteps(ig)=0
182 END IF
183 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
184 IF (step_counter(ng).le.(windowsteps(ig)+extra)) THEN
185 IF (telescoping(ng)) THEN
186 my_nsteps(ig)=1
187 step_counter(ng)=step_counter(ng)+1
189 DO il=nl+1,nestlayers ! When a parent steps,
190 gn=gridnumber(ig,il) ! set all its telescoped
191 IF (telescoping(gn)) THEN ! children counters to
192 refinestepscounter(gn)=0 ! zero
193 END IF
194 END DO
195 ELSE
196 my_nsteps(ig)=refinesteps(ng)
199 & refinesteps(ng)
200 END IF
201 ELSE
202 my_nsteps(ig)=0
203 END IF
204 END IF
205#else
206 my_nsteps(ig)=max(my_nsteps(ig), windowsteps(ig)+extra)
207 step_counter(ng)=step_counter(ng)+windowsteps(ig)+extra
208#endif
209
210#if defined ADJOINT && !defined NESTING
211!
212! Advancing model backwards: Adjoint model.
213!
214 ELSE IF (model.eq.iadm) THEN
215 my_nsteps(ig)=max(my_nsteps(ig), windowsteps(ig)+extra)
216 step_counter(ng)=step_counter(ng)+windowsteps(ig)+extra
217#endif
218 END IF
219 END DO
220!
221! Insure that the steps per time-window are the same.
222!
223 IF (gridsinlayer(nl).gt.1) THEN
224 DO ig=2,gridsinlayer(nl)
225 IF (windowsteps(ig).ne.windowsteps(ig-1)) THEN
226 ngm1=gridnumber(ig-1,nl)
227 ng =gridnumber(ig ,nl)
228 IF (master) THEN
229 WRITE (stdout,10) nl, ngm1, dt(ngm1), ng, dt(ng)
230 10 FORMAT (/,' NTIMESTEPS - timestep size are not the ', &
231 & ' same in nesting layer: ',i2, &
232 & 2(/,14x,'Grid ',i2.2,3x,'dt = ',f11.3))
233 END IF
234 exit_flag=5
235 END IF
236 END DO
237 END IF
238!
239! Set number of time-steps to execute. Choose minimum values.
240!
241 nsteps=my_nsteps(1)
242 rsteps=windowsteps(1)+extra
243 DO ig=2,gridsinlayer(nl)
244 nsteps=min(nsteps, my_nsteps(ig))
245 rsteps=min(rsteps, windowsteps(ig)+extra)
246 END DO
247 END IF
248
249 RETURN
250 END SUBROUTINE ntimesteps
251
252#if defined NESTING && defined ADJOINT
253!
254 SUBROUTINE nlm_step_sequence (RunInterval, icount)
255!
256!git $Id$
257!================================================== Hernan G. Arango ===
258! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
259! Licensed under a MIT/X style license !
260! See License_ROMS.md !
261!=======================================================================
262! !
263! Given a time interval window (seconds), this routine stores all the !
264! parameters needed to timestep the adjoint model backwards in nested !
265! applications. It is an elegant way to have the correct parameters !
266! in complex setups. !
267! !
268! StepInfo(:,1) NL, current nested layer number !
269! StepInfo(:,2) Nsteps, number of timesteps for grids in NL !
270! StepInfo(:,3) Rsteps, number of timesteps for simulation window !
271! StepInfo(:,ng+3) step_counter(ng), timestep counter for each grid !
272! !
273! On Input: !
274! !
275! RunInterval Time interval window (seconds) to advance forward !
276! or backwards the primitive equations (scalar) !
277! !
278! On Input: !
279! !
280! icount Number of stored values in "StepInfo" (integer) !
281! StepInfo(icount,Ngrid+3) !
282! !
283! It calls repetitively routine "ntimestep". !
284! !
285!=======================================================================
286!
287 USE mod_param
288 USE mod_parallel
289 USE mod_iounits
290 USE mod_nesting
291 USE mod_scalars
292!
293! Imported variable declarations.
294!
295 integer, intent(out) :: icount
296!
297 real(dp), intent(in) :: RunInterval
298!
299! Local variable declarations.
300!
301 logical :: DoNestLayer, Time_Step
302 integer :: NL, Nsteps, Rsteps
303 integer :: dg, ig, istep, ng
304!
305!=======================================================================
306! Compute the NLM sequence of nested timesteps parameters.
307!=======================================================================
308!
309! If applicable, allocate and intialize parameter sequence array. It
310! can only be allocated in "ROMS_run" phase when the RunInterval value
311! is known, and niot before. That value is specified by the Master or
312! Coupler subroutine. Notice that we use the (rs-1) factor to discount
313! the coarser grid step that other grids need to execute to arrive to
314! the same simulation time.
315!
316 IF (.not.allocated(stepinfo)) THEN
317 dg=1
318 rs=1
319 DO ng=1,ngrids
320 IF (refinescale(ng).gt.0) THEN
321 rs=rs*refinescale(ng)
322 END IF
323 END DO
324 nstepinfo=int((runinterval+0.5_r8*dt(dg))/dt(dg))*(rs-1)+ &
325 & ngrids
326 allocate ( stepinfo(nstepinfo,ngrids+3) )
327 stepinfo=0
328 END IF
329!
330! Initialize.
331!
332 time_step=.true.
333 donestlayer=.true.
334 icount=0
335!
336! Compute and store nonlinear model sequence of NL, Nstep, Rstep, and
337! step_counter parameters. The sequence of values will be reversed
338! when timestepping the adjoint model backwards.
339!
340 DO WHILE (time_step)
341 nl=0
342 DO WHILE (donestlayer)
343 CALL ntimesteps (inlm, runinterval, nl, nsteps, rsteps)
344 icount=icount+1
345 IF (icount.le.nstepinfo) THEN
346 stepinfo(icount,1)=nl
347 stepinfo(icount,2)=nsteps
348 stepinfo(icount,3)=rsteps
349 DO ng=1,ngrids
350 stepinfo(icount,ng+3)=step_counter(ng)
351 END DO
352!
353 IF ((nl.le.0).or.(nl.gt.nestlayers)) EXIT
354!
355 DO istep=1,nsteps
356 DO ig=1,gridsinlayer(nl)
357 ng=gridnumber(ig,nl)
358 IF (step_counter(ng).eq.rsteps) time_step=.false.
359 END DO
360 END DO
361 IF (.not.time_step.and.nl.eq.nestlayers) EXIT
362 ELSE
363 IF (master) WRITE(stdout,10) ', nStepInfo = ', nstepinfo, &
364 & icount
365 exit_flag=5
366 RETURN
367 END IF
368 END DO
369 END DO
370!
371 10 FORMAT (/,' NLM_STEP_SEQUENCE - too small dimension parameter ', &
372 & a,i0,3x,i0,/,21x,'Recompute dimensions for ''StepInfo''')
373!
374 RETURN
375 END SUBROUTINE nlm_step_sequence
376#endif
integer stdout
integer, dimension(:), allocatable refinestepscounter
integer, dimension(:), allocatable refinesteps
logical, dimension(:), allocatable telescoping
logical master
integer, parameter inlm
Definition mod_param.F:662
integer, parameter irpm
Definition mod_param.F:664
integer, dimension(:,:), allocatable gridnumber
Definition mod_param.F:127
integer, parameter iadm
Definition mod_param.F:665
integer nestlayers
Definition mod_param.F:118
integer ngrids
Definition mod_param.F:113
integer, parameter itlm
Definition mod_param.F:663
integer, dimension(:), allocatable gridsinlayer
Definition mod_param.F:122
integer, dimension(:), allocatable ntimes
real(dp), dimension(:), allocatable dt
integer exit_flag
logical, dimension(:,:), allocatable compositegrid
real(dp), dimension(:), allocatable time
logical, dimension(:), allocatable refinedgrid
integer, dimension(:), allocatable refinescale
integer, dimension(:), allocatable step_counter
real(dp), dimension(:), allocatable initime
subroutine nlm_step_sequence(runinterval, icount)
Definition ntimestep.F:255
subroutine ntimesteps(model, runinterval, nl, nsteps, rsteps)
Definition ntimestep.F:3