ROMS
Loading...
Searching...
No Matches
optobs_roms.h
Go to the documentation of this file.
1 MODULE roms_kernel_mod
2!
3!git $Id$
4!================================================== Hernan G. Arango ===
5! Copyright (c) 2002-2025 The ROMS Group W. G. Zhang !
6! Licensed under a MIT/X style license !
7! See License_ROMS.md !
8!=======================================================================
9! !
10! ROMS Optimal Observation Driver: !
11! !
12! This driver computes the adjoint sensitivity of a function or !
13! index, J, in terms of space and/or time integrals of the model !
14! state, S(zeta,u,v,T,...). Small changes, dS, in S will lead to !
15! changes dJ in J: !
16! !
17! dJ = (dJ/dzeta) dzeta + (dJ/du) du + (dJ/dv) dv + (dJ/dt) dT ... !
18! !
19! and !
20! !
21! dJ/dS = transpose(R) S !
22! !
23! where transpose(R) is the adjoint propagator. It implies that !
24! the sensitivity for ALL variables, parameters, and space-time !
25! points can be computed from a single integration of the adjoint !
26! model. !
27! !
28! These routines control the initialization, time-stepping, and !
29! finalization of ROMS model following ESMF conventions: !
30! !
31! ROMS_initialize !
32! ROMS_run !
33! ROMS_finalize !
34! !
35! Reference: !
36! !
37! !
38!=======================================================================
39!
40 USE mod_param
41 USE mod_parallel
42 USE mod_arrays
43 USE mod_fourdvar
44 USE mod_iounits
45 USE mod_ncparam
46 USE mod_ocean
47 USE mod_scalars
48 USE mod_stepping
49!
52 USE def_norm_mod, ONLY : def_norm
53 USE inp_par_mod, ONLY : inp_par
54 USE get_state_mod, ONLY : get_state
55#ifdef MCT_LIB
56# ifdef ATM_COUPLING
57 USE mct_coupler_mod, ONLY : initialize_ocn2atm_coupling
58# endif
59# ifdef WAV_COUPLING
60 USE mct_coupler_mod, ONLY : initialize_ocn2wav_coupling
61# endif
62#endif
65 USE strings_mod, ONLY : founderror
66 USE tl_def_ini_mod, ONLY : tl_def_ini
67 USE wrt_rst_mod, ONLY : wrt_rst
68#if defined BALANCE_OPERATOR && defined ZETA_ELLIPTIC
70#endif
71!
72 implicit none
73!
74 PUBLIC :: roms_initialize
75 PUBLIC :: roms_run
76 PUBLIC :: roms_finalize
77!
78 CONTAINS
79!
80 SUBROUTINE roms_initialize (first, mpiCOMM)
81!
82!=======================================================================
83! !
84! This routine allocates and initializes ROMS state variables !
85! and internal and external parameters. !
86! !
87!=======================================================================
88!
89! Imported variable declarations.
90!
91 logical, intent(inout) :: first
92!
93 integer, intent(in), optional :: mpiCOMM
94!
95! Local variable declarations.
96!
97 logical :: allocate_vars = .true.
98!
99#ifdef DISTRIBUTE
100 integer :: MyError, MySize
101#endif
102 integer :: STDrec, Tindex
103 integer :: chunk_size, ng, thread
104#ifdef _OPENMP
105 integer :: my_threadnum
106#endif
107!
108 character (len=*), parameter :: MyFile = &
109 & __FILE__//", ROMS_initialize"
110
111#ifdef DISTRIBUTE
112!
113!-----------------------------------------------------------------------
114! Set distribute-memory (mpi) world communictor.
115!-----------------------------------------------------------------------
116!
117 IF (PRESENT(mpicomm)) THEN
118 ocn_comm_world=mpicomm
119 ELSE
120 ocn_comm_world=mpi_comm_world
121 END IF
122 CALL mpi_comm_rank (ocn_comm_world, myrank, myerror)
123 CALL mpi_comm_size (ocn_comm_world, mysize, myerror)
124#endif
125!
126!-----------------------------------------------------------------------
127! On first pass, initialize model parameters a variables for all
128! nested/composed grids. Notice that the logical switch "first"
129! is used to allow multiple calls to this routine during ensemble
130! configurations.
131!-----------------------------------------------------------------------
132!
133 IF (first) THEN
134 first=.false.
135!
136! Initialize parallel control switches. These scalars switches are
137! independent from standard input parameters.
138!
140!
141! Set the ROMS standard output unit to write verbose execution info.
142! Notice that the default standard out unit in Fortran is 6.
143!
144! In some applications like coupling or disjointed mpi-communications,
145! it is advantageous to write standard output to a specific filename
146! instead of the default Fortran standard output unit 6. If that is
147! the case, it opens such formatted file for writing.
148!
149 IF (set_stdoutunit) THEN
151 set_stdoutunit=.false.
152 END IF
153!
154! Read in model tunable parameters from standard input. Allocate and
155! initialize variables in several modules after the number of nested
156! grids and dimension parameters are known.
157!
158 CALL inp_par (iadm)
159 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
160!
161! Set domain decomposition tile partition range. This range is
162! computed only once since the "first_tile" and "last_tile" values
163! are private for each parallel thread/node.
164!
165#if defined _OPENMP
166 mythread=my_threadnum()
167#elif defined DISTRIBUTE
169#else
170 mythread=0
171#endif
172 DO ng=1,ngrids
173 chunk_size=(ntilex(ng)*ntilee(ng)+numthreads-1)/numthreads
174 first_tile(ng)=mythread*chunk_size
175 last_tile(ng)=first_tile(ng)+chunk_size-1
176 END DO
177!
178! Initialize internal wall clocks. Notice that the timings does not
179! includes processing standard input because several parameters are
180! needed to allocate clock variables.
181!
182 IF (master) THEN
183 WRITE (stdout,10)
184 10 FORMAT (/,' Process Information:',/)
185 END IF
186!
187 DO ng=1,ngrids
188 DO thread=thread_range
189 CALL wclock_on (ng, iadm, 0, __line__, myfile)
190 END DO
191 END DO
192!
193! Allocate and initialize all model state arrays.
194!
195 CALL roms_allocate_arrays (allocate_vars)
197 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
198
199 END IF
200
201#if defined MCT_LIB && (defined ATM_COUPLING || defined WAV_COUPLING)
202!
203!-----------------------------------------------------------------------
204! Initialize coupling streams between model(s).
205!-----------------------------------------------------------------------
206!
207 DO ng=1,ngrids
208# ifdef ATM_COUPLING
209 CALL initialize_ocn2atm_coupling (ng, myrank)
210# endif
211# ifdef WAV_COUPLING
212 CALL initialize_ocn2wav_coupling (ng, myrank)
213# endif
214 END DO
215#endif
216!
217!-----------------------------------------------------------------------
218! Read in standard deviation factors for error covariance.
219!-----------------------------------------------------------------------
220!
221! Initial conditions standard deviation. They are loaded in Tindex=1
222! of the e_var(...,Tindex) state variables.
223!
224 stdrec=1
225 tindex=1
226 DO ng=1,ngrids
227 CALL get_state (ng, 10, 10, std(1,ng), stdrec, tindex)
228 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
229 END DO
230!
231! Model error standard deviation. They are loaded in Tindex=2
232! of the e_var(...,Tindex) state variables.
233!
234 stdrec=1
235 tindex=2
236 DO ng=1,ngrids
237 IF (nsa.eq.2) THEN
238 CALL get_state (ng, 11, 11, std(2,ng), stdrec, tindex)
239 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
240 END IF
241 END DO
242
243#ifdef ADJUST_BOUNDARY
244!
245! Open boundary conditions standard deviation.
246!
247 stdrec=1
248 tindex=1
249 DO ng=1,ngrids
250 CALL get_state (ng, 12, 12, std(3,ng), stdrec, tindex)
251 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
252 END DO
253#endif
254#if defined ADJUST_WSTRESS || defined ADJUST_STFLUX
255!
256! Surface forcing standard deviation.
257!
258 stdrec=1
259 tindex=1
260 DO ng=1,ngrids
261 CALL get_state (ng, 13, 13, std(4,ng), stdrec, tindex)
262 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
263 END DO
264#endif
265!
266 RETURN
267 END SUBROUTINE roms_initialize
268!
269 SUBROUTINE roms_run (RunInterval)
270!
271!=======================================================================
272! !
273! This routine computes the adjoint sensitivity analysis, dJ/dS, !
274! to the specified functional J. The sensitivity masking arrays !
275! Rscope, Uscope, and Vscope are used to evaluate the functional !
276! in the desired spatial area. !
277! !
278!=======================================================================
279!
280! Imported variable declarations
281!
282 real(dp), intent(in) :: RunInterval ! seconds
283!
284! Local variable declarations.
285!
286 logical :: Lposterior
287!
288 integer :: ng, tile
289 integer :: Lbck, Lini, Rec, Rec1, Rec2
290 integer :: NRMrec
291!
292 real (r8) :: str_day, end_day
293!
294 character (len=6) :: driver
295
296 character (len=*), parameter :: MyFile = &
297 & __FILE__//", ROMS_run"
298!
299!=======================================================================
300! Run model for all nested grids, if any.
301!=======================================================================
302!
303! Initialize relevant parameters.
304!
305 DO ng=1,ngrids
306 lold(ng)=1 ! old minimization time index
307 lnew(ng)=2 ! new minimization time index
308 lreadfwd(ng)=.true.
309 END DO
310 lini=1 ! NLM initial conditions record in INI
311 lbck=1 ! background record in INI
312 rec1=1
313 rec2=2
314 driver='optobs'
315
316#ifdef FORWARD_FLUXES
317!
318! Set the BLK structure to contain the nonlinear model surface fluxes
319! needed by the tangent linear and adjoint models. Also, set switches
320! to process that structure in routine "check_multifile". Notice that
321! it is possible to split the solution into multiple NetCDF files to
322! reduce their size.
323!
324! The switch LreadFRC is deactivated because all the atmospheric
325! forcing, including shortwave radiation, is read from the NLM
326! surface fluxes or is assigned during ESM coupling. Such fluxes
327! are available from the QCK structure. There is no need for reading
328! and processing from the FRC structure input forcing-files.
329!
330 CALL edit_multifile ('QCK2BLK')
331 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
332 DO ng=1,ngrids
333 lreadblk(ng)=.true.
334 lreadfrc(ng)=.false.
335 END DO
336#endif
337#ifdef BALANCE_OPERATOR
338!
339! Set NLM model background trajectory to process in the balance
340! operator.
341!
342 DO ng=1,ngrids
343 ini(ng)%name=fwd(ng)%name
344 END DO
345
346# ifdef ZETA_ELLIPTIC
347!
348! Compute the reference zeta and biconjugate gradient arrays
349! required for the balance of free surface.
350!
351 IF (balance(isfsur)) THEN
352 DO ng=1,ngrids
353 CALL get_state (ng, inlm, 2, ini(ng), lini, lini)
354 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
355!
356 DO tile=first_tile(ng),last_tile(ng),+1
357 CALL balance_ref (ng, tile, lini)
358 CALL biconj (ng, tile, inlm, lini)
359 END DO
360 wrtzetaref(ng)=.true.
361 END DO
362 END IF
363# endif
364#endif
365!
366!-----------------------------------------------------------------------
367! Initialize adjoint model and define sensitivity functional.
368!-----------------------------------------------------------------------
369!
370 DO ng=1,ngrids
371 lstiffness=.false.
372 CALL ad_initial (ng)
373 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
374 END DO
375!
376! Activate adjoint output.
377!
378 DO ng=1,ngrids
379 ldefadj(ng)=.true.
380 lwrtadj(ng)=.true.
381 lcycleadj(ng)=.false.
382 END DO
383!
384! Define tangent linear initial conditions file.
385!
386 DO ng=1,ngrids
387 ldefitl(ng)=.true.
388 CALL tl_def_ini (ng)
389 ldefitl(ng)=.false.
390 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
391 END DO
392!
393!-----------------------------------------------------------------------
394! Compute or read in background-error correlations normalization
395! factors.
396!-----------------------------------------------------------------------
397!
398! If computing, write out factors to NetCDF. This is an expensive
399! computation and needs to be computed once for a particular
400! application grid.
401!
402 DO ng=1,ngrids
403 IF (any(lwrtnrm(:,ng))) THEN
404 CALL def_norm (ng, inlm, 1)
405 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
406
407 IF (nsa.eq.2) THEN
408 CALL def_norm (ng, inlm, 2)
409 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
410 END IF
411
412#ifdef ADJUST_BOUNDARY
413 CALL def_norm (ng, inlm, 3)
414 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
415#endif
416
417#if defined ADJUST_WSTRESS || defined ADJUST_STFLUX
418 CALL def_norm (ng, inlm, 4)
419 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
420#endif
421
422 DO tile=first_tile(ng),last_tile(ng),+1
423 CALL normalization (ng, tile, 2)
424 END DO
425 ldefnrm(1:4,ng)=.false.
426 lwrtnrm(1:4,ng)=.false.
427 ELSE
428 nrmrec=1
429 CALL get_state (ng, 14, 14, nrm(1,ng), nrmrec, 1)
430 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
431
432 IF (nsa.eq.2) THEN
433 CALL get_state (ng, 15, 15, nrm(2,ng), nrmrec, 2)
434 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
435 END IF
436
437#ifdef ADJUST_BOUNDARY
438 CALL get_state (ng, 16, 16, nrm(3,ng), nrmrec, 1)
439 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
440#endif
441
442#if defined ADJUST_WSTRESS || defined ADJUST_STFLUX
443 CALL get_state (ng, 17, 17, nrm(3,ng), nrmrec, 1)
444 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
445#endif
446 END IF
447 END DO
448!
449!------------------------------------------------------------------------
450! Compute the gradient or index, dJ/dS, of the sensitivity functional.
451!------------------------------------------------------------------------
452!
453 DO ng=1,ngrids
454 str_day=tdays(ng)
455 end_day=str_day-ntimes(ng)*dt(ng)*sec2day
456 IF ((dstrs(ng).eq.0.0_r8).and.(dends(ng).eq.0.0_r8)) THEN
457 dstrs(ng)=end_day
458 dends(ng)=str_day
459 END IF
460 IF (master) THEN
461 WRITE (stdout,10) 'AD', dends(ng), dstrs(ng)
462 END IF
463 IF ((dstrs(ng).gt.str_day).or.(dstrs(ng).lt.end_day)) THEN
464 IF (master) WRITE (stdout,20) 'DstrS = ', dstrs(ng), &
465 & end_day, str_day
466 exit_flag=7
467 RETURN
468 END IF
469 IF ((dends(ng).gt.str_day).or.(dends(ng).lt.end_day)) THEN
470 IF (master) WRITE (stdout,20) 'DendS = ', dends(ng), &
471 & end_day, str_day
472 exit_flag=7
473 RETURN
474 END IF
475 END DO
476!
477! Time-step adjoint model.
478!
479 DO ng=1,ngrids
480 IF (master) THEN
481 WRITE (stdout,30) 'AD', ng, ntstart(ng), ntend(ng)
482 END IF
483 END DO
484!
485#ifdef SOLVE3D
486 CALL ad_main3d (runinterval)
487#else
488 CALL ad_main2d (runinterval)
489#endif
490 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
491!
492! Convolve adjoint trajectory with error covariances and write
493! TLM initial conditions.
494!
495 lposterior=.false.
496 CALL error_covariance (itlm, driver, -1, -1, &
497 & lbck, lini, lold, lnew, &
498 & rec1, rec2, lposterior)
499 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
500!
501!-----------------------------------------------------------------------
502! Run tangent linear model for all nested grids, if any.
503!-----------------------------------------------------------------------
504!
505! Initialize tangent linear model.
506!
507 DO ng=1,ngrids
508 CALL tl_initial (ng)
509 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
510 END DO
511!
512! Activate tangent linear output.
513!
514 DO ng=1,ngrids
515 ldeftlm(ng)=.true.
516 lwrttlm(ng)=.true.
517 lcycletlm(ng)=.false.
518 END DO
519!
520! Time-step tangent linear model.
521!
522 DO ng=1,ngrids
523 IF (master) THEN
524 WRITE (stdout,30) 'TL', ng, ntstart(ng), ntend(ng)
525 END IF
526 END DO
527!
528#ifdef SOLVE3D
529 CALL tl_main3d (runinterval)
530#else
531 CALL tl_main2d (runinterval)
532#endif
533 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
534!
535 10 FORMAT (/,'AD ROMS: adjoint forcing time range: ', &
536 & f12.4,' - ',f12.4 ,/)
537 20 FORMAT (/,' Out of range adjoint forcing time, ',a,f12.4,/, &
538 & ' It must be between ',f12.4,' and ',f12.4)
539 30 FORMAT (/,1x,a,1x,'ROMS: started time-stepping:', &
540 & ' (Grid: ',i2.2,' TimeSteps: ',i0,' - ',i0,')',/)
541!
542 RETURN
543 END SUBROUTINE roms_run
544!
545 SUBROUTINE roms_finalize
546!
547!=======================================================================
548! !
549! This routine terminates ROMS nonlinear and adjoint models !
550! execution. !
551! !
552!=======================================================================
553!
554! Local variable declarations.
555!
556 integer :: Fcount, ng, thread
557!
558 character (len=*), parameter :: MyFile = &
559 & __FILE__//", ROMS_finalize"
560!
561!-----------------------------------------------------------------------
562! If blowing-up, save latest model state into RESTART NetCDF file.
563!-----------------------------------------------------------------------
564!
565! If cycling restart records, write solution into the next record.
566!
567 DO ng=1,ngrids
568 IF (lwrtrst(ng).and.(exit_flag.eq.1)) THEN
569 IF (master) WRITE (stdout,10)
570 10 FORMAT (/,' Blowing-up: Saving latest model state into ', &
571 & ' RESTART file',/)
572 fcount=rst(ng)%load
573 IF (lcyclerst(ng).and.(rst(ng)%Nrec(fcount).ge.2)) THEN
574 rst(ng)%Rindex=2
575 lcyclerst(ng)=.false.
576 END IF
579#ifdef DISTRIBUTE
580 CALL wrt_rst (ng, myrank)
581#else
582 CALL wrt_rst (ng, -1)
583#endif
584 END IF
585 END DO
586!
587!-----------------------------------------------------------------------
588! Stop model and time profiling clocks, report memory requirements, and
589! close output NetCDF files.
590!-----------------------------------------------------------------------
591!
592! Stop time clocks.
593!
594 IF (master) THEN
595 WRITE (stdout,20)
596 20 FORMAT (/,'Elapsed wall CPU time for each process (seconds):',/)
597 END IF
598!
599 DO ng=1,ngrids
600 DO thread=thread_range
601 CALL wclock_off (ng, iadm, 0, __line__, myfile)
602 END DO
603 END DO
604!
605! Report dynamic memory and automatic memory requirements.
606!
607 CALL memory
608!
609! Close IO files.
610!
611 DO ng=1,ngrids
612 CALL close_inp (ng, iadm)
613 END DO
614 CALL close_out
615!
616 RETURN
617 END SUBROUTINE roms_finalize
618
619 END MODULE roms_kernel_mod
subroutine ad_initial(ng)
Definition ad_initial.F:4
subroutine ad_main2d
Definition ad_main2d.F:586
subroutine ad_main3d(runinterval)
Definition ad_main3d.F:4
subroutine edit_multifile(task)
subroutine memory
Definition memory.F:3
subroutine, public close_out
Definition close_io.F:175
subroutine, public close_inp(ng, model)
Definition close_io.F:92
subroutine, public error_covariance(model, driver, outloop, innloop, rbck, rini, rold, rnew, rec1, rec2, lposterior)
Definition convolve.F:215
subroutine, public def_norm(ng, model, ifile)
Definition def_norm.F:53
subroutine, public get_state(ng, model, msg, s, inirec, tindex)
Definition get_state.F:90
subroutine, public inp_par(model)
Definition inp_par.F:56
subroutine, public roms_initialize_arrays
Definition mod_arrays.F:351
subroutine, public roms_allocate_arrays(allocate_vars)
Definition mod_arrays.F:114
logical, dimension(:), allocatable wrtzetaref
type(t_io), dimension(:,:), allocatable std
type(t_io), dimension(:,:), allocatable nrm
type(t_io), dimension(:), allocatable fwd
type(t_io), dimension(:), allocatable rst
type(t_io), dimension(:), allocatable ini
integer stdout
integer isfsur
subroutine, public initialize_parallel
integer numthreads
integer, dimension(:), allocatable first_tile
integer mythread
logical master
integer, dimension(:), allocatable last_tile
integer ocn_comm_world
integer, parameter inlm
Definition mod_param.F:662
integer, dimension(:), allocatable ntilex
Definition mod_param.F:685
integer, parameter iadm
Definition mod_param.F:665
integer ngrids
Definition mod_param.F:113
integer, dimension(:), allocatable ntilee
Definition mod_param.F:686
integer, parameter itlm
Definition mod_param.F:663
integer nsa
Definition mod_param.F:612
logical, dimension(:,:), allocatable lwrtnrm
integer, dimension(:), allocatable ntimes
logical, dimension(:), allocatable ldefitl
real(dp), dimension(:), allocatable dt
logical lstiffness
integer blowup
logical, dimension(:), allocatable balance
logical, dimension(:), allocatable lreadfrc
real(dp), dimension(:), allocatable tdays
real(r8), dimension(:), allocatable dends
logical, dimension(:), allocatable ldefadj
logical, dimension(:), allocatable lcycleadj
logical, dimension(:), allocatable lwrtadj
real(dp), parameter sec2day
integer, dimension(:), allocatable ntend
logical, dimension(:), allocatable lcycletlm
logical, dimension(:,:), allocatable ldefnrm
integer exit_flag
real(r8), dimension(:), allocatable dstrs
logical, dimension(:), allocatable lwrttlm
logical, dimension(:), allocatable lwrtrst
logical, dimension(:), allocatable ldeftlm
integer, dimension(:), allocatable ntstart
logical, dimension(:), allocatable lreadfwd
integer noerror
logical, dimension(:), allocatable lcyclerst
logical, dimension(:), allocatable lreadblk
integer, dimension(:), allocatable lold
integer, dimension(:), allocatable lnew
subroutine, public normalization(ng, tile, ifac)
subroutine, public roms_finalize
Definition ad_roms.h:283
subroutine, public roms_run(runinterval)
Definition ad_roms.h:239
subroutine, public roms_initialize(first, mpicomm)
Definition ad_roms.h:52
integer function, public stdout_unit(mymaster)
Definition stdout_mod.F:48
logical, save set_stdoutunit
Definition stdout_mod.F:41
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52
subroutine, public tl_def_ini(ng)
Definition tl_def_ini.F:43
subroutine, public wrt_rst(ng, tile)
Definition wrt_rst.F:63
subroutine, public biconj(ng, tile, model, lbck)
subroutine, public balance_ref(ng, tile, lbck)
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
subroutine tl_initial(ng)
Definition tl_initial.F:4
subroutine tl_main2d
Definition tl_main2d.F:429
subroutine tl_main3d(runinterval)
Definition tl_main3d.F:4