ROMS
Loading...
Searching...
No Matches
tl_main2d.F
Go to the documentation of this file.
1#include "cppdefs.h"
2#if defined TANGENT && !defined SOLVE3D
3 SUBROUTINE tl_main2d (RunInterval)
4!
5!git $Id$
6!================================================== Hernan G. Arango ===
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md !
10!=======================================================================
11! !
12! This routine is the main driver for tangent linear ROMS when !
13! configurated as shallow water (barotropic ) ocean model only. It !
14! advances forward the tangent linear model for all nested grids, !
15! if any, by the specified time interval (seconds), RunInterval. !
16! !
17!=======================================================================
18!
19 USE mod_param
20 USE mod_parallel
21# if defined MODEL_COUPLING && defined MCT_LIB
22 USE mod_coupler
23# endif
24 USE mod_iounits
25 USE mod_scalars
26 USE mod_stepping
27!
28 USE dateclock_mod, ONLY : time_string
29# ifdef TLM_CHECK
31# endif
32# ifdef TIDE_GENERATING_FORCES
33 USE equilibrium_tide_mod, ONLY : equilibrium_tide
34# endif
35# if defined ATM_COUPLING_NOT_YET && defined MCT_LIB
36 USE mct_coupler_mod, ONLY : ocn2atm_coupling
37# endif
38# if defined WAV_COUPLING_NOT_YET && defined MCT_LIB
39 USE mct_coupler_mod, ONLY : ocn2wav_coupling
40# endif
41 USE strings_mod, ONLY : founderror
42 USE tl_diag_mod, ONLY : tl_diag
43# if defined WEAK_CONSTRAINT || defined FORCING_SV
44 USE tl_forcing_mod, ONLY : tl_forcing
45# endif
46# ifdef ADJUST_WSTRESS
48# endif
50# ifdef ADJUST_BOUNDARY
52# endif
53# ifdef NEARSHORE_MELLOR_NOT_YET
54!! USE tl_radiation_stress_mod, ONLY : tl_radiation_stress
55# endif
56# ifdef TL_AVERAGES
57 USE tl_set_avg_mod, ONLY : tl_set_avg
58# endif
59# if defined SSH_TIDES_NOT_YET || defined UV_TIDES_NOT_YET
60!! USE tl_set_tides_mod, ONLY : tl_set_tides
61# endif
62 USE tl_set_vbc_mod, ONLY : tl_set_vbc
63 USE tl_step2d_mod, ONLY : tl_step2d
64# ifdef FLOATS_NOT_YET
65!! USE tl_step_floats_mod, ONLY : tl_step_floats
66# endif
67!
68 implicit none
69!
70! Imported variable declarations.
71!
72 real(dp), intent(in) :: RunInterval
73!
74! Local variable declarations.
75!
76 integer :: ng, tile
77 integer :: next_indx1
78# ifdef FLOATS_NOT_YET
79 integer :: Lend, Lstr, chunk_size
80# endif
81!
82 real(r8) :: MaxDT, my_StepTime
83!
84 character (len=*), parameter :: MyFile = &
85 & __FILE__
86!
87!=======================================================================
88! Time-step tangent linear vertically integrated equations.
89!=======================================================================
90!
91 my_steptime=0.0_r8
92 maxdt=maxval(dt)
93
94 step_loop : DO WHILE (my_steptime.le.(runinterval+0.5_r8*maxdt))
95
96 my_steptime=my_steptime+maxdt
97!
98! Set time clock.
99!
100 DO ng=1,ngrids
101 iic(ng)=iic(ng)+1
102!$OMP MASTER
103 time(ng)=time(ng)+dt(ng)
104 tdays(ng)=time(ng)*sec2day
105 CALL time_string (time(ng), time_code(ng))
106!$OMP END MASTER
107!$OMP BARRIER
108 END DO
109!
110!-----------------------------------------------------------------------
111! Read in required data, if any, from input NetCDF files.
112!-----------------------------------------------------------------------
113!
114 DO ng=1,ngrids
115!$OMP MASTER
116 CALL tl_get_data (ng)
117!$OMP END MASTER
118!$OMP BARRIER
119 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
120 END DO
121!
122!-----------------------------------------------------------------------
123! If applicable, process input data: time interpolate between data
124! snapshots.
125!-----------------------------------------------------------------------
126!
127 DO ng=1,ngrids
128 DO tile=first_tile(ng),last_tile(ng),+1
129 CALL tl_set_data (ng, tile)
130 END DO
131!$OMP BARRIER
132 END DO
133 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
134
135# if (defined WEAK_CONSTRAINT || defined FORCING_SV) && \
136 !defined SP4DVAR
137!
138!-----------------------------------------------------------------------
139! If appropriate, add convolved adjoint solution impulse forcing to
140! the representer model solution. Notice that the forcing is only
141! needed after finishing all inner loops. The forcing is continuous.
142! That is, it is time interpolated at every time-step from available
143! snapshots (FrequentImpulse=TRUE).
144!-----------------------------------------------------------------------
145!
146 DO ng=1,ngrids
147# ifdef WEAK_CONSTRAINT
148# ifdef WEAK_NOINTERP
149 IF ((iic(ng).gt.1).and.(iic(ng).ne.ntend(ng)+1).and. &
150 & (mod(iic(ng)-1,nadj(ng)).eq.0)) THEN
151 IF (master) THEN
152 WRITE (stdout,*) ' FORCING TLM at iic = ',iic(ng)
153 END IF
154# endif
155 IF (frequentimpulse(ng)) THEN
156 DO tile=first_tile(ng),last_tile(ng),+1
157 CALL tl_forcing (ng, tile, kstp(ng), nstp(ng))
158 END DO
159!$OMP BARRIER
160 END IF
161# ifdef WEAK_NOINTERP
162 END IF
163# endif
164# else
165 DO tile=first_tile(ng),last_tile(ng),+1
166 CALL tl_forcing (ng, tile, kstp(ng), nstp(ng))
167 END DO
168!$OMP BARRIER
169# endif
170 END DO
171# endif
172!
173!-----------------------------------------------------------------------
174! If not a restart, initialize all time levels and compute other
175! initial fields.
176!-----------------------------------------------------------------------
177!
178 DO ng=1,ngrids
179 IF (iic(ng).eq.ntstart(ng)) THEN
180!
181! Initialize free-surface.
182!
183 DO tile=first_tile(ng),last_tile(ng),+1
184 CALL tl_ini_zeta (ng, tile, itlm)
185 END DO
186!$OMP BARRIER
187!
188! Initialize other state variables.
189!
190 DO tile=last_tile(ng),first_tile(ng),-1
191 CALL tl_ini_fields (ng, tile, itlm)
192 END DO
193!$OMP BARRIER
194 END IF
195 END DO
196!
197!-----------------------------------------------------------------------
198! Compute and report diagnostics. If appropriate, accumulate time-
199! averaged output data which needs a irreversible loop in shared-memory
200! jobs.
201!-----------------------------------------------------------------------
202!
203 DO ng=1,ngrids
204 DO tile=first_tile(ng),last_tile(ng),+1 ! irreversible
205# ifdef TL_AVERAGES
206 CALL tl_set_avg (ng, tile)
207# endif
208# ifdef DIAGNOSTICS
209!! CALL set_diags (ng, tile)
210# endif
211# ifdef TIDE_GENERATING_FORCES
212 CALL equilibrium_tide (ng, tile, itlm)
213# endif
214 CALL tl_diag (ng, tile)
215# ifdef TLM_CHECK
216 CALL tl_dotproduct (ng, tile, lnew(ng))
217# endif
218 END DO
219!$OMP BARRIER
220 END DO
221
222# if defined ATM_COUPLING_NOT_YET && defined MCT_LIB
223!
224!-----------------------------------------------------------------------
225! Couple to atmospheric model every CoupleSteps(Iatmos) timesteps: get
226! air/sea fluxes.
227!-----------------------------------------------------------------------
228!
229 DO ng=1,ngrids
230 IF ((iic(ng).ne.ntstart(ng)).and. &
231 & mod(iic(ng)-1,couplesteps(iatmos,ng)).eq.0) THEN
232 DO tile=last_tile(ng),first_tile(ng),-1
233 CALL ocn2atm_coupling (ng, tile)
234 END DO
235!$OMP BARRIER
236 END IF
237 END DO
238# endif
239
240# ifdef ADJUST_BOUNDARY
241!
242!-----------------------------------------------------------------------
243! Interpolate open boundary increments and adjust open boundaries.
244! Skip the last output timestep.
245!-----------------------------------------------------------------------
246!
247 DO ng=1,ngrids
248 IF (iic(ng).lt.(ntend(ng)+1)) THEN
249 DO tile=first_tile(ng),last_tile(ng),+1
250 CALL tl_obc_adjust (ng, tile, lbinp(ng))
251 END DO
252!$OMP BARRIER
253 END IF
254 END DO
255# endif
256
257# ifdef ADJUST_WSTRESS
258!
259!-----------------------------------------------------------------------
260! Interpolate surface forcing increments and adjust surface forcing.
261! Skip the last output timestep.
262!-----------------------------------------------------------------------
263!
264 DO ng=1,ngrids
265 IF (iic(ng).lt.(ntend(ng)+1)) THEN
266 DO tile=first_tile(ng),last_tile(ng),+1
267 CALL tl_frc_adjust (ng, tile, lfinp(ng))
268 END DO
269!$OMP BARRIER
270 END IF
271 END DO
272# endif
273
274# if defined WAV_COUPLING_NOT_YET && defined MCT_LIB
275!
276!-----------------------------------------------------------------------
277! Couple to waves model every CoupleSteps(Iwaves) timesteps: get
278! waves/sea fluxes.
279!-----------------------------------------------------------------------
280!
281 DO ng=1,ngrids
282 IF ((iic(ng).ne.ntstart(ng)).and. &
283 & mod(iic(ng)-1,couplesteps(iwaves,ng)).eq.0) THEN
284 DO tile=first_tile(ng),last_tile(ng),+1
285 CALL ocn2wav_coupling (ng, tile)
286 END DO
287!$OMP BARRIER
288 END IF
289 END DO
290# endif
291
292# ifdef NEARSHORE_MELLOR_NOT_YET
293!
294!-----------------------------------------------------------------------
295! Compute radiation stress terms.
296!-----------------------------------------------------------------------
297!
298 DO ng=1,ngrids
299 DO tile=last_tile(ng),first_tile(ng),-1
300 CALL tl_radiation_stress (ng, tile)
301 END DO
302!$OMP BARRIER
303 END DO
304# endif
305!
306!-----------------------------------------------------------------------
307! Set vertical boundary conditions. Process tidal forcing.
308!-----------------------------------------------------------------------
309!
310 DO ng=1,ngrids
311 DO tile=first_tile(ng),last_tile(ng),+1
312 CALL tl_set_vbc (ng, tile)
313# if defined SSH_TIDES_NOT_YET || defined UV_TIDES_NOT_YET
314 CALL tl_set_tides (ng, tile)
315# endif
316 END DO
317!$OMP BARRIER
318 END DO
319!
320!-----------------------------------------------------------------------
321! If appropriate, write out fields into output NetCDF files. Notice
322! that IO data is written in delayed and serial mode. Exit if last
323! time step.
324!-----------------------------------------------------------------------
325!
326 DO ng=1,ngrids
327!$OMP MASTER
328 CALL tl_output (ng)
329!$OMP END MASTER
330!$OMP BARRIER
331 IF ((founderror(exit_flag, noerror, __line__, myfile)).or. &
332 & ((iic(ng).eq.(ntend(ng)+1)).and.(ng.eq.ngrids))) RETURN
333 END DO
334!
335!-----------------------------------------------------------------------
336! Solve the vertically integrated primitive equations for the
337! free-surface and momentum components.
338!-----------------------------------------------------------------------
339!
340! Set time indices for predictor step. The PREDICTOR_2D_STEP switch
341! it is assumed to be false before the first time-step.
342!
343 DO ng=1,ngrids
344 iif(ng)=1
345 nfast(ng)=1
346 next_indx1=3-indx1(ng)
347 IF (.not.predictor_2d_step(ng)) THEN
348 predictor_2d_step(ng)=.true.
349 IF (first_2d_step) THEN
350 kstp(ng)=indx1(ng)
351 ELSE
352 kstp(ng)=3-indx1(ng)
353 END IF
354 knew(ng)=3
355 krhs(ng)=indx1(ng)
356 END IF
357!
358! Predictor step - Advance barotropic equations using 2D time-step
359! ============== predictor scheme.
360!
361 DO tile=last_tile(ng),first_tile(ng),-1
362 CALL tl_step2d (ng, tile)
363 END DO
364!$OMP BARRIER
365 END DO
366!
367! Set time indices for corrector step.
368!
369 DO ng=1,ngrids
370 IF (predictor_2d_step(ng)) THEN
371 predictor_2d_step(ng)=.false.
372 knew(ng)=next_indx1
373 kstp(ng)=3-knew(ng)
374 krhs(ng)=3
375 IF (iif(ng).lt.(nfast(ng)+1)) indx1(ng)=next_indx1
376 END IF
377!
378! Corrector step - Apply 2D time-step corrector scheme. Notice that
379! ============== there is not need for a corrector step during the
380! auxiliary (nfast+1) time-step.
381!
382 IF (iif(ng).lt.(nfast(ng)+1)) THEN
383 DO tile=first_tile(ng),last_tile(ng),+1
384 CALL tl_step2d (ng, tile)
385 END DO
386!$OMP BARRIER
387 END IF
388 END DO
389
390# ifdef FLOATS_NOT_YET
391!
392!-----------------------------------------------------------------------
393! Compute Lagrangian drifters trajectories: Split all the drifters
394! between all the computational threads, except in distributed-memory
395! and serial configurations. In distributed-memory, the parallel node
396! containing the drifter is selected internally since the state
397! variables do not have a global scope.
398!-----------------------------------------------------------------------
399!
400 DO ng=1,ngrids
401 IF (lfloats(ng)) THEN
402# ifdef _OPENMP
403 chunk_size=(nfloats(ng)+numthreads-1)/numthreads
404 lstr=1+my_thread*chunk_size
405 lend=min(nfloats(ng),lstr+chunk_size-1)
406# else
407 lstr=1
408 lend=nfloats(ng)
409# endif
410 CALL tl_step_floats (ng, lstr, lend)
411!$OMP BARRIER
412!
413! Shift floats time indices.
414!
415 nfp1(ng)=mod(nfp1(ng)+1,nft+1)
416 nf(ng)=mod(nf(ng)+1,nft+1)
417 nfm1(ng)=mod(nfm1(ng)+1,nft+1)
418 nfm2(ng)=mod(nfm2(ng)+1,nft+1)
419 nfm3(ng)=mod(nfm3(ng)+1,nft+1)
420 END IF
421 END DO
422# endif
423 END DO step_loop
424!
425 RETURN
426 END SUBROUTINE tl_main2d
427#else
428 SUBROUTINE tl_main2d
429 RETURN
430 END SUBROUTINE tl_main2d
431#endif
subroutine, public time_string(mytime, date_string)
Definition dateclock.F:1272
subroutine, public tl_dotproduct(ng, tile, linp)
integer, dimension(:,:), allocatable couplesteps
integer iatmos
integer iwaves
integer stdout
integer numthreads
integer, dimension(:), allocatable first_tile
logical master
integer, dimension(:), allocatable last_tile
integer, dimension(:), allocatable nfloats
Definition mod_param.F:543
integer, parameter nft
Definition mod_param.F:539
integer ngrids
Definition mod_param.F:113
integer, parameter itlm
Definition mod_param.F:663
logical, dimension(:), allocatable lfloats
integer, dimension(:), allocatable iic
real(dp), dimension(:), allocatable dt
logical, dimension(:), allocatable predictor_2d_step
real(dp), dimension(:), allocatable tdays
logical, dimension(:), allocatable frequentimpulse
integer, dimension(:), allocatable nfast
real(dp), parameter sec2day
integer, dimension(:), allocatable ntend
character(len=22), dimension(:), allocatable time_code
integer exit_flag
integer, dimension(:), allocatable indx1
real(dp), dimension(:), allocatable time
integer, dimension(:), allocatable ntstart
integer, dimension(:), allocatable nadj
integer, dimension(:), allocatable iif
integer noerror
integer, dimension(:), allocatable nfm2
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable knew
integer, dimension(:), allocatable nfm1
integer, dimension(:), allocatable lfinp
integer, dimension(:), allocatable lbinp
integer, dimension(:), allocatable nf
integer, dimension(:), allocatable nfm3
integer, dimension(:), allocatable nfp1
integer, dimension(:), allocatable lnew
integer, dimension(:), allocatable krhs
integer, dimension(:), allocatable nstp
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52
subroutine, public tl_diag(ng, tile)
Definition tl_diag.F:25
subroutine, public tl_forcing(ng, tile, kfrc, nfrc)
Definition tl_forcing.F:26
subroutine, public tl_frc_adjust(ng, tile, linp)
subroutine, public tl_ini_fields(ng, tile, model)
subroutine, public tl_ini_zeta(ng, tile, model)
subroutine, public tl_obc_adjust(ng, tile, linp)
subroutine, public tl_set_avg(ng, tile)
Definition tl_set_avg.F:27
subroutine, public tl_set_vbc(ng, tile)
Definition tl_set_vbc.F:33
subroutine, public tl_step2d(ng, tile)
subroutine tl_get_data(ng)
Definition tl_get_data.F:4
subroutine tl_main2d
Definition tl_main2d.F:429
subroutine tl_output(ng)
Definition tl_output.F:4
subroutine tl_set_data(ng, tile)
Definition tl_set_data.F:4