ROMS
Loading...
Searching...
No Matches
obs_sen_rbl4dvar_forecast.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 Andrew M. Moore !
6! Licensed under a MIT/X style license !
7! See License_ROMS.md !
8!=======================================================================
9! !
10! ROMS Strong/Weak Constraint 4-Dimensional Variational Data !
11! Assimilation and Observation Sensitivity Driver: Restricted !
12! B-preconditioned Lanczos (RBL4D-Var). !
13! !
14! This driver is used for the dual formulation (observation space), !
15! strong or weak constraint 4D-Var where errors may be considered !
16! in both model and observations. !
17! !
18! It also computes the sensitivity of the assimilation system to !
19! each observation. It measures the degree to which each observation !
20! contributes to the uncertainty in the estimate. This analysis can !
21! be used to determine the type of measurements that need to be made, !
22! where to observe, and when. !
23! !
24! The routines in this driver control the initialization, time- !
25! stepping, and finalization of ROMS model following ESMF !
26! conventions: !
27! !
28! ROMS_initialize !
29! ROMS_run !
30! ROMS_finalize !
31! !
32! References: !
33! !
34! Moore, A.M., H.G. Arango, G. Broquet, B.S. Powell, A.T. Weaver, !
35! and J. Zavala-Garay, 2011: The Regional Ocean Modeling System !
36! (ROMS) 4-dimensional variational data assimilations systems, !
37! Part I - System overview and formulation, Prog. Oceanogr., 91, !
38! 34-49, doi:10.1016/j.pocean.2011.05.004. !
39! !
40! Moore, A.M., H.G. Arango, G. Broquet, C. Edward, M. Veneziani, !
41! B. Powell, D. Foley, J.D. Doyle, D. Costa, and P. Robinson, !
42! 2011: The Regional Ocean Modeling System (ROMS) 4-dimensional !
43! variational data assimilations systems, Part II - Performance !
44! and application to the California Current System, Prog. !
45! Oceanogr., 91, 50-73, doi:10.1016/j.pocean.2011.05.003. !
46! !
47! Moore, A.M., H.G. Arango, G. Broquet, C. Edward, M. Veneziani, !
48! B. Powell, D. Foley, J.D. Doyle, D. Costa, and P. Robinson, !
49! 2011: The Regional Ocean Modeling System (ROMS) 4-dimensional !
50! variational data assimilations systems, Part III - Observation !
51! impact and observation sensitivity in the California Current !
52! System, Prog. Oceanogr., 91, 74-94, !
53! doi:10.1016/j.pocean.2011.05.005. !
54! !
55!=======================================================================
56!
57 USE mod_param
58 USE mod_parallel
59 USE mod_arrays
60 USE mod_fourdvar
61 USE mod_iounits
62 USE mod_ncparam
63 USE mod_netcdf
64#if defined PIO_LIB && defined DISTRIBUTE
66#endif
67 USE mod_scalars
68 USE mod_stepping
69!
70#ifdef ADJUST_BOUNDARY
72#endif
74 USE mod_ocean, ONLY : initialize_ocean
75!
76 USE ad_wrt_his_mod, ONLY : ad_wrt_his
78# ifndef RPCG
79 USE congrad_mod, ONLY : congrad
80# endif
83 USE def_mod_mod, ONLY : def_mod
84 USE def_norm_mod, ONLY : def_norm
85 USE get_state_mod, ONLY : get_state
86 USE ini_adjust_mod, ONLY : load_adtotl
87 USE ini_adjust_mod, ONLY : load_tltoad
88 USE inp_par_mod, ONLY : inp_par
89#ifdef MCT_LIB
90# ifdef ATM_COUPLING
91 USE mct_coupler_mod, ONLY : initialize_ocn2atm_coupling
92# endif
93# ifdef WAV_COUPLING
94 USE mct_coupler_mod, ONLY : initialize_ocn2wav_coupling
95# endif
96#endif
98#ifdef RPCG
100#endif
103 USE strings_mod, ONLY : founderror, uppercase
104 USE tl_def_ini_mod, ONLY : tl_def_ini
105 USE wrt_ini_mod, ONLY : wrt_ini
106 USE wrt_rst_mod, ONLY : wrt_rst
107#if defined BALANCE_OPERATOR && defined ZETA_ELLIPTIC
109#endif
110!
111 implicit none
112!
113 PRIVATE
114 PUBLIC :: roms_initialize
115 PUBLIC :: roms_run
116 PUBLIC :: roms_finalize
117!
118 CONTAINS
119!
120 SUBROUTINE roms_initialize (first, mpiCOMM)
121!
122!=======================================================================
123! !
124! This routine allocates and initializes ROMS state variables !
125! and internal and external parameters. !
126! !
127!=======================================================================
128!
129! Imported variable declarations.
130!
131 logical, intent(inout) :: first
132!
133 integer, intent(in), optional :: mpiCOMM
134!
135! Local variable declarations.
136!
137 logical :: allocate_vars = .true.
138!
139#ifdef DISTRIBUTE
140 integer :: MyError, MySize
141#endif
142 integer :: STDrec, Tindex
143 integer :: chunk_size, ng, thread
144#ifdef _OPENMP
145 integer :: my_threadnum
146#endif
147!
148 character (len=*), parameter :: MyFile = &
149 & __FILE__//", ROMS_initialize"
150
151#ifdef DISTRIBUTE
152!
153!-----------------------------------------------------------------------
154! Set distribute-memory (mpi) world communictor.
155!-----------------------------------------------------------------------
156!
157 IF (PRESENT(mpicomm)) THEN
158 ocn_comm_world=mpicomm
159 ELSE
160 ocn_comm_world=mpi_comm_world
161 END IF
162 CALL mpi_comm_rank (ocn_comm_world, myrank, myerror)
163 CALL mpi_comm_size (ocn_comm_world, mysize, myerror)
164#endif
165!
166!-----------------------------------------------------------------------
167! On first pass, initialize model parameters a variables for all
168! nested/composed grids. Notice that the logical switch "first"
169! is used to allow multiple calls to this routine during ensemble
170! configurations.
171!-----------------------------------------------------------------------
172!
173 IF (first) THEN
174 first=.false.
175!
176! Initialize parallel control switches. These scalars switches are
177! independent from standard input parameters.
178!
180!
181! Set the ROMS standard output unit to write verbose execution info.
182! Notice that the default standard out unit in Fortran is 6.
183!
184! In some applications like coupling or disjointed mpi-communications,
185! it is advantageous to write standard output to a specific filename
186! instead of the default Fortran standard output unit 6. If that is
187! the case, it opens such formatted file for writing.
188!
189 IF (set_stdoutunit) THEN
191 set_stdoutunit=.false.
192 END IF
193!
194! Read in model tunable parameters from standard input. Allocate and
195! initialize variables in several modules after the number of nested
196! grids and dimension parameters are known.
197!
198 CALL inp_par (inlm)
199 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
200!
201! Set domain decomposition tile partition range. This range is
202! computed only once since the "first_tile" and "last_tile" values
203! are private for each parallel thread/node.
204!
205#if defined _OPENMP
206 mythread=my_threadnum()
207#elif defined DISTRIBUTE
209#else
210 mythread=0
211#endif
212 DO ng=1,ngrids
213 chunk_size=(ntilex(ng)*ntilee(ng)+numthreads-1)/numthreads
214 first_tile(ng)=mythread*chunk_size
215 last_tile(ng)=first_tile(ng)+chunk_size-1
216 END DO
217!
218! Initialize internal wall clocks. Notice that the timings does not
219! includes processing standard input because several parameters are
220! needed to allocate clock variables.
221!
222 IF (master) THEN
223 WRITE (stdout,10)
224 10 FORMAT (/,' Process Information:',/)
225 END IF
226!
227 DO ng=1,ngrids
228 DO thread=thread_range
229 CALL wclock_on (ng, inlm, 0, __line__, myfile)
230 END DO
231 END DO
232!
233! Allocate and initialize modules variables.
234!
235 CALL roms_allocate_arrays (allocate_vars)
237 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
238
239 END IF
240
241#if defined MCT_LIB && (defined ATM_COUPLING || defined WAV_COUPLING)
242!
243!-----------------------------------------------------------------------
244! Initialize coupling streams between model(s).
245!-----------------------------------------------------------------------
246!
247 DO ng=1,ngrids
248# ifdef ATM_COUPLING
249 CALL initialize_ocn2atm_coupling (ng, myrank)
250# endif
251# ifdef WAV_COUPLING
252 CALL initialize_ocn2wav_coupling (ng, myrank)
253# endif
254 END DO
255#endif
256
257#ifndef OBS_SPACE
258!
259!-----------------------------------------------------------------------
260! If the required vectors and arrays from congrad from a previous
261! run of the assimilation cycle are available, read them here from
262! LCZ(ng)%name NetCDF file.
263!-----------------------------------------------------------------------
264!
265 sourcefile=myfile
266 DO ng=1,ngrids
267 SELECT CASE (lcz(ng)%IOtype)
268 CASE (io_nf90)
269 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
270 & 'cg_beta', cg_beta, &
271 & start = (/1,1/), &
272 & total = (/ninner+1,nouter/))
273 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
274!
275 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
276 & 'cg_delta', cg_delta, &
277 & start = (/1,1/), &
278 & total = (/ninner,nouter/))
279 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
280!
281 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
282 & 'cg_Gnorm_v', cg_gnorm_v, &
283 & start = (/1/), &
284 & total = (/nouter/))
285 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
286!
287 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
288 & 'cg_dla', cg_dla, &
289 & start = (/1,1/), &
290 & total = (/ninner,nouter/))
291 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
292!
293 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
294 & 'cg_QG', cg_qg, &
295 & start = (/1,1/), &
296 & total = (/ninner+1,nouter/))
297 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
298!
299 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
300 & 'zgrad0', zgrad0, &
301 & start = (/1,1/), &
302 & total = (/ndatum(ng)+1,nouter/))
303 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
304!
305 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
306 & 'zcglwk', zcglwk, &
307 & start = (/1,1,1/), &
308 & total = (/ndatum(ng)+1,ninner+1, &
309 & nouter/))
310 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
311!
312 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
313 & 'TLmodVal_S', tlmodval_s, &
314 & start = (/1,1,1/), &
315 & total = (/ndatum(ng),ninner,nouter/), &
316 & broadcast = .false.) ! Master only
317 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
318
319# if defined RPCG && !defined OBS_IMPACT
320!
321 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
322 & 'Hbk', hbk, &
323 & start = (/1,1/), &
324 & total = (/ndatum(ng),nouter/))
325 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
326!
327 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
328 & 'Jb0', jb0, &
329 & start = (/1/), &
330 & total = (/nouter+1/))
331 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
332!
333 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
334 & 'vcglwk', vcglwk, &
335 & start = (/1,1,1/), &
336 & total = (/ndatum(ng)+1,ninner+1,
337 & nouter/))
338 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
339# endif
340
341# if defined PIO_LIB && defined DISTRIBUTE
342 CASE (io_pio)
343 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
344 & 'cg_beta', cg_beta, &
345 & start = (/1,1/), &
346 & total = (/ninner+1,nouter/))
347 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
348!
349 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
350 & 'cg_delta', cg_delta, &
351 & start = (/1,1/), &
352 & total = (/ninner,nouter/))
353 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
354!
355 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
356 & 'cg_Gnorm_v', cg_gnorm_v, &
357 & start = (/1/), &
358 & total = (/nouter/))
359 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
360!
361 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
362 & 'cg_dla', cg_dla, &
363 & start = (/1,1/), &
364 & total = (/ninner,nouter/))
365 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
366!
367 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
368 & 'cg_QG', cg_qg, &
369 & start = (/1,1/), &
370 & total = (/ninner+1,nouter/))
371 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
372!
373 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
374 & 'zgrad0', zgrad0, &
375 & start = (/1,1/), &
376 & total = (/ndatum(ng)+1,nouter/))
377 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
378!
379 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
380 & 'zcglwk', zcglwk, &
381 & start = (/1,1,1/), &
382 & total = (/ndatum(ng)+1,ninner+1, &
383 & nouter/))
384 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
385!
386 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
387 & 'TLmodVal_S', tlmodval_s, &
388 & start = (/1,1,1/), &
389 & total = (/ndatum(ng),ninner, &
390 & nouter/))
391 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
392
393# if defined RPCG && !defined OBS_IMPACT
394 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
395 & 'Hbk', hbk, &
396 & start = (/1,1/), &
397 & total = (/ndatum(ng),nouter/))
398 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
399!
400 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
401 & 'Jb0', jb0, &
402 & start = (/1/), &
403 & total = (/nouter+1/))
404 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
405!
406 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
407 & 'vcglwk', vcglwk, &
408 & start = (/1,1,1/), &
409 & total = (/ndatum(ng)+1,ninner+1, &
410 & nouter/))
411 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
412# endif
413# endif
414 END SELECT
415 END DO
416!
417!-----------------------------------------------------------------------
418! If skiping runing nonlinear model, read in observation screening and
419! quality control flag.
420!-----------------------------------------------------------------------
421!
422 sourcefile=myfile
423 wrtobsscale(1:ngrids)=.false.
424 DO ng=1,ngrids
425 SELECT CASE (lcz(ng)%IOtype)
426 CASE (io_nf90)
427 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
428 & vname(1,idobss), obsscale)
429 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
430
431# if defined PIO_LIB && defined DISTRIBUTE
432 CASE (io_pio)
433 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
434 & vname(1,idobss), obsscale)
435 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
436# endif
437 END SELECT
438 END DO
439#endif
440!
441!-----------------------------------------------------------------------
442! Read in standard deviation factors for error covariance.
443!-----------------------------------------------------------------------
444!
445! Initial conditions standard deviation. They are loaded in Tindex=1
446! of the e_var(...,Tindex) state variables.
447!
448 stdrec=1
449 tindex=1
450 DO ng=1,ngrids
451 CALL get_state (ng, 10, 10, std(1,ng), stdrec, tindex)
452 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
453 END DO
454!
455! Model error standard deviation. They are loaded in Tindex=2
456! of the e_var(...,Tindex) state variables.
457!
458 stdrec=1
459 tindex=2
460 IF (nsa.eq.2) THEN
461 DO ng=1,ngrids
462 CALL get_state (ng, 11, 11, std(2,ng), stdrec, tindex)
463 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
464 END DO
465 END IF
466
467#ifdef ADJUST_BOUNDARY
468!
469! Open boundary conditions standard deviation.
470!
471 stdrec=1
472 tindex=1
473 DO ng=1,ngrids
474 CALL get_state (ng, 12, 12, std(3,ng), stdrec, tindex)
475 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
476 END DO
477#endif
478#if defined ADJUST_WSTRESS || defined ADJUST_STFLUX
479!
480! Surface forcing standard deviation.
481!
482 stdrec=1
483 tindex=1
484 DO ng=1,ngrids
485 CALL get_state (ng, 13, 13, std(4,ng), stdrec, tindex)
486 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
487 END DO
488#endif
489!
490 RETURN
491 END SUBROUTINE roms_initialize
492!
493 SUBROUTINE roms_run (RunInterval)
494!
495!=======================================================================
496! !
497! This routine time-steps ROMS nonlinear, tangent linear and !
498! adjoint models. !
499! !
500!=======================================================================
501!
502! Imported variable declarations
503!
504 real(dp), intent(inout) :: runinterval ! seconds
505!
506! Local variable declarations.
507!
508 logical :: lcgini, linner, lposterior, add
509!
510 integer :: my_inner, my_outer
511 integer :: lbck, lini, rec1, rec2, impord
512 integer :: i, lstr, ng, status, tile
513 integer :: fcount, nrmrec
514
515 integer, dimension(Ngrids) :: indxsave
516 integer, dimension(Ngrids) :: nrec
517!
518 real(r8) :: str_day, end_day, dstarts, rtime
519!
520 character (len=1) :: chara, charb
521#ifdef OBS_SPACE
522 character (len=1) :: charc
523#endif
524 character (len=25) :: driver
525 character (len=20) :: string
526
527 character (len=*), parameter :: myfile = &
528 & __FILE__//", ROMS_run"
529!
530!=======================================================================
531! Run model for all nested grids, if any.
532!=======================================================================
533!
534! Initialize relevant parameters.
535!
536 DO ng=1,ngrids
537#if defined ADJUST_STFLUX || defined ADJUST_WSTRESS
538 lfinp(ng)=1 ! forcing index for input
539 lfout(ng)=1 ! forcing index for output history files
540#endif
541#ifdef ADJUST_BOUNDARY
542 lbinp(ng)=1 ! boundary index for input
543 lbout(ng)=1 ! boundary index for output history files
544#endif
545 lold(ng)=1 ! old minimization time index
546 lnew(ng)=2 ! new minimization time index
547 END DO
548 lini=1 ! NLM initial conditions record in INI
549 lbck=2 ! background record in INI
550 rec1=1
551 rec2=2
552 nrun=1
553 outer=0
554 inner=0
555 erstr=1
557 driver='obs_sen_rbl4dvar_forecast'
558 chara='A'
559 charb='B'
560#ifdef OBS_SPACE
561 charc='C'
562#endif
563!
564! Choose the desired order approximation of the forecast impact
565! calculations:
566!
567! ImpOrd: 1 (1st order)
568! 2 (2nd order) or
569! 3 (3rd order)
570!
571! ImpOrd=1 is the least accurate and can be in error by a factor of 2.
572!
573 impord=3
574!
575! For this driver DSTART is assumed to be the start time of the
576! analysis cycle that yields the analysis for the subsequent forecast
577! cycle. The number of time steps in the analysis cycle is NTIMES_ANA,
578! and the number of time steps in the forecast cycle is NTIMES_FCT.
579!
580! Initially we will be running the adjoint model over the forecast
581! time interval so we temporarily reset dstart and RunInterval to
582! reflect the forecast interval start and run time.
583!
585 rtime=0.0_r8
586 DO ng=1,ngrids
587 rtime=max(rtime, dt(ng)*ntimes_ana(ng))
588 END DO
589 dstart=dstart+rtime*sec2day
590!
591 rtime=0.0_r8
592 DO ng=1,ngrids
593 rtime=max(rtime, dt(ng)*ntimes_fct(ng))
594 ntimes(ng)=ntimes_fct(ng)
595 END DO
596 runinterval=rtime
597!
598!-----------------------------------------------------------------------
599! Configure weak constraint RBL4D-Var algorithm.
600!-----------------------------------------------------------------------
601!
602! Initialize the switch to gather weak constraint forcing.
603!
604 DO ng=1,ngrids
605 wrtforce(ng)=.false.
606 END DO
607!
608! Initialize and set nonlinear model initial conditions.
609!
610 DO ng=1,ngrids
611 lreadfwd(ng)=.false.
612 lreadblk(ng)=.false.
613 wrtnlmod(ng)=.true.
614 wrtrpmod(ng)=.false.
615 wrttlmod(ng)=.false.
616 lsen4dvar(ng)=.false.
617#ifdef OBS_SPACE
618 lobspace(ng)=.false.
619# ifndef OBS_IMPACT
620 ladjvar(ng)=.false.
621# endif
622#endif
623 END DO
624
625 CALL initial
626 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
627!
628! Save nonlinear initial conditions (currently in time index 1,
629! background) into record "Lbck" of INI(ng)%name NetCDF file. The
630! record "Lbck" becomes the background state record and the record
631! "Lini" becomes current nonlinear initial conditions.
632!
633 DO ng=1,ngrids
634 ini(ng)%Rindex=1
635 fcount=ini(ng)%load
636 ini(ng)%Nrec(fcount)=1
637#ifdef DISTRIBUTE
638 CALL wrt_ini (ng, myrank, 1)
639#else
640 CALL wrt_ini (ng, -1, 1)
641#endif
642 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
643 END DO
644!
645! The FWD structure is used several times and contains the nonnlinear
646! background trajectories (RBL4DVAR, FCTA, and FCTB) used to linearize
647! the tangent linear and adjoint models. These trajectories can be
648! split into multi-files. If so, the user needs to specified such
649! multiple files in the standard input script (roms.in). Since the HIS
650! structure is not used here, copy the original FWD containing the
651! trajectory information loaded during configuration in "inp_par",
652! which has the values for the regular RBL4D-Var at specified outer
653! loop into the HIS structure.
654!
655 CALL edit_multifile ('FWD2HIS')
656 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
657!
658!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
659! Model-error covariance normalization and stardard deviation factors.
660!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
661!
662! Compute or read in the error covariance normalization factors.
663! If computing, write out factors to NetCDF. This is an expensive
664! computation that needs to be computed only once for a particular
665! application grid and decorrelation scales.
666!
667 DO ng=1,ngrids
668 IF (any(lwrtnrm(:,ng))) THEN
669 CALL def_norm (ng, inlm, 1)
670 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
671
672 IF (nsa.eq.2) THEN
673 CALL def_norm (ng, inlm, 2)
674 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
675 END IF
676
677#ifdef ADJUST_BOUNDARY
678 CALL def_norm (ng, inlm, 3)
679 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
680#endif
681
682#if defined ADJUST_WSTRESS || defined ADJUST_STFLUX
683 CALL def_norm (ng, inlm, 4)
684 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
685#endif
686
687 DO tile=first_tile(ng),last_tile(ng),+1
688 CALL normalization (ng, tile, 2)
689 END DO
690
691 ldefnrm(1:4,ng)=.false.
692 lwrtnrm(1:4,ng)=.false.
693 ELSE
694 nrmrec=1
695 CALL get_state (ng, 14, 14, nrm(1,ng), nrmrec, 1)
696 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
697
698 IF (nsa.eq.2) THEN
699 CALL get_state (ng, 15, 15, nrm(2,ng), nrmrec, 2)
700 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
701 END IF
702
703#ifdef ADJUST_BOUNDARY
704 CALL get_state (ng, 16, 16, nrm(3,ng), nrmrec, 1)
705 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
706#endif
707
708#if defined ADJUST_WSTRESS || defined ADJUST_STFLUX
709 CALL get_state (ng, 17, 17, nrm(4,ng), nrmrec, 1)
710 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
711#endif
712 END IF
713 END DO
714
715#if defined BALANCE_OPERATOR && defined ZETA_ELLIPTIC
716!
717! Compute the reference zeta and biconjugate gradient arrays
718! required for the balance of free surface.
719!
720 IF (balance(isfsur)) THEN
721 DO ng=1,ngrids
722 CALL get_state (ng, inlm, 2, ini(ng), lini, lini)
723 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
724!
725 DO tile=first_tile(ng),last_tile(ng),+1
726 CALL balance_ref (ng, tile, lini)
727 CALL biconj (ng, tile, inlm, lini)
728 END DO
729 wrtzetaref(ng)=.true.
730 END DO
731 END IF
732#endif
733!
734! Define tangent linear initial conditions file.
735!
736 DO ng=1,ngrids
737 ldefitl(ng)=.true.
738 CALL tl_def_ini (ng)
739 ldefitl(ng)=.false.
740 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
741 END DO
742!
743! Define impulse forcing NetCDF file.
744!
745 DO ng=1,ngrids
746 ldeftlf(ng)=.true.
747 CALL def_impulse (ng)
748 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
749 END DO
750
751#ifndef OBS_SPACE
752!
753! Define output 4DVAR NetCDF file containing all processed data
754! at observation locations.
755!
756 DO ng=1,ngrids
757 ldefmod(ng)=.true.
758 CALL def_mod (ng)
759 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
760 END DO
761!
762! Write out outer loop beeing processed.
763!
764 sourcefile=myfile
765 DO ng=1,ngrids
766 SELECT CASE (dav(ng)%IOtype)
767 CASE (io_nf90)
768 CALL netcdf_put_ivar (ng, inlm, dav(ng)%name, &
769 & 'Nimpact', nimpact, &
770 & (/0/), (/0/), &
771 & ncid = dav(ng)%ncid)
772 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
773
774# if defined PIO_LIB && defined DISTRIBUTE
775 CASE (io_pio)
776 CALL pio_netcdf_put_ivar (ng, inlm, dav(ng)%name, &
777 & 'Nimpact', nimpact, &
778 & (/0/), (/0/), &
779 & piofile = dav(ng)%pioFile)
780 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
781# endif
782 END SELECT
783 END DO
784#endif
785!
786!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
787! Run the adjoint model forced by forecast minus analysis difference
788! using appropriate solution trajectories.
789!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
790!
791! Reset the start and end times for the adjoint forcing.
792!
793 DO ng=1,ngrids
794 str_day=tdays(ng)+ntimes(ng)*dt(ng)*sec2day
795 end_day=tdays(ng)
796 IF ((dstrs(ng).eq.0.0_r8).and.(dends(ng).eq.0.0_r8)) THEN
797 dstrs(ng)=end_day
798 dends(ng)=str_day
799 END IF
800 IF (master) THEN
801 WRITE (stdout,70) 'AD', dends(ng), dstrs(ng)
802 END IF
803 END DO
804!
805! Set basic state trajectory and adjoint forcing file.
806!
807! ads_xxxx_a.nc if obs_space is off
808! ads_xxxx_b.nc
809!
810! obs_xxxx_a.nc if obs_space is on
811! obs_xxxx_b.nc
812!
813 sourcefile=myfile
814 DO ng=1,ngrids
815#ifdef OBS_SPACE
816 CALL close_file (ng, inlm, obs(ng), obs(ng)%name)
817 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
818
819 IF (impord.ne.2) THEN
820 WRITE (obs(ng)%name,90) trim(oifb(ng)%head)
821 ELSE
822 WRITE (obs(ng)%name,90) trim(oifa(ng)%head)
823 END IF
824#else
825 CALL close_file (ng, inlm, ads(ng), ads(ng)%name)
826 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
827
828 IF (impord.ne.2) THEN
829 WRITE (ads(ng)%name,90) trim(foib(ng)%head)
830 ELSE
831 WRITE (ads(ng)%name,90) trim(foia(ng)%head)
832 END IF
833#endif
834 END DO
835!
836! Set structure for the nonlinear forward trajectory to be processed
837! by the tangent linear and adjoint models. Also, set switches to
838! process the FWD structure in routine "check_multifile". Notice that
839! it is possible to split solution into multiple NetCDF files to reduce
840! their size.
841!
842 CALL edit_multifile ('FCTB2FWD')
843 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
844 DO ng=1,ngrids
845 lreadfwd(ng)=.true.
846 END DO
847!
848! Set structure for the nonlinear surface fluxes to be processed by
849! by the tangent linear and adjoint models. Also, set switches to
850! process the BLK structure in routine "check_multifile". Notice that
851! it is possible to split solution into multiple NetCDF files to reduce
852! their size.
853!
854 CALL edit_multifile ('FCTB2BLK')
855 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
856 DO ng=1,ngrids
857 lreadblk(ng)=.true.
858 END DO
859!
860 IF (master) THEN
861 WRITE (stdout,50)
862 END IF
863!
864! Initialize the adjoint model: initialize using dI/dxf is appropriate.
865!
866! TODO: There is problem with the "haveObsMeta" flag. It is always
867! FALSE on the first call to OBS_READ, so "ObsMeta" is zero
868! for the radials. The reason this is happening is because that
869! NetCDF file had already been opened somewhere prior to calling
870! OBS_INITIAL in TL_INITIAL (and AD_INITIAL when appropriate),
871! so the obs arrays were not getting initialized properly by
872! these routines. For now, the solution is to close the NetCDF
873! file before each call to TL_INITIAL and AD_INITIAL.
874!
875 DO ng=1,ngrids
876 lstiffness=.false.
877#ifdef OBS_SPACE
878 lsen4dvar(ng)=.false.
879 lsenfct(ng)=.true.
880 lobspace(ng)=.true.
881# ifndef OBS_IMPACT
882 ladjvar(ng)=.false.
883# endif
884#else
885 lsen4dvar(ng)=.true.
886 lsenfct(ng)=.false.
887#endif
888 CALL close_file (ng, iadm, obs(ng), obs(ng)%name)
889 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
890!
891 CALL ad_initial (ng)
892 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
893 END DO
894!
895! Set adjoint history NetCDF parameters. Define adjoint history
896! file once to avoid opening to many files.
897!
898 DO ng=1,ngrids
899 wrtforce=.true.
900 IF (adm(ng)%ncid.ne.-1) ldefadj(ng)=.false.
901 fcount=adm(ng)%load
902 adm(ng)%Nrec(fcount)=0
903 adm(ng)%Rindex=0
904 END DO
905
906! Time-step adjoint model backwards.
907!
908 DO ng=1,ngrids
909 IF (master) THEN
910 WRITE (stdout,20) 'AD', ng, ntstart(ng), ntend(ng)
911 END IF
912 END DO
913!
914#ifdef SOLVE3D
915 CALL ad_main3d (runinterval)
916#else
917 CALL ad_main2d (runinterval)
918#endif
919 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
920!
921! Write out last weak-constraint forcing (WRTforce is still .TRUE.)
922! record into the adjoint history file. Note that the weak-constraint
923! forcing is delayed by nADJ time-steps.
924!
925 DO ng=1,ngrids
926#ifdef DISTRIBUTE
927 CALL ad_wrt_his (ng, myrank)
928#else
929 CALL ad_wrt_his (ng, -1)
930#endif
931 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
932 END DO
933!
934! Write out adjoint initial condition record into the adjoint
935! history file.
936!
937 DO ng=1,ngrids
938 wrtforce(ng)=.false.
939 lwrtstate2d(ng)=.false.
940#ifdef DISTRIBUTE
941 CALL ad_wrt_his (ng, myrank)
942#else
943 CALL ad_wrt_his (ng, -1)
944#endif
945 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
946 END DO
947!
948! Retrieve adjoint solution into tl_var(1).
949!
950 DO ng=1,ngrids
951 CALL get_state (ng, itlm, 4, adm(ng), adm(ng)%Rindex, rec1)
952 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
953 END DO
954!
955! Set basic state trajectory and adjoint forcing file.
956!
957 IF (impord.gt.1) THEN
958 sourcefile=myfile
959 DO ng=1,ngrids
960#ifdef OBS_SPACE
961 CALL close_file (ng, inlm, obs(ng), obs(ng)%name)
962 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
963
964 IF (impord.eq.2) THEN
965 WRITE (obs(ng)%name,90) trim(oifb(ng)%head)
966 ELSE
967 WRITE (obs(ng)%name,90) trim(oifa(ng)%head)
968 END IF
969#else
970 CALL close_file (ng, inlm, ads(ng), ads(ng)%name)
971 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
972
973 IF (impord.eq.2) THEN
974 WRITE (ads(ng)%name,90) trim(foib(ng)%head)
975 ELSE
976 WRITE (ads(ng)%name,90) trim(foia(ng)%head)
977 END IF
978#endif
979 END DO
980!
981! Set structure for the nonlinear forward trajectory to be processed
982! by the tangent linear and adjoint models. Also, set switches to
983! process the FWD structure in routine "check_multifile". Notice that
984! it is possible to split solution into multiple NetCDF files to reduce
985! their size.
986!
987 CALL edit_multifile ('FCTA2FWD')
988 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
989 DO ng=1,ngrids
990 lreadfwd(ng)=.true.
991 END DO
992!
993! Set structure for the nonlinear surface fluxes to be processed by
994! by the tangent linear and adjoint models. Also, set switches to
995! process the BLK structure in routine "check_multifile". Notice that
996! it is possible to split solution into multiple NetCDF files to reduce
997! their size.
998!
999 CALL edit_multifile ('FCTA2BLK')
1000 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1001 DO ng=1,ngrids
1002 lreadblk(ng)=.true.
1003 END DO
1004!
1005 IF (master) THEN
1006 WRITE (stdout,50)
1007 END IF
1008!
1009! Initialize the adjoint model: initialize using dI/dxf is appropriate.
1010!
1011 DO ng=1,ngrids
1012 lstiffness=.false.
1013#ifdef OBS_SPACE
1014 lsen4dvar(ng)=.false.
1015 lsenfct(ng)=.true.
1016 lobspace(ng)=.true.
1017# ifndef OBS_IMPACT
1018 ladjvar(ng)=.false.
1019# endif
1020#else
1021 lsen4dvar(ng)=.true.
1022 lsenfct(ng)=.false.
1023#endif
1024 CALL close_file (ng, iadm, obs(ng), obs(ng)%name)
1025 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1026!
1027 CALL ad_initial (ng)
1028 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1029 END DO
1030!
1031! Set adjoint history NetCDF parameters. Define adjoint history
1032! file once to avoid opening to many files.
1033!
1034 DO ng=1,ngrids
1035 wrtforce=.true.
1036 IF (adm(ng)%ncid.ne.-1) ldefadj(ng)=.false.
1037 fcount=adm(ng)%load
1038 adm(ng)%Nrec(fcount)=0
1039 adm(ng)%Rindex=0
1040 END DO
1041!
1042! Time-step adjoint model backwards.
1043!
1044 DO ng=1,ngrids
1045 IF (master) THEN
1046 WRITE (stdout,20) 'AD', ng, ntstart(ng), ntend(ng)
1047 END IF
1048 END DO
1049!
1050#ifdef SOLVE3D
1051 CALL ad_main3d (runinterval)
1052#else
1053 CALL ad_main2d (runinterval)
1054#endif
1055 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1056!
1057! Write out last weak-constraint forcing (WRTforce is still .TRUE.)
1058! record into the adjoint history file. Note that the weak-constraint
1059! forcing is delayed by nADJ time-steps.
1060!
1061 DO ng=1,ngrids
1062#ifdef DISTRIBUTE
1063 CALL ad_wrt_his (ng, myrank)
1064#else
1065 CALL ad_wrt_his (ng, -1)
1066#endif
1067 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1068 END DO
1069!
1070! Write out adjoint initial condition record into the adjoint
1071! history file.
1072!
1073 DO ng=1,ngrids
1074 wrtforce(ng)=.false.
1075 lwrtstate2d(ng)=.false.
1076#ifdef DISTRIBUTE
1077 CALL ad_wrt_his (ng, myrank)
1078#else
1079 CALL ad_wrt_his (ng, -1)
1080#endif
1081 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1082 END DO
1083!
1084! Retrieve adjoint solution.
1085!
1086 DO ng=1,ngrids
1087 CALL get_state (ng, iadm, 4, adm(ng), adm(ng)%Rindex, rec1)
1088 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1089 END DO
1090!
1091! Add the retrieved adjoint solution to the previous solution saved
1092! in tl_var(1).
1093!
1094 add=.true.
1095 DO ng=1,ngrids
1096 DO tile=first_tile(ng),last_tile(ng),+1
1097 CALL load_adtotl (ng, tile, rec1, rec1, add)
1098 END DO
1099 END DO
1100!
1101 END IF
1102!
1103!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1104!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1105! Adjoint of RBL4D-Var to compute the observation sensitivity.
1106!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1107!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1108!
1109! Reset the start and end times for the adjoint forcing.
1110!
1111 DO ng=1,ngrids
1112 str_day=tdays(ng)+ntimes(ng)*dt(ng)*sec2day
1113 end_day=tdays(ng)
1114 IF ((dstrs(ng).eq.0.0_r8).and.(dends(ng).eq.0.0_r8)) THEN
1115 dstrs(ng)=end_day
1116 dends(ng)=str_day
1117 END IF
1118 IF (master) THEN
1119 WRITE (stdout,70) 'AD', dends(ng), dstrs(ng)
1120 END IF
1121 END DO
1122!
1123! WARNING: ONLY one outer loop can be used for this application.
1124! ======= For more than 1 outer-loop, we require the second
1125! derivative of each model operator (i.e. the tangent linear
1126! of the tangent linear operator).
1127!
1128 ad_outer_loop : DO my_outer=nimpact,nimpact
1129 outer=my_outer
1130 inner=0
1131!
1132!-----------------------------------------------------------------------
1133! Run the adjoint model initialized and forced by dI/dx where I is the
1134! chosen function of the analysis/forecast state x.
1135!-----------------------------------------------------------------------
1136!
1137! Reset DSTART and RunIterval for the analysis cycle interval using
1138! NTIMES_ANA.
1139!
1141!
1142 rtime=0.0_r8
1143 DO ng=1,ngrids
1144 rtime=max(rtime, dt(ng)*ntimes_ana(ng))
1145 ntimes(ng)=ntimes_ana(ng)
1146 END DO
1147 runinterval=rtime
1148!
1149#if defined ADJUST_WSTRESS || defined ADJUST_STFLUX || \
1150 defined adjust_boundary
1151!
1152! Call initial again to reset the forcing and obc adjustment time
1153! counters.
1154!
1155 DO ng=1,ngrids
1156 lreadfwd(ng)=.false.
1157 lreadblk(ng)=.false.
1158 END DO
1159!
1160 CALL initial
1161#endif
1162!
1163! Set basic state trajectory.
1164!
1165 DO ng=1,ngrids
1166 WRITE (fwd(ng)%name,10) trim(his(ng)%head), outer-1
1167 END DO
1168!
1169! Set structure for the nonlinear forward trajectory (from regular
1170! RBL4DVAR) to be processed by the tangent linear and adjoint models.
1171! Also, set switches to process the FWD structure in routine
1172! "check_multifile". Notice that it is possible to split solution
1173! into multiple NetCDF files to reduce their size.
1174!
1175 CALL edit_multifile ('HIS2FWD')
1176 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1177 DO ng=1,ngrids
1178 lreadfwd(ng)=.true.
1179 END DO
1180
1181#ifdef FORWARD_FLUXES
1182!
1183! Set the BLK structure to contain the nonlinear model surface fluxes
1184! needed by the tangent linear and adjoint models. Also, set switches
1185! to process that structure in routine "check_multifile". Notice that
1186! it is possible to split the solution into multiple NetCDF files to
1187! reduce their size.
1188!
1189! The switch LreadFRC is deactivated because all the atmospheric
1190! forcing, including shortwave radiation, is read from the NLM
1191! surface fluxes or is assigned during ESM coupling. Such fluxes
1192! are available from the QCK structure. There is no need for reading
1193! and processing from the FRC structure input forcing-files.
1194!
1195 CALL edit_multifile ('QCK2BLK')
1196 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1197 DO ng=1,ngrids
1198 lreadblk(ng)=.true.
1199 lreadfrc(ng)=.false.
1200 lreadqck(ng)=.false.
1201 END DO
1202#endif
1203!
1204 IF (master) THEN
1205 WRITE (stdout,50)
1206 END IF
1207!
1208! Initialize the adjoint model: initialize using dI/dxf is
1209! appropriate.
1210!
1211 DO ng=1,ngrids
1212 lstiffness=.false.
1213!!AMM
1214#ifdef OBS_SPACE
1215 lsen4dvar(ng)=.true.
1216 lsenfct(ng)=.false.
1217 lobspace(ng)=.false.
1218# ifndef OBS_IMPACT
1219 ladjvar(ng)=.false.
1220# endif
1221#else
1222 lsen4dvar(ng)=.false.
1223 lsenfct(ng)=.false.
1224#endif
1225!!AMM
1226 CALL close_file (ng, iadm, obs(ng), obs(ng)%name)
1227 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1228!
1229 CALL ad_initial (ng)
1230 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1231 END DO
1232!
1233! Copy tl_var(1) solution to ad_var(1).
1234! It is important to clear the tl surface forcing and boundary
1235! arrays first since only the initial fields must be used to
1236! initialize the next run of the adjoint.
1237!
1238 add=.false.
1239 DO ng=1,ngrids
1240 DO tile=first_tile(ng),last_tile(ng),+1
1241#if defined ADJUST_WSTRESS || defined ADJUST_STFLUX
1242 CALL initialize_forces (ng, tile, itlm)
1243#endif
1244#ifdef ADJUST_BOUNDARY
1245 CALL initialize_boundary (ng, tile, itlm)
1246#endif
1247 CALL load_tltoad (ng, tile, rec1, rec1, add)
1248 END DO
1249 END DO
1250
1251#if defined ADJUST_WSTRESS || defined ADJUST_STFLUX || \
1252 defined adjust_boundary
1253!
1254! Clear the adjoint forcing and open boundary arrays
1255!
1256 DO ng=1,ngrids
1257 DO tile=first_tile(ng),last_tile(ng),+1
1258 CALL initialize_forces (ng, tile, iadm)
1259# ifdef ADJUST_BOUNDARY
1260 CALL initialize_boundary (ng, tile, iadm)
1261# endif
1262 END DO
1263 END DO
1264#endif
1265!
1266! Set adjoint history NetCDF parameters. Define adjoint history
1267! file one to avoid opening to many files.
1268!
1269 DO ng=1,ngrids
1270 wrtforce=.true.
1271 IF (adm(ng)%ncid.ne.-1) ldefadj(ng)=.false.
1272 fcount=adm(ng)%load
1273 adm(ng)%Nrec(fcount)=0
1274 adm(ng)%Rindex=0
1275 END DO
1276!
1277! Time-step adjoint model backwards.
1278! ??? What do we do in the case of model error? Save forcing for TLM?
1279!
1280 DO ng=1,ngrids
1281 IF (master) THEN
1282 WRITE (stdout,20) 'AD', ng, ntstart(ng), ntend(ng)
1283 END IF
1284 END DO
1285!
1286#ifdef SOLVE3D
1287 CALL ad_main3d (runinterval)
1288#else
1289 CALL ad_main2d (runinterval)
1290#endif
1291 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1292!
1293! Write out last weak-constraint forcing (WRTforce is still .TRUE.)
1294! record into the adjoint history file. Note that the weak-constraint
1295! forcing is delayed by nADJ time-steps.
1296!
1297 DO ng=1,ngrids
1298#ifdef DISTRIBUTE
1299 CALL ad_wrt_his (ng, myrank)
1300#else
1301 CALL ad_wrt_his (ng, -1)
1302#endif
1303 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1304 END DO
1305!
1306! Write out adjoint initial condition record into the adjoint
1307! history file.
1308!
1309 DO ng=1,ngrids
1310 wrtforce(ng)=.false.
1311 lwrtstate2d(ng)=.false.
1312#ifdef DISTRIBUTE
1313 CALL ad_wrt_his (ng, myrank)
1314#else
1315 CALL ad_wrt_his (ng, -1)
1316#endif
1317 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1318 END DO
1319!
1320! Convolve adjoint trajectory with error covariances.
1321!
1322 lposterior=.false.
1323 CALL error_covariance (itlm, driver, outer, inner, &
1324 & lbck, lini, lold, lnew, &
1325 & rec1, rec2, lposterior)
1326 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1327!
1328! Convert the current adjoint solution in ADM(ng)%name to impulse
1329! forcing. Write out impulse forcing into TLF(ng)%name NetCDF file.
1330! To facilitate the forcing to the TLM and RPM, the forcing is
1331! processed and written in increasing time coordinates (recall that
1332! the adjoint solution in ADM(ng)%name is backwards in time).
1333!
1334! AMM: Do not know what to do in the weak constraint case yet.
1335!
1336!! IF (Master) THEN
1337!! WRITE (stdout,40) outer, inner
1338!! END IF
1339!! DO ng=1,Ngrids
1340!! TLF(ng)%Rindex=0
1341#ifdef DISTRIBUTE
1342!! CALL wrt_impulse (ng, MyRank, iADM, ADM(ng)%name)
1343#else
1344!! CALL wrt_impulse (ng, -1, iADM, ADM(ng)%name)
1345#endif
1346!! IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN
1347!! END DO
1348!
1349!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
1350! Integrate tangent linear model forced by the convolved adjoint
1351! trajectory.
1352!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
1353!
1354#ifdef OBS_SPACE
1355!
1356! Reinitialize fourdvar arrays using the analysis interval observation
1357! file.
1358!
1359 DO ng=1,ngrids
1360 CALL close_file (ng, inlm, obs(ng), obs(ng)%name)
1361 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1362
1363 WRITE (obs(ng)%name,80) trim(obs(ng)%head), charc
1364 END DO
1365!
1368!
1369!-----------------------------------------------------------------------
1370! If the required vectors and arrays from congrad from a previous
1371! run of the assimilation cycle are available, read them here from
1372! LCZ(ng)%name NetCDF file.
1373!-----------------------------------------------------------------------
1374!
1375 sourcefile=myfile
1376 DO ng=1,ngrids
1377 SELECT CASE (lcz(ng)%IOtype)
1378 CASE (io_nf90)
1379 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1380 & 'cg_beta', cg_beta)
1382 & __line__, myfile)) RETURN
1383
1384 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1385 & 'cg_delta', cg_delta)
1387 & __line__, myfile)) RETURN
1388
1389 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1390 & 'cg_Gnorm_v', cg_gnorm_v)
1392 & __line__, myfile)) RETURN
1393
1394 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1395 & 'cg_dla', cg_dla)
1397 & __line__, myfile)) RETURN
1398
1399 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1400 & 'cg_QG', cg_qg)
1402 & __line__, myfile)) RETURN
1403
1404 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1405 & 'zgrad0', zgrad0)
1407 & __line__, myfile)) RETURN
1408
1409 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1410 & 'zcglwk', zcglwk)
1412 & __line__, myfile)) RETURN
1413
1414 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1415 & 'TLmodVal_S', tlmodval_s, &
1416 & broadcast = .false.) ! Master only
1418 & __line__, myfile)) RETURN
1419
1420# ifdef RPCG
1421 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1422 & 'Hbk', hbk)
1424 & __line__, myfile)) RETURN
1425
1426 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1427 & 'Jb0', jb0)
1429 & __line__, myfile)) RETURN
1430
1431 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1432 & 'vcglwk', vcglwk)
1434 & __line__, myfile)) RETURN
1435# endif
1436
1437# if defined PIO_LIB && defined DISTRIBUTE
1438 CASE (io_pio)
1439 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1440 & 'cg_beta', cg_beta)
1442 & __line__, myfile)) RETURN
1443
1444 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1445 & 'cg_delta', cg_delta)
1447 & __line__, myfile)) RETURN
1448
1449 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1450 & 'cg_Gnorm_v', cg_gnorm_v)
1452 & __line__, myfile)) RETURN
1453
1454 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1455 & 'cg_dla', cg_dla)
1457 & __line__, myfile)) RETURN
1458
1459 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1460 & 'cg_QG', cg_qg)
1462 & __line__, myfile)) RETURN
1463
1464 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1465 & 'zgrad0', zgrad0)
1467 & __line__, myfile)) RETURN
1468
1469 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1470 & 'zcglwk', zcglwk)
1472 & __line__, myfile)) RETURN
1473
1474 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1475 & 'TLmodVal_S', tlmodval_s)
1477 & __line__, myfile)) RETURN
1478
1479# ifdef RPCG
1480 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1481 & 'Hbk', hbk)
1483 & __line__, myfile)) RETURN
1484
1485 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1486 & 'Jb0', jb0)
1488 & __line__, myfile)) RETURN
1489
1490 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1491 & 'vcglwk', vcglwk)
1493 & __line__, myfile)) RETURN
1494# endif
1495# endif
1496 END SELECT
1497 END DO
1498!
1499!-----------------------------------------------------------------------
1500! If skiping runing nonlinear model, read in observation screening and
1501! quality control flag.
1502!-----------------------------------------------------------------------
1503!
1504 sourcefile=myfile
1505 wrtobsscale(1:ngrids)=.false.
1506 DO ng=1,ngrids
1507 SELECT CASE (lcz(ng)%IOtype)
1508 CASE (io_nf90)
1509 CALL netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1510 & vname(1,idobss), obsscale)
1511 IF (founderror(exit_flag, &
1512 & noerror, __line__, myfile)) RETURN
1513
1514# if defined PIO_LIB && defined DISTRIBUTE
1515 CASE (io_pio)
1516 CALL pio_netcdf_get_fvar (ng, itlm, lcz(ng)%name, &
1517 & vname(1,idobss), obsscale)
1518 IF (founderror(exit_flag, &
1519 & noerror, __line__, myfile)) RETURN
1520# endif
1521 END SELECT
1522 END DO
1523!
1524! Define output 4DVAR NetCDF file containing all processed data
1525! at observation locations.
1526!
1527 DO ng=1,ngrids
1528 ldefmod(ng)=.true.
1529 CALL def_mod (ng)
1530 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1531 END DO
1532#endif /* OBS_SPACE */
1533!
1534 DO ng=1,ngrids
1535 wrtnlmod(ng)=.false.
1536 wrttlmod(ng)=.true.
1537 lwrttlm(ng)=.false.
1538 END DO
1539!
1540! Clear tangent linear forcing arrays before entering inner-loop.
1541! This is very important.
1542!
1543 DO ng=1,ngrids
1544 DO tile=first_tile(ng),last_tile(ng),+1
1545 CALL initialize_forces (ng, tile, itlm)
1546#ifdef ADJUST_BOUNDARY
1547 CALL initialize_boundary (ng, tile, itlm)
1548#endif
1549 END DO
1550 END DO
1551!
1552! Set basic state trajectory.
1553!
1554 DO ng=1,ngrids
1555 WRITE (fwd(ng)%name,10) trim(his(ng)%head), outer-1
1556 lstr=len_trim(fwd(ng)%name)
1557 fwd(ng)%base=fwd(ng)%name(1:lstr-3)
1558 END DO
1559!
1560! If weak constraint, the impulses are time-interpolated at each
1561! time-steps.
1562!
1563 DO ng=1,ngrids
1564 IF (frcrec(ng).gt.3) THEN
1565 frequentimpulse(ng)=.true.
1566 END IF
1567 END DO
1568!
1569! Initialize tangent linear model from ITL(ng)%name, record Rec1.
1570!
1571 DO ng=1,ngrids
1572 itl(ng)%Rindex=rec1
1573 CALL close_file (ng, itlm, obs(ng), obs(ng)%name)
1574 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1575!
1576 CALL tl_initial (ng)
1577 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1578 END DO
1579!
1580! Run tangent linear model forward and force with convolved adjoint
1581! trajectory impulses. Compute (HMBM'H')_n * PSI at observation points
1582! which are used in the conjugate gradient algorithm.
1583!
1584 DO ng=1,ngrids
1585 IF (master) THEN
1586 WRITE (stdout,20) 'TL', ng, ntstart(ng), ntend(ng)
1587 END IF
1588 END DO
1589!
1590#ifdef SOLVE3D
1591 CALL tl_main3d (runinterval)
1592#else
1593 CALL tl_main2d (runinterval)
1594#endif
1595 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1596!
1597 DO ng=1,ngrids
1598 wrtnlmod(ng)=.false.
1599 wrttlmod(ng)=.false.
1600 END DO
1601
1602#ifdef OBS_IMPACT
1603!
1604! Compute observation impact to the data assimilation system.
1605!
1606 DO ng=1,ngrids
1607# ifdef RPCG
1608 CALL rep_matrix (ng, itlm, outer, ninner-1)
1609# else
1610 CALL rep_matrix (ng, itlm, outer, ninner)
1611# endif
1612 END DO
1613#else
1614!
1615! Set basic state trajectory for adjoint inner-loops.
1616!
1617 DO ng=1,ngrids
1618 WRITE (fwd(ng)%name,10) trim(his(ng)%head), outer-1
1619 lstr=len_trim(fwd(ng)%name)
1620 fwd(ng)%base=fwd(ng)%name(1:lstr-3)
1621 END DO
1622!
1623! Clear tangent linear forcing arrays before entering inner-loop.
1624!
1625 DO ng=1,ngrids
1626 DO tile=first_tile(ng),last_tile(ng),+1
1627 CALL initialize_forces (ng, tile, itlm)
1628# ifdef ADJUST_BOUNDARY
1629 CALL initialize_boundary (ng, tile, itlm)
1630# endif
1631 END DO
1632 END DO
1633!
1634# ifdef RPCG
1635 ad_inner_loop : DO my_inner=ninner,0,-1
1636# else
1637 ad_inner_loop : DO my_inner=ninner,1,-1
1638# endif
1639 inner=my_inner
1640
1641# ifdef RPCG
1642!
1643! Retrieve NLmodVal when inner=0 for use as BCKmodVal.
1644!
1645 IF (inner.eq.0) THEN
1646 DO ng=1,ngrids
1647 SELECT CASE (dav(ng)%IOtype)
1648 CASE (io_nf90)
1649 CALL netcdf_get_fvar (ng, itlm, dav(ng)%name, &
1650 & 'NLmodel_value', nlmodval)
1652 & __line__, myfile)) RETURN
1653
1654# if defined PIO_LIB && defined DISTRIBUTE
1655 CASE (io_pio)
1656 CALL pio_netcdf_get_fvar (ng, itlm, dav(ng)%name, &
1657 & 'NLmodel_value', nlmodval)
1659 & __line__, myfile)) RETURN
1660# endif
1661 END SELECT
1662 END DO
1663 END IF
1664 IF (inner.ne.ninner) THEN
1665 linner=.true.
1666 ELSE
1667 linner=.false.
1668 END IF
1669# endif
1670 IF (master) THEN
1671 WRITE (stdout,60) 'Adjoint of', uppercase('rbl4dvar'), &
1672 & outer, inner
1673 END IF
1674# ifdef RPCG
1675!
1676 inner_compute : IF (linner) THEN
1677!
1678# else
1679!
1680! Call adjoint conjugate gradient algorithm.
1681!
1682 lcgini=.false.
1683 DO ng=1,ngrids
1684 CALL ad_congrad (ng, itlm, outer, inner, ninner, lcgini)
1685 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1686 END DO
1687# endif
1688!
1689! Initialize the adjoint model from rest.
1690!
1691 DO ng=1,ngrids
1692 lsen4dvar(ng)=.false.
1693 lsenfct(ng)=.true.
1694# ifdef OBS_SPACE
1695 lobspace(ng)=.true.
1696# ifndef OBS_IMPACT
1697 ladjvar(ng)=.true.
1698# endif
1699# endif
1700 CALL close_file (ng, iadm, obs(ng), obs(ng)%name)
1701 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1702!
1703 CALL ad_initial (ng)
1704 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1705 wrtmisfit(ng)=.false.
1706 END DO
1707!
1708! Set adjoint history NetCDF parameters. Define adjoint history
1709! file only once to avoid opening too many files.
1710!
1711 DO ng=1,ngrids
1712 wrtforce(ng)=.true.
1713 IF (adm(ng)%ncid.ne.-1) ldefadj(ng)=.false.
1714 fcount=adm(ng)%load
1715 adm(ng)%Nrec(fcount)=0
1716 adm(ng)%Rindex=0
1717 END DO
1718!
1719! Time-step adjoint model backwards forced with current PSI vector.
1720!
1721 DO ng=1,ngrids
1722 IF (master) THEN
1723 WRITE (stdout,20) 'AD', ng, ntstart(ng), ntend(ng)
1724 END IF
1725 END DO
1726!
1727# ifdef SOLVE3D
1728 CALL ad_main3d (runinterval)
1729# else
1730 CALL ad_main2d (runinterval)
1731# endif
1732 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1733!
1734! Write out last weak-constraint forcing (WRTforce is still .TRUE.)
1735! record into the adjoint history file. Note that the weak-constraint
1736! forcing is delayed by nADJ time-steps.
1737!
1738 DO ng=1,ngrids
1739# ifdef DISTRIBUTE
1740 CALL ad_wrt_his (ng, myrank)
1741# else
1742 CALL ad_wrt_his (ng, -1)
1743# endif
1744 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1745 END DO
1746!
1747! Write out adjoint initial condition record into the adjoint
1748! history file.
1749!
1750 DO ng=1,ngrids
1751 wrtforce(ng)=.false.
1752 lwrtstate2d(ng)=.false.
1753# ifdef DISTRIBUTE
1754 CALL ad_wrt_his (ng, myrank)
1755# else
1756 CALL ad_wrt_his (ng, -1)
1757# endif
1758 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1759 END DO
1760!
1761! Convolve adjoint trajectory with error covariances.
1762!
1763 lposterior=.false.
1764 CALL error_covariance (itlm, driver, outer, inner, &
1765 & lbck, lini, lold, lnew, &
1766 & rec1, rec2, lposterior)
1767 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1768!
1769! Convert the current adjoint solution in ADM(ng)%name to impulse
1770! forcing. Write out impulse forcing into TLF(ng)%name NetCDF file.
1771! To facilitate the forcing to the TLM and RPM, the forcing is
1772! processed and written in increasing time coordinates (recall that
1773! the adjoint solution in ADM(ng)%name is backwards in time).
1774!
1775!!
1776!! AMM: Do not know what to do in the weak constraint case yet.
1777!!
1778!! IF (Master) THEN
1779!! WRITE (stdout,40) outer, inner
1780!! END IF
1781!! DO ng=1,Ngrids
1782!! TLF(ng)%Rindex=0
1783# ifdef DISTRIBUTE
1784!! CALL wrt_impulse (ng, MyRank, iADM, ADM(ng)%name)
1785# else
1786!! CALL wrt_impulse (ng, -1, iADM, ADM(ng)%name)
1787# endif
1788!! IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN
1789!! END DO
1790!
1791!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
1792! Integrate tangent linear model.
1793!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
1794!
1795! Initialize tangent linear model from initial impulse which is now
1796! stored in file ITL(ng)%name.
1797!
1798 DO ng=1,ngrids
1799 wrtnlmod(ng)=.false.
1800 wrttlmod(ng)=.true.
1801 END DO
1802!
1803! If weak constraint, the impulses are time-interpolated at each
1804! time-steps.
1805!
1806 DO ng=1,ngrids
1807 IF (frcrec(ng).gt.3) THEN
1808 frequentimpulse(ng)=.true.
1809 END IF
1810 END DO
1811!
1812! Initialize tangent linear model from ITL(ng)%name, record 1.
1813!
1814 DO ng=1,ngrids
1815 itl(ng)%Rindex=rec1
1816 CALL close_file (ng, itlm, obs(ng), obs(ng)%name)
1817!
1818 CALL tl_initial (ng)
1819 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1820 END DO
1821!
1822! Set tangent linear history NetCDF parameters. Define tangent linear
1823! history file at the beggining of each inner loop to avoid opening
1824! too many NetCDF files.
1825!
1826 DO ng=1,ngrids
1827 IF (inner.gt.ninner) ldeftlm(ng)=.false.
1828 fcount=tlm(ng)%load
1829 tlm(ng)%Nrec(fcount)=0
1830 tlm(ng)%Rindex=0
1831 END DO
1832!
1833! Run tangent linear model forward and force with convolved adjoint
1834! trajectory impulses.
1835!
1836 DO ng=1,ngrids
1837 IF (master) THEN
1838 WRITE (stdout,20) 'TL', ng, ntstart(ng), ntend(ng)
1839 END IF
1840 END DO
1841!
1842# ifdef SOLVE3D
1843 CALL tl_main3d (runinterval)
1844# else
1845 CALL tl_main2d (runinterval)
1846# endif
1847 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1848
1849 DO ng=1,ngrids
1850 wrtnlmod(ng)=.false.
1851 wrttlmod(ng)=.false.
1852 END DO
1853# ifdef RPCG
1854 END IF inner_compute
1855!
1856 DO ng=1,ngrids
1857 CALL ad_rpcg_lanczos (ng, irpm, outer, inner, ninner, &
1858 & lcgini)
1859 END DO
1860# endif
1861
1862 END DO ad_inner_loop
1863# ifndef RPCG
1864!
1865! Call adjoint conjugate gradient algorithm.
1866!
1867 inner=0
1868 lcgini=.true.
1869 DO ng=1,ngrids
1870 CALL ad_congrad (ng, itlm, outer, inner, ninner, lcgini)
1871 END DO
1872# endif
1873
1874#endif /* !OBS_IMPACT */
1875
1876#ifdef OBS_IMPACT
1877!
1878! Write out total observation impact.
1879!
1880 sourcefile=myfile
1881 DO ng=1,ngrids
1882 SELECT CASE (dav(ng)%IOtype)
1883 CASE (io_nf90)
1884 CALL netcdf_put_fvar (ng, inlm, dav(ng)%name, &
1885 & 'ObsImpact_total', ad_obsval, &
1886# ifdef IMPACT_INNER
1887 & (/1,1/), (/mobs,ninner/), &
1888# else
1889 & (/1/), (/mobs/), &
1890# endif
1891 & ncid = dav(ng)%ncid)
1893 & __line__, myfile)) RETURN
1894!
1895 CALL netcdf_sync (ng, inlm, dav(ng)%name, dav(ng)%ncid)
1897 & __line__, myfile)) RETURN
1898
1899# if defined PIO_LIB && defined DISTRIBUTE
1900 CASE (io_pio)
1901 CALL pio_netcdf_put_fvar (ng, inlm, dav(ng)%name, &
1902 & 'ObsImpact_total', ad_obsval, &
1903# ifdef IMPACT_INNER
1904 & (/1,1/), (/mobs,ninner/), &
1905# else
1906 & (/1/), (/mobs/), &
1907# endif
1908 & piofile = dav(ng)%pioFile)
1910 & __line__, myfile)) RETURN
1911!
1912 CALL pio_netcdf_sync (ng, inlm, dav(ng)%name, &
1913 & dav(ng)%pioFile)
1915 & __line__, myfile)) RETURN
1916# endif
1917 END SELECT
1918 END DO
1919#else
1920!
1921! Write out observation sensitivity.
1922!
1923 sourcefile=myfile
1924 DO ng=1,ngrids
1925 SELECT CASE (dav(ng)%IOtype)
1926 CASE (io_nf90)
1927 CALL netcdf_put_fvar (ng, itlm, dav(ng)%name, &
1928 & 'ObsSens_total', ad_obsval, &
1929# ifdef IMPACT_INNER
1930 & (/1,1/), (/mobs,ninner/), &
1931# else
1932 & (/1/), (/mobs/), &
1933# endif
1934 & ncid = dav(ng)%ncid)
1936 & __line__, myfile)) RETURN
1937!
1938 CALL netcdf_sync (ng, inlm, dav(ng)%name, dav(ng)%ncid)
1940 & __line__, myfile)) RETURN
1941
1942# if defined PIO_LIB && defined DISTRIBUTE
1943 CASE (io_pio)
1944 CALL pio_netcdf_put_fvar (ng, itlm, dav(ng)%name, &
1945 & 'ObsSens_total', ad_obsval, &
1946# ifdef IMPACT_INNER
1947 & (/1,1/), (/mobs,ninner/), &
1948# else
1949 & (/1/), (/mobs/), &
1950# endif
1951 & piofile = dav(ng)%pioFile)
1953 & __line__, myfile)) RETURN
1954!
1955 CALL pio_netcdf_sync (ng, inlm, dav(ng)%name, &
1956 & dav(ng)%pioFile)
1958 & __line__, myfile)) RETURN
1959# endif
1960 END SELECT
1961 END DO
1962#endif
1963!
1964! Close tangent linear NetCDF file.
1965!
1966 sourcefile=myfile
1967 DO ng=1,ngrids
1968 CALL close_file (ng, itlm, tlm(ng), tlm(ng)%name)
1969 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1970 END DO
1971
1972#if defined OBS_IMPACT && defined OBS_IMPACT_SPLIT
1973!
1974!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
1975! Integrate tangent linear model with initial condition increments
1976! only to compute the observation impact associated with the initial
1977! conditions.
1978!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
1979!
1980 DO ng=1,ngrids
1981 wrtnlmod(ng)=.false.
1982 wrttlmod(ng)=.true.
1983 lwrttlm(ng)=.false.
1984 END DO
1985!
1986! Clear tangent linear forcing arrays before entering inner-loop.
1987! This is very important since these arrays are non-zero after
1988! running the representer model and must be zero when running the
1989! tangent linear model.
1990!
1991 DO ng=1,ngrids
1992 DO tile=first_tile(ng),last_tile(ng),+1
1993 CALL initialize_forces (ng, tile, itlm)
1994# ifdef ADJUST_BOUNDARY
1995 CALL initialize_boundary (ng, tile, itlm)
1996# endif
1997 END DO
1998 END DO
1999!
2000! Set basic state trajectory.
2001!
2002 DO ng=1,ngrids
2003 WRITE (fwd(ng)%name,10) trim(his(ng)%head), outer-1
2004 lstr=len_trim(fwd(ng)%name)
2005 fwd(ng)%base=fwd(ng)%name(1:lstr-3)
2006 END DO
2007!
2008! If weak constraint, the impulses are time-interpolated at each
2009! time-steps.
2010!
2011 DO ng=1,ngrids
2012 IF (frcrec(ng).gt.3) THEN
2013 frequentimpulse(ng)=.true.
2014 END IF
2015 END DO
2016!
2017! Initialize tangent linear model from ITL(ng)%name, record 1.
2018!
2019 DO ng=1,ngrids
2020 itl(ng)%Rindex=rec1
2021 CALL close_file (ng, itlm, obs(ng), obs(ng)%name)
2022 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2023!
2024 CALL tl_initial (ng)
2025 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2026 END DO
2027!
2028! Clear tangent linear forcing arrays and boundary arrays
2029! before the obs impact initial condition calculation.
2030!
2031 DO ng=1,ngrids
2032 DO tile=first_tile(ng),last_tile(ng),+1
2033 CALL initialize_forces (ng, tile, itlm)
2034# ifdef ADJUST_BOUNDARY
2035 CALL initialize_boundary (ng, tile, itlm)
2036# endif
2037 END DO
2038 END DO
2039!
2040! Run tangent linear model forward and force with convolved adjoint
2041! trajectory impulses. Compute (HMBM'H')_n * PSI at observation points
2042! which are used in the conjugate gradient algorithm.
2043!
2044 DO ng=1,ngrids
2045 IF (master) THEN
2046 WRITE (stdout,20) 'TL', ng, ntstart(ng), ntend(ng)
2047 END IF
2048 END DO
2049!
2050# ifdef SOLVE3D
2051 CALL tl_main3d (runinterval)
2052# else
2053 CALL tl_main2d (runinterval)
2054# endif
2055 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2056!
2057 DO ng=1,ngrids
2058 wrtnlmod(ng)=.false.
2059 wrttlmod(ng)=.false.
2060 END DO
2061!
2062! Compute observation impact to the data assimilation system.
2063!
2064 DO ng=1,ngrids
2065# ifdef RPCG
2066 CALL rep_matrix (ng, itlm, outer, ninner-1)
2067# else
2068 CALL rep_matrix (ng, itlm, outer, ninner)
2069# endif
2070 END DO
2071!
2072! Write out observation sentivity.
2073!
2074 sourcefile=myfile
2075 DO ng=1,ngrids
2076 SELECT CASE (dav(ng)%IOtype)
2077 CASE (io_nf90)
2078 CALL netcdf_put_fvar (ng, itlm, dav(ng)%name, &
2079 & 'ObsImpact_IC', ad_obsval, &
2080# ifdef IMPACT_INNER
2081 & (/1,1/), (/mobs,ninner/), &
2082# else
2083 & (/1/), (/mobs/), &
2084# endif
2085 & ncid = dav(ng)%ncid)
2087 & __line__, myfile)) RETURN
2088!
2089 CALL netcdf_sync (ng, inlm, dav(ng)%name, dav(ng)%ncid)
2091 & __line__, myfile)) RETURN
2092
2093# if defined PIO_LIB && defined DISTRIBUTE
2094 CASE (io_pio)
2095 CALL pio_netcdf_put_fvar (ng, itlm, dav(ng)%name, &
2096 & 'ObsImpact_IC', ad_obsval, &
2097# ifdef IMPACT_INNER
2098 & (/1,1/), (/mobs,ninner/), &
2099# else
2100 & (/1/), (/mobs/), &
2101# endif
2102 & piofile = dav(ng)%pioFile)
2104 & __line__, myfile)) RETURN
2105!
2106 CALL pio_netcdf_sync (ng, inlm, dav(ng)%name, &
2107 & dav(ng)%pioFile)
2109 & __line__, myfile)) RETURN
2110# endif
2111 END SELECT
2112 END DO
2113
2114# if defined ADJUST_WSTRESS || defined ADJUST_STFLUX
2115!
2116!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
2117! Integrate tangent linear model with surface forcing increments
2118! only to compute the observation impact associated with the surface
2119! forcing.
2120!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
2121!
2122 DO ng=1,ngrids
2123 wrtnlmod(ng)=.false.
2124 wrttlmod(ng)=.true.
2125 lwrttlm(ng)=.false.
2126 END DO
2127!
2128! Clear tangent linear forcing arrays before entering inner-loop.
2129! This is very important since these arrays are non-zero after
2130! running the representer model and must be zero when running the
2131! tangent linear model.
2132!
2133 DO ng=1,ngrids
2134 DO tile=first_tile(ng),last_tile(ng),+1
2135 CALL initialize_forces (ng, tile, itlm)
2136# ifdef ADJUST_BOUNDARY
2137 CALL initialize_boundary (ng, tile, itlm)
2138# endif
2139 END DO
2140 END DO
2141!
2142! Set basic state trajectory.
2143!
2144 DO ng=1,ngrids
2145 WRITE (fwd(ng)%name,10) trim(his(ng)%head), outer-1
2146 lstr=len_trim(fwd(ng)%name)
2147 fwd(ng)%base=fwd(ng)%name(1:lstr-3)
2148 END DO
2149!
2150! If weak constraint, the impulses are time-interpolated at each
2151! time-steps.
2152!
2153 DO ng=1,ngrids
2154 IF (frcrec(ng).gt.3) THEN
2155 frequentimpulse(ng)=.true.
2156 END IF
2157 END DO
2158!
2159! Initialize tangent linear model from ITL(ng)%name, record 1.
2160!
2161 DO ng=1,ngrids
2162 itl(ng)%Rindex=rec1
2163 CALL close_file (ng, itlm, obs(ng), obs(ng)%name)
2164 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2165!
2166 CALL tl_initial (ng)
2167 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2168 END DO
2169!
2170! Clear tangent initial condition arrays and boundary arrays
2171! before the obs impact initial condition calculation.
2172!
2173 DO ng=1,ngrids
2174 DO tile=first_tile(ng),last_tile(ng),+1
2175 CALL initialize_ocean (ng, tile, itlm)
2176# ifdef ADJUST_BOUNDARY
2177 CALL initialize_boundary (ng, tile, itlm)
2178# endif
2179 END DO
2180 END DO
2181!
2182! Run tangent linear model forward and force with convolved adjoint
2183! trajectory impulses. Compute (HMBM'H')_n * PSI at observation points
2184! which are used in the conjugate gradient algorithm.
2185!
2186 DO ng=1,ngrids
2187 IF (master) THEN
2188 WRITE (stdout,20) 'TL', ng, ntstart(ng), ntend(ng)
2189 END IF
2190 END DO
2191!
2192# ifdef SOLVE3D
2193 CALL tl_main3d (runinterval)
2194# else
2195 CALL tl_main2d (runinterval)
2196# endif
2197 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2198!
2199 DO ng=1,ngrids
2200 wrtnlmod(ng)=.false.
2201 wrttlmod(ng)=.false.
2202 END DO
2203!
2204! Compute observation impact to the data assimilation system.
2205!
2206 DO ng=1,ngrids
2207# ifdef RPCG
2208 CALL rep_matrix (ng, itlm, outer, ninner-1)
2209# else
2210 CALL rep_matrix (ng, itlm, outer, ninner)
2211# endif
2212 END DO
2213!
2214! Write out observation sentivity.
2215!
2216 sourcefile=myfile
2217 DO ng=1,ngrids
2218 SELECT CASE (dav(ng)%IOtype)
2219 CASE (io_nf90)
2220 CALL netcdf_put_fvar (ng, itlm, dav(ng)%name, &
2221 & 'ObsImpact_FC', ad_obsval, &
2222# ifdef IMPACT_INNER
2223 & (/1,1/), (/mobs,ninner/), &
2224# else
2225 & (/1/), (/mobs/), &
2226# endif
2227 & ncid = dav(ng)%ncid)
2229 & __line__, myfile)) RETURN
2230!
2231 CALL netcdf_sync (ng, inlm, dav(ng)%name, dav(ng)%ncid)
2233 & __line__, myfile)) RETURN
2234
2235# if defined PIO_LIB && defined DISTRIBUTE
2236 CASE (io_pio)
2237 CALL pio_netcdf_put_fvar (ng, itlm, dav(ng)%name, &
2238 & 'ObsImpact_FC', ad_obsval, &
2239# ifdef IMPACT_INNER
2240 & (/1,1/), (/mobs,ninner/), &
2241# else
2242 & (/1/), (/mobs/), &
2243# endif
2244 & piofile = dav(ng)%pioFile)
2246 & __line__, myfile)) RETURN
2247!
2248 CALL pio_netcdf_sync (ng, inlm, dav(ng)%name, &
2249 & dav(ng)%pioFile)
2251 & __line__, myfile)) RETURN
2252# endif
2253 END SELECT
2254 END DO
2255# endif
2256
2257# if defined ADJUST_BOUNDARY
2258!
2259!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
2260! Integrate tangent linear model with boundary condition increments
2261! only to compute the observation impact associated with the boundary
2262! conditions.
2263!:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
2264!
2265 DO ng=1,ngrids
2266 wrtnlmod(ng)=.false.
2267 wrttlmod(ng)=.true.
2268 lwrttlm(ng)=.false.
2269 END DO
2270!
2271! Clear tangent linear forcing arrays before entering inner-loop.
2272! This is very important since these arrays are non-zero after
2273! running the representer model and must be zero when running the
2274! tangent linear model.
2275!
2276 DO ng=1,ngrids
2277 DO tile=first_tile(ng),last_tile(ng),+1
2278 CALL initialize_forces (ng, tile, itlm)
2279# ifdef ADJUST_BOUNDARY
2280 CALL initialize_boundary (ng, tile, itlm)
2281# endif
2282 END DO
2283 END DO
2284!
2285! Set basic state trajectory.
2286!
2287 DO ng=1,ngrids
2288 WRITE (fwd(ng)%name,10) trim(his(ng)%head), outer-1
2289 lstr=len_trim(fwd(ng)%name)
2290 fwd(ng)%base=fwd(ng)%name(1:lstr-3)
2291 END DO
2292!
2293! If weak constraint, the impulses are time-interpolated at each
2294! time-steps.
2295!
2296 DO ng=1,ngrids
2297 IF (frcrec(ng).gt.3) THEN
2298 frequentimpulse(ng)=.true.
2299 END IF
2300 END DO
2301!
2302! Initialize tangent linear model from ITL(ng)%name, record Rec1.
2303!
2304 DO ng=1,ngrids
2305 itl(ng)%Rindex=rec1
2306 CALL close_file (ng, itlm, obs(ng), obs(ng)%name)
2307 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2308!
2309 CALL tl_initial (ng)
2310 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2311 END DO
2312!
2313! Clear tangent linear initial condition and forcing arrays
2314! before the obs impact initial condition calculation.
2315!
2316 DO ng=1,ngrids
2317 DO tile=first_tile(ng),last_tile(ng),+1
2318 CALL initialize_ocean (ng, tile, itlm)
2319 CALL initialize_forces (ng, tile, itlm)
2320 END DO
2321 END DO
2322!
2323! Run tangent linear model forward and force with convolved adjoint
2324! trajectory impulses. Compute (HMBM'H')_n * PSI at observation points
2325! which are used in the conjugate gradient algorithm.
2326!
2327 DO ng=1,ngrids
2328 IF (master) THEN
2329 WRITE (stdout,20) 'TL', ng, ntstart(ng), ntend(ng)
2330 END IF
2331 END DO
2332!
2333# ifdef SOLVE3D
2334 CALL tl_main3d (runinterval)
2335# else
2336 CALL tl_main2d (runinterval)
2337# endif
2338 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2339!
2340 DO ng=1,ngrids
2341 wrtnlmod(ng)=.false.
2342 wrttlmod(ng)=.false.
2343 END DO
2344!
2345! Compute observation impact to the data assimilation system.
2346!
2347 DO ng=1,ngrids
2348# ifdef RPCG
2349 CALL rep_matrix (ng, itlm, outer, ninner-1)
2350# else
2351 CALL rep_matrix (ng, itlm, outer, ninner)
2352# endif
2353 END DO
2354!
2355! Write out observation sentivity.
2356!
2357 sourcefile=myfile
2358 DO ng=1,ngrids
2359 SELECT CASE (dav(ng)%IOtype)
2360 CASE (io_nf90)
2361 CALL netcdf_put_fvar (ng, itlm, dav(ng)%name, &
2362 & 'ObsImpact_BC', ad_obsval, &
2363# ifdef IMPACT_INNER
2364 & (/1,1/), (/mobs,ninner/), &
2365# else
2366 & (/1/), (/mobs/), &
2367# endif
2368 & ncid = dav(ng)%ncid)
2370 & __line__, myfile)) RETURN
2371!
2372 CALL netcdf_sync (ng, inlm, dav(ng)%name, dav(ng)%ncid)
2374 & __line__, myfile)) RETURN
2375
2376# if defined PIO_LIB && defined DISTRIBUTE
2377 CASE (io_pio)
2378 CALL pio_netcdf_put_fvar (ng, itlm, dav(ng)%name, &
2379 & 'ObsImpact_BC', ad_obsval, &
2380# ifdef IMPACT_INNER
2381 & (/1,1/), (/mobs,ninner/), &
2382# else
2383 & (/1/), (/mobs/), &
2384# endif
2385 & piofile = dav(ng)%pioFile)
2387 & __line__, myfile)) RETURN
2388!
2389 CALL pio_netcdf_sync (ng, inlm, dav(ng)%name, &
2390 & dav(ng)%pioFile)
2392 & __line__, myfile)) RETURN
2393# endif
2394 END SELECT
2395 END DO
2396# endif
2397#endif /* OBS_IMPACT_SPLIT */
2398!
2399! Close current forward NetCDF file.
2400!
2401 sourcefile=myfile
2402 DO ng=1,ngrids
2403 CALL close_file (ng, inlm, fwd(ng), fwd(ng)%name)
2404 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2405!
2406 IF (his(ng)%IOtype.eq.io_nf90) THEN
2407 his(ng)%ncid=-1
2408#if defined PIO_LIB && defined DISTRIBUTE
2409 ELSE IF (his(ng)%IOtype.eq.io_pio) THEN
2410 his(ng)%pioFile%fh=-1
2411#endif
2412 END IF
2413 END DO
2414
2415 END DO ad_outer_loop
2416!
2417 10 FORMAT (a,'_outer',i0,'.nc')
2418 20 FORMAT (/,1x,a,1x,'ROMS: started time-stepping:', &
2419 & ' (Grid: ',i2.2,' TimeSteps: ',i8.8,' - ',i8.8,')',/)
2420 30 FORMAT (' (',i3.3,',',i3.3,'): ',a,' data penalty, Jdata = ', &
2421 & 1p,e17.10,0p,t68,a)
2422 40 FORMAT (/,' Converting Convolved Adjoint Trajectory to', &
2423 & ' Impulses: Outer = ',i3.3,' Inner = ',i3.3,/)
2424 50 FORMAT (/,'ROMS: Started adjoint Sensitivity calculation', &
2425 & ' ...',/)
2426 60 FORMAT (/,'ROMS: ',a,1x,a,', Outer = ',i3.3, &
2427 & ' Inner = ',i3.3,/)
2428 70 FORMAT (/,1x,a,1x,'ROMS: adjoint forcing time range: ', &
2429 & f12.4,' - ',f12.4 ,/)
2430 80 FORMAT (a,'_',a,'.nc')
2431 90 FORMAT (a,'.nc')
2432!
2433 RETURN
2434 END SUBROUTINE roms_run
2435!
2436 SUBROUTINE roms_finalize
2437!
2438!=======================================================================
2439! !
2440! This routine terminates ROMS nonlinear, tangent linear, and !
2441! adjoint models execution. !
2442! !
2443!=======================================================================
2444!
2445! Local variable declarations.
2446!
2447 integer :: fcount, ng, thread
2448!
2449 character (len=*), parameter :: myfile = &
2450 & __FILE__//", ROMS_finalize"
2451!
2452!-----------------------------------------------------------------------
2453! Read and write observation variables for completeness.
2454!-----------------------------------------------------------------------
2455!
2456 DO ng=1,ngrids
2457#ifdef DISTRIBUTE
2458 CALL stats_modobs (ng, myrank)
2459#else
2460 CALL stats_modobs (ng, -1)
2461#endif
2462 END DO
2463!
2464!-----------------------------------------------------------------------
2465! If blowing-up, save latest model state into RESTART NetCDF file.
2466!-----------------------------------------------------------------------
2467!
2468! If cycling restart records, write solution into record 3.
2469!
2470 IF (exit_flag.eq.1) THEN
2471 DO ng=1,ngrids
2472 IF (lwrtrst(ng)) THEN
2473 IF (master) WRITE (stdout,10)
2474 10 FORMAT (/,' Blowing-up: Saving latest model state into ', &
2475 & ' RESTART file',/)
2476 fcount=rst(ng)%load
2477 IF (lcyclerst(ng).and.(rst(ng)%Nrec(fcount).ge.2)) THEN
2478 rst(ng)%Rindex=2
2479 lcyclerst(ng)=.false.
2480 END IF
2483#ifdef DISTRIBUTE
2484 CALL wrt_rst (ng, myrank)
2485#else
2486 CALL wrt_rst (ng, -1)
2487#endif
2488 END IF
2489 END DO
2490 END IF
2491!
2492!-----------------------------------------------------------------------
2493! Stop model and time profiling clocks. Close output NetCDF files.
2494!-----------------------------------------------------------------------
2495!
2496! Stop time clocks.
2497!
2498 IF (master) THEN
2499 WRITE (stdout,20)
2500 20 FORMAT (/,'Elapsed wall CPU time for each process (seconds):',/)
2501 END IF
2502!
2503 DO ng=1,ngrids
2504 DO thread=thread_range
2505 CALL wclock_off (ng, inlm, 0, __line__, myfile)
2506 END DO
2507 END DO
2508!
2509! Report dynamic memory and automatic memory requirements.
2510!
2511 CALL memory
2512!
2513! Close IO files.
2514!
2515 DO ng=1,ngrids
2516 CALL close_inp (ng, inlm)
2517 END DO
2518 CALL close_out
2519!
2520 RETURN
2521 END SUBROUTINE roms_finalize
2522
2523 END MODULE roms_kernel_mod
subroutine ad_congrad
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 initial
Definition initial.F:3
subroutine memory
Definition memory.F:3
subroutine, public ad_wrt_his(ng, tile)
Definition ad_wrt_his.F:66
subroutine, public close_out
Definition close_io.F:175
subroutine, public close_file(ng, model, s, ncname, lupdate)
Definition close_io.F:43
subroutine, public close_inp(ng, model)
Definition close_io.F:92
subroutine, public congrad(ng, model, outloop, innloop, ninnloop, lcgini)
Definition congrad.F:163
subroutine, public error_covariance(model, driver, outloop, innloop, rbck, rini, rold, rnew, rec1, rec2, lposterior)
Definition convolve.F:215
subroutine, public def_impulse(ng)
Definition def_impulse.F:49
subroutine, public def_mod(ng)
Definition def_mod.F:49
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 load_adtotl(ng, tile, linp, lout, add)
Definition ini_adjust.F:814
subroutine, public load_tltoad(ng, tile, linp, lout, add)
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
subroutine, public initialize_boundary(ng, tile, model)
subroutine, public initialize_forces(ng, tile, model)
real(dp), dimension(:), allocatable cg_gnorm_v
real(r8), dimension(:,:,:), allocatable vcglwk
real(dp), dimension(:,:), allocatable cg_beta
subroutine, public deallocate_fourdvar
integer, dimension(:), allocatable ndatum
real(r8), dimension(:,:), allocatable cg_dla
real(dp), dimension(:,:), allocatable cg_qg
logical, dimension(:), allocatable wrtmisfit
integer, dimension(:), allocatable ntimes_fct
logical, dimension(:), allocatable wrtrpmod
logical, dimension(:), allocatable lsenfct
logical, dimension(:), allocatable lsen4dvar
real(r8), dimension(:,:), allocatable ad_obsval
logical, dimension(:), allocatable wrttlmod
logical, dimension(:), allocatable wrtnlmod
logical, dimension(:), allocatable wrtforce
real(r8), dimension(:), allocatable obsscale
real(r8), dimension(:), allocatable jb0
integer nimpact
subroutine, public initialize_fourdvar
real(r8), dimension(:,:,:), allocatable zcglwk
logical, dimension(:), allocatable wrtobsscale
real(r8), dimension(:), allocatable nlmodval
real(dp), dimension(:,:), allocatable cg_delta
logical, dimension(:), allocatable lobspace
integer, dimension(:), allocatable ntimes_ana
logical, dimension(:), allocatable wrtzetaref
real(r8), dimension(:,:), allocatable zgrad0
type(t_io), dimension(:), allocatable ads
type(t_io), dimension(:), allocatable lcz
type(t_io), dimension(:,:), allocatable std
type(t_io), dimension(:), allocatable his
type(t_io), dimension(:), allocatable oifa
type(t_io), dimension(:,:), allocatable nrm
type(t_io), dimension(:), allocatable adm
type(t_io), dimension(:), allocatable obs
type(t_io), dimension(:), allocatable oifb
type(t_io), dimension(:), allocatable tlm
type(t_io), dimension(:), allocatable itl
type(t_io), dimension(:), allocatable dav
type(t_io), dimension(:), allocatable fwd
type(t_io), dimension(:), allocatable rst
type(t_io), dimension(:), allocatable ini
integer stdout
character(len=256) sourcefile
integer idobss
integer, parameter io_nf90
Definition mod_ncparam.F:95
integer, parameter io_pio
Definition mod_ncparam.F:96
integer isfsur
character(len=maxlen), dimension(6, 0:nv) vname
subroutine, public netcdf_sync(ng, model, ncname, ncid)
subroutine, public initialize_ocean(ng, tile, model)
Definition mod_ocean.F:1526
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, parameter irpm
Definition mod_param.F:664
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
subroutine, public pio_netcdf_sync(ng, model, ncname, piofile)
integer ninner
logical, dimension(:,:), allocatable lwrtnrm
logical, dimension(:), allocatable lreadqck
integer nouter
integer, dimension(:), allocatable ntimes
logical, dimension(:), allocatable ldefitl
real(dp), dimension(:), allocatable dt
logical lstiffness
integer blowup
logical, dimension(:), allocatable balance
logical, dimension(:), allocatable ldeftlf
integer erend
integer, dimension(:), allocatable frcrec
logical, dimension(:), allocatable lreadfrc
real(dp), dimension(:), allocatable tdays
real(r8), dimension(:), allocatable dends
logical, dimension(:), allocatable ldefadj
real(dp) dstart
logical, dimension(:), allocatable frequentimpulse
real(dp), parameter sec2day
integer, dimension(:), allocatable ntend
logical, dimension(:), allocatable ldefmod
logical, dimension(:,:), allocatable ldefnrm
integer exit_flag
logical, dimension(:), allocatable lwrtstate2d
integer erstr
real(r8), dimension(:), allocatable dstrs
logical, dimension(:), allocatable lwrttlm
logical, dimension(:), allocatable lwrtrst
real(dp) dstarts
logical, dimension(:), allocatable ldeftlm
integer, dimension(:), allocatable ntstart
integer nrun
logical, dimension(:), allocatable lreadfwd
integer inner
integer noerror
logical, dimension(:), allocatable lcyclerst
integer outer
logical, dimension(:), allocatable lreadblk
integer, dimension(:), allocatable lold
integer, dimension(:), allocatable lbout
integer, dimension(:), allocatable lfinp
integer, dimension(:), allocatable lbinp
integer, dimension(:), allocatable lnew
integer, dimension(:), allocatable lfout
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
subroutine, public rpcg_lanczos(ng, model, outloop, innloop, ninnloop, lcgini)
subroutine, public stats_modobs(ng, tile)
integer function, public stdout_unit(mymaster)
Definition stdout_mod.F:48
logical, save set_stdoutunit
Definition stdout_mod.F:41
character(len(sinp)) function, public uppercase(sinp)
Definition strings.F:582
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_ini(ng, tile, tindex, outrec)
Definition wrt_ini.F:89
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)
subroutine rep_matrix(ng, model, outloop, ninnloop)
Definition rep_matrix.F:4
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