ROMS
Loading...
Searching...
No Matches
wrt_hessian_mod Module Reference

Functions/Subroutines

subroutine, public wrt_hessian (ng, tile, kout, nout)
 
subroutine, private wrt_hessian_nf90 (ng, tile, kout, nout, lbij, ubij, lbi, ubi, lbj, ubj)
 
subroutine, private wrt_hessian_pio (ng, tile, kout, nout, lbij, ubij, lbi, ubi, lbj, ubj)
 

Function/Subroutine Documentation

◆ wrt_hessian()

subroutine, public wrt_hessian_mod::wrt_hessian ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) kout,
integer, intent(in) nout )

Definition at line 61 of file wrt_hessian.F.

62!***********************************************************************
63!
64! Imported variable declarations.
65!
66 integer, intent(in) :: ng, tile, kout, nout
67!
68! Local variable declarations.
69!
70# ifdef ADJUST_BOUNDARY
71 integer :: LBij, UBij
72# endif
73 integer :: LBi, UBi, LBj, UBj
74!
75 character (len=*), parameter :: MyFile = &
76 & __FILE__
77!
78!-----------------------------------------------------------------------
79! Write out time-averaged fields according to IO type.
80!-----------------------------------------------------------------------
81!
82# ifdef ADJUST_BOUNDARY
83 lbij=bounds(ng)%LBij
84 ubij=bounds(ng)%UBij
85# endif
86 lbi=bounds(ng)%LBi(tile)
87 ubi=bounds(ng)%UBi(tile)
88 lbj=bounds(ng)%LBj(tile)
89 ubj=bounds(ng)%UBj(tile)
90!
91 SELECT CASE (hss(ng)%IOtype)
92 CASE (io_nf90)
93 CALL wrt_hessian_nf90 (ng, tile, kout, nout, &
94# ifdef ADJUST_BOUNDARY
95 & lbij, ubij, &
96# endif
97 & lbi, ubi, lbj, ubj)
98
99# if defined PIO_LIB && defined DISTRIBUTE
100 CASE (io_pio)
101 CALL wrt_hessian_pio (ng, tile, kout, nout, &
102# ifdef ADJUST_BOUNDARY
103 & lbij, ubij, &
104# endif
105 & lbi, ubi, lbj, ubj)
106# endif
107 CASE DEFAULT
108 IF (master) WRITE (stdout,10) hss(ng)%IOtype
109 exit_flag=3
110 END SELECT
111 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
112!
113 10 FORMAT (' WRT_HESSIAN - Illegal output file type, io_type = ',i0, &
114 & /,15x,'Check KeyWord ''OUT_LIB'' in ''roms.in''.')
115!
116 RETURN

References mod_param::bounds, mod_scalars::exit_flag, strings_mod::founderror(), mod_iounits::hss, mod_ncparam::io_nf90, mod_ncparam::io_pio, mod_parallel::master, mod_scalars::noerror, mod_iounits::stdout, wrt_hessian_nf90(), and wrt_hessian_pio().

Referenced by convolve_mod::error_covariance(), cgradient_mod::hessian_evecs(), r4dvar_mod::increment(), rbl4dvar_mod::increment(), and posterior_mod::posterior_eofs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ wrt_hessian_nf90()

subroutine, private wrt_hessian_mod::wrt_hessian_nf90 ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) kout,
integer, intent(in) nout,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )
private

Definition at line 120 of file wrt_hessian.F.

125!***********************************************************************
126!
127 USE mod_netcdf
128!
129! Imported variable declarations.
130!
131 integer, intent(in) :: ng, tile, kout, nout
132# ifdef ADJUST_BOUNDARY
133 integer, intent(in) :: LBij, UBij
134# endif
135 integer, intent(in) :: LBi, UBi, LBj, UBj
136!
137! Local variable declarations.
138!
139 integer :: Fcount, i, j, gfactor, gtype, status
140# ifdef SOLVE3D
141 integer :: itrc, k
142# endif
143!
144 real(dp) :: scale
145!
146 character (len=*), parameter :: MyFile = &
147 & __FILE__//", wrt_hessian_nf90"
148!
149 sourcefile=myfile
150!
151!-----------------------------------------------------------------------
152! Write out Hessian eigenvectors.
153!-----------------------------------------------------------------------
154!
155 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
156!
157! Set grid type factor to write full (gfactor=1) fields or water
158! points (gfactor=-1) fields only.
159!
160# if defined WRITE_WATER && defined MASKING
161 gfactor=-1
162# else
163 gfactor=1
164# endif
165!
166! Set time record index.
167!
168 hss(ng)%Rindex=hss(ng)%Rindex+1
169 fcount=hss(ng)%Fcount
170 hss(ng)%Nrec(fcount)=hss(ng)%Nrec(fcount)+1
171!
172! Report.
173!
174# ifdef SOLVE3D
175# ifdef NESTING
176 IF (master) WRITE (stdout,10) kout, nout, hss(ng)%Rindex, ng
177# else
178 IF (master) WRITE (stdout,10) kout, nout, hss(ng)%Rindex
179# endif
180# else
181# ifdef NESTING
182 IF (master) WRITE (stdout,10) kout, hss(ng)%Rindex, ng
183# else
184 IF (master) WRITE (stdout,10) kout, hss(ng)%Rindex
185# endif
186# endif
187!
188! Write out model time (s).
189!
190 CALL netcdf_put_fvar (ng, iadm, hss(ng)%name, &
191 & trim(vname(1,idtime)), time(ng:), &
192 & (/hss(ng)%Rindex/), (/1/), &
193 & ncid = hss(ng)%ncid, &
194 & varid = hss(ng)%Vid(idtime))
195 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
196!
197! Write out free-surface.
198!
199 scale=1.0_dp
200 gtype=gfactor*r2dvar
201 status=nf_fwrite2d(ng, iadm, hss(ng)%ncid, idfsur, &
202 & hss(ng)%Vid(idfsur), &
203 & hss(ng)%Rindex, gtype, &
204 & lbi, ubi, lbj, ubj, scale, &
205# ifdef MASKING
206 & grid(ng) % rmask, &
207# endif
208 & ocean(ng)% ad_zeta(:,:,kout))
209 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
210 IF (master) THEN
211 WRITE (stdout,20) trim(vname(1,idfsur)), hss(ng)%Rindex
212 END IF
213 exit_flag=3
214 ioerror=status
215 RETURN
216 END IF
217
218# ifdef ADJUST_BOUNDARY
219!
220! Write out free-surface open boundaries.
221!
222 IF (any(lobc(:,isfsur,ng))) THEN
223 scale=1.0_dp
224 status=nf_fwrite2d_bry(ng, iadm, hss(ng)%name, hss(ng)%ncid, &
225 & vname(1,idsbry(isfsur)), &
226 & hss(ng)%Vid(idsbry(isfsur)), &
227 & hss(ng)%Rindex, r2dvar, &
228 & lbij, ubij, nbrec(ng), scale, &
229 & boundary(ng) % ad_zeta_obc(lbij:,:,:, &
230 & kout))
231 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
232 IF (master) THEN
233 WRITE (stdout,20) trim(vname(1,idsbry(isfsur))), &
234 & hss(ng)%Rindex
235 END IF
236 exit_flag=3
237 ioerror=status
238 RETURN
239 END IF
240 END IF
241# endif
242!
243! Write out 2D U-momentum component.
244!
245 scale=1.0_dp
246 gtype=gfactor*u2dvar
247 status=nf_fwrite2d(ng, iadm, hss(ng)%ncid, idubar, &
248 & hss(ng)%Vid(idubar), &
249 & hss(ng)%Rindex, gtype, &
250 & lbi, ubi, lbj, ubj, scale, &
251# ifdef MASKING
252 & grid(ng) % umask_full, &
253# endif
254 & ocean(ng) % ad_ubar(:,:,kout))
255 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
256 IF (master) THEN
257 WRITE (stdout,20) trim(vname(1,idubar)), hss(ng)%Rindex
258 END IF
259 exit_flag=3
260 ioerror=status
261 RETURN
262 END IF
263
264# ifdef ADJUST_BOUNDARY
265!
266! Write out 2D U-momentum component open boundaries.
267!
268 IF (any(lobc(:,isubar,ng))) THEN
269 scale=1.0_dp
270 status=nf_fwrite2d_bry(ng, iadm, hss(ng)%name, hss(ng)%ncid, &
271 & vname(1,idsbry(isubar)), &
272 & hss(ng)%Vid(idsbry(isubar)), &
273 & hss(ng)%Rindex, u2dvar, &
274 & lbij, ubij, nbrec(ng), scale, &
275 & boundary(ng) % ad_ubar_obc(lbij:,:,:, &
276 & kout))
277 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
278 IF (master) THEN
279 WRITE (stdout,20) trim(vname(1,idsbry(isubar))), &
280 & hss(ng)%Rindex
281 END IF
282 exit_flag=3
283 ioerror=status
284 RETURN
285 END IF
286 END IF
287# endif
288!
289! Write out 2D V-momentum component.
290!
291 scale=1.0_dp
292 gtype=gfactor*v2dvar
293 status=nf_fwrite2d(ng, iadm, hss(ng)%ncid, idvbar, &
294 & hss(ng)%Vid(idvbar), &
295 & hss(ng)%Rindex, gtype, &
296 & lbi, ubi, lbj, ubj, scale, &
297# ifdef MASKING
298 & grid(ng) % vmask_full, &
299# endif
300 & ocean(ng) % ad_vbar(:,:,kout))
301 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
302 IF (master) THEN
303 WRITE (stdout,20) trim(vname(1,idvbar)), hss(ng)%Rindex
304 END IF
305 exit_flag=3
306 ioerror=status
307 RETURN
308 END IF
309
310# ifdef ADJUST_BOUNDARY
311!
312! Write out 2D V-momentum component open boundaries.
313!
314 IF (any(lobc(:,isvbar,ng))) THEN
315 scale=1.0_dp
316 status=nf_fwrite2d_bry(ng, iadm, hss(ng)%name, hss(ng)%ncid, &
317 & vname(1,idsbry(isvbar)), &
318 & hss(ng)%Vid(idsbry(isvbar)), &
319 & hss(ng)%Rindex, v2dvar, &
320 & lbij, ubij, nbrec(ng), scale, &
321 & boundary(ng) % ad_vbar_obc(lbij:,:,:, &
322 & kout))
323 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
324 IF (master) THEN
325 WRITE (stdout,20) trim(vname(1,idsbry(isvbar))), &
326 & hss(ng)%Rindex
327 END IF
328 exit_flag=3
329 ioerror=status
330 RETURN
331 END IF
332 END IF
333# endif
334
335# ifdef SOLVE3D
336!
337! Write out 3D U-momentum component.
338!
339 scale=1.0_dp
340 gtype=gfactor*u3dvar
341 status=nf_fwrite3d(ng, iadm, hss(ng)%ncid, iduvel, &
342 & hss(ng)%Vid(iduvel), &
343 & hss(ng)%Rindex, gtype, &
344 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
345# ifdef MASKING
346 & grid(ng) % umask_full, &
347# endif
348 & ocean(ng) % ad_u(:,:,:,nout))
349 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
350 IF (master) THEN
351 WRITE (stdout,20) trim(vname(1,iduvel)), hss(ng)%Rindex
352 END IF
353 exit_flag=3
354 ioerror=status
355 RETURN
356 END IF
357
358# ifdef ADJUST_BOUNDARY
359!
360! Write out 3D U-momentum component open boundaries.
361!
362 IF (any(lobc(:,isuvel,ng))) THEN
363 scale=1.0_dp
364 status=nf_fwrite3d_bry(ng, iadm, hss(ng)%name, hss(ng)%ncid, &
365 & vname(1,idsbry(isuvel)), &
366 & hss(ng)%Vid(idsbry(isuvel)), &
367 & hss(ng)%Rindex, u3dvar, &
368 & lbij, ubij, 1, n(ng), nbrec(ng), scale, &
369 & boundary(ng) % ad_u_obc(lbij:,:,:,:, &
370 & nout))
371 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
372 IF (master) THEN
373 WRITE (stdout,20) trim(vname(1,idsbry(isuvel))), &
374 & hss(ng)%Rindex
375 END IF
376 exit_flag=3
377 ioerror=status
378 RETURN
379 END IF
380 END IF
381# endif
382!
383! Write out 3D V-momentum component.
384!
385 scale=1.0_dp
386 gtype=gfactor*v3dvar
387 status=nf_fwrite3d(ng, iadm, hss(ng)%ncid, idvvel, &
388 & hss(ng)%Vid(idvvel), &
389 & hss(ng)%Rindex, gtype, &
390 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
391# ifdef MASKING
392 & grid(ng) % vmask_full, &
393# endif
394 & ocean(ng) % ad_v(:,:,:,nout))
395 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
396 IF (master) THEN
397 WRITE (stdout,20) trim(vname(1,idvvel)), hss(ng)%Rindex
398 END IF
399 exit_flag=3
400 ioerror=status
401 RETURN
402 END IF
403
404# ifdef ADJUST_BOUNDARY
405!
406! Write out 3D V-momentum component open boundaries.
407!
408 IF (any(lobc(:,isvvel,ng))) THEN
409 scale=1.0_dp
410 status=nf_fwrite3d_bry(ng, iadm, hss(ng)%name, hss(ng)%ncid, &
411 & vname(1,idsbry(isvvel)), &
412 & hss(ng)%Vid(idsbry(isvvel)), &
413 & hss(ng)%Rindex, v3dvar, &
414 & lbij, ubij, 1, n(ng), nbrec(ng), scale, &
415 & boundary(ng) % ad_v_obc(lbij:,:,:,:, &
416 & nout))
417 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
418 IF (master) THEN
419 WRITE (stdout,20) trim(vname(1,idsbry(isvvel))), &
420 & hss(ng)%Rindex
421 END IF
422 exit_flag=3
423 ioerror=status
424 RETURN
425 END IF
426 END IF
427# endif
428!
429! Write out tracer type variables.
430!
431 DO itrc=1,nt(ng)
432 scale=1.0_dp
433 gtype=gfactor*r3dvar
434 status=nf_fwrite3d(ng, iadm, hss(ng)%ncid, idtvar(itrc), &
435 & hss(ng)%Tid(itrc), &
436 & hss(ng)%Rindex, gtype, &
437 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
438# ifdef MASKING
439 & grid(ng) % rmask, &
440# endif
441 & ocean(ng) % ad_t(:,:,:,nout,itrc))
442 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
443 IF (master) THEN
444 WRITE (stdout,20) trim(vname(1,idtvar(itrc))), &
445 & hss(ng)%Rindex
446 END IF
447 exit_flag=3
448 ioerror=status
449 RETURN
450 END IF
451 END DO
452
453# ifdef ADJUST_BOUNDARY
454!
455! Write out tracers open boundaries.
456!
457 DO itrc=1,nt(ng)
458 IF (any(lobc(:,istvar(itrc),ng))) THEN
459 scale=1.0_dp
460 status=nf_fwrite3d_bry(ng, iadm, hss(ng)%name, hss(ng)%ncid, &
461 & vname(1,idsbry(istvar(itrc))), &
462 & hss(ng)%Vid(idsbry(istvar(itrc))), &
463 & hss(ng)%Rindex, r3dvar, &
464 & lbij, ubij, 1, n(ng), nbrec(ng), &
465 & scale, &
466 & boundary(ng) % ad_t_obc(lbij:,:,:,:, &
467 & nout,itrc))
468 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
469 IF (master) THEN
470 WRITE (stdout,20) trim(vname(1,idsbry(istvar(itrc)))), &
471 & hss(ng)%Rindex
472 END IF
473 exit_flag=3
474 ioerror=status
475 RETURN
476 END IF
477 END IF
478 END DO
479# endif
480
481# ifdef ADJUST_STFLUX
482!
483! Write out surface net tracers fluxes. Notice that fluxes have their
484! own fixed time-dimension (of size Nfrec) to allow 4DVAR adjustments
485! at other times in addition to initialization time.
486!
487 DO itrc=1,nt(ng)
488 IF (lstflux(itrc,ng)) THEN
489 scale=1.0_dp
490 gtype=gfactor*r3dvar
491 status=nf_fwrite3d(ng, iadm, hss(ng)%ncid, idtsur(itrc), &
492 & hss(ng)%Vid(idtsur(itrc)), &
493 & hss(ng)%Rindex, gtype, &
494 & lbi, ubi, lbj, ubj, 1, nfrec(ng), scale, &
495# ifdef MASKING
496 & grid(ng) % rmask, &
497# endif
498 & forces(ng) % ad_tflux(:,:,:,kout,itrc))
499 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
500 IF (master) THEN
501 WRITE (stdout,20) trim(vname(1,idtsur(itrc))), &
502 & hss(ng)%Rindex
503 END IF
504 exit_flag=3
505 ioerror=status
506 RETURN
507 END IF
508 END IF
509 END DO
510# endif
511# endif
512# ifdef ADJUST_WSTRESS
513!
514! Write out surface U-momentum stress. Notice that the stress has its
515! own fixed time-dimension (of size Nfrec) to allow 4DVAR adjustments
516! at other times in addition to initialization time.
517!
518 scale=1.0_dp
519 gtype=gfactor*u3dvar
520 status=nf_fwrite3d(ng, iadm, hss(ng)%ncid, idusms, &
521 & hss(ng)%Vid(idusms), &
522 & hss(ng)%Rindex, gtype, &
523 & lbi, ubi, lbj, ubj, 1, nfrec(ng), scale, &
524# ifdef MASKING
525 & grid(ng) % umask, &
526# endif
527 & forces(ng) % ad_ustr(:,:,:,kout))
528 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
529 IF (master) THEN
530 WRITE (stdout,20) trim(vname(1,idusms)), hss(ng)%Rindex
531 END IF
532 exit_flag=3
533 ioerror=status
534 RETURN
535 END IF
536!
537! Write out surface V-momentum stress.
538!
539 scale=1.0_dp
540 gtype=gfactor*v3dvar
541 status=nf_fwrite3d(ng, iadm, hss(ng)%ncid, idvsms, &
542 & hss(ng)%Vid(idvsms), &
543 & hss(ng)%Rindex, gtype, &
544 & lbi, ubi, lbj, ubj, 1, nfrec(ng), scale, &
545# ifdef MASKING
546 & grid(ng) % vmask, &
547# endif
548 & forces(ng) % ad_vstr(:,:,:,kout))
549 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
550 IF (master) THEN
551 WRITE (stdout,20) trim(vname(1,idvsms)), hss(ng)%Rindex
552 END IF
553 exit_flag=3
554 ioerror=status
555 RETURN
556 END IF
557# endif
558!
559!-----------------------------------------------------------------------
560! Synchronize Hessian eigenvectors NetCDF file to disk to allow other
561! processes to access data immediately after it is written.
562!-----------------------------------------------------------------------
563!
564 CALL netcdf_sync (ng, iadm, hss(ng)%name, hss(ng)%ncid)
565 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
566!
567 10 FORMAT (2x,'WRT_HESSIAN_NF90 - writing Hessian', t42, &
568# ifdef SOLVE3D
569# ifdef NESTING
570 & 'fields (Index=',i1,',',i1,') in record = ',i0,t92,i2.2)
571# else
572 & 'fields (Index=',i1,',',i1,') in record = ',i0)
573# endif
574# else
575# ifdef NESTING
576 & 'fields (Index=',i1,') in record = ',i0,t92,i2.2)
577# else
578 & 'fields (Index=',i1,') in record = ',i0)
579# endif
580# endif
581 20 FORMAT (/,' WRT_HESSIAN_NF90 - error while writing variable: ',a, &
582 & /,20x,'into Hessian NetCDF file for time record: ',i0)
583!
584 RETURN
subroutine, public netcdf_sync(ng, model, ncname, ncid)

References mod_boundary::boundary, mod_scalars::exit_flag, mod_forces::forces, strings_mod::founderror(), mod_grid::grid, mod_iounits::hss, mod_param::iadm, mod_ncparam::idfsur, mod_ncparam::idsbry, mod_ncparam::idtime, mod_ncparam::idtsur, mod_ncparam::idtvar, mod_ncparam::idubar, mod_ncparam::idusms, mod_ncparam::iduvel, mod_ncparam::idvbar, mod_ncparam::idvsms, mod_ncparam::idvvel, mod_iounits::ioerror, mod_ncparam::isfsur, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::lobc, mod_scalars::lstflux, mod_parallel::master, mod_param::n, mod_scalars::nbrec, mod_netcdf::netcdf_sync(), mod_scalars::nfrec, mod_scalars::noerror, mod_param::nt, mod_ocean::ocean, mod_param::r2dvar, mod_param::r3dvar, mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::time, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, mod_param::v3dvar, and mod_ncparam::vname.

Referenced by wrt_hessian().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ wrt_hessian_pio()

subroutine, private wrt_hessian_mod::wrt_hessian_pio ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) kout,
integer, intent(in) nout,
integer, intent(in) lbij,
integer, intent(in) ubij,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )
private

Definition at line 590 of file wrt_hessian.F.

595!***********************************************************************
596!
598!
599! Imported variable declarations.
600!
601 integer, intent(in) :: ng, tile, kout, nout
602# ifdef ADJUST_BOUNDARY
603 integer, intent(in) :: LBij, UBij
604# endif
605 integer, intent(in) :: LBi, UBi, LBj, UBj
606!
607! Local variable declarations.
608!
609 integer :: Fcount, i, ifield, j, status
610# ifdef SOLVE3D
611 integer :: itrc, k
612# endif
613!
614 real(dp) :: scale
615!
616 character (len=*), parameter :: MyFile = &
617 & __FILE__//", wrt_hessian_pio"
618!
619 TYPE (IO_desc_t), pointer :: ioDesc
620!
621 sourcefile=myfile
622!
623!-----------------------------------------------------------------------
624! Write out Hessian eigenvectors.
625!-----------------------------------------------------------------------
626!
627 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
628!
629! Set time record index.
630!
631 hss(ng)%Rindex=hss(ng)%Rindex+1
632 fcount=hss(ng)%Fcount
633 hss(ng)%Nrec(fcount)=hss(ng)%Nrec(fcount)+1
634!
635! Report.
636!
637# ifdef SOLVE3D
638# ifdef NESTING
639 IF (master) WRITE (stdout,10) kout, nout, hss(ng)%Rindex, ng
640# else
641 IF (master) WRITE (stdout,10) kout, nout, hss(ng)%Rindex
642# endif
643# else
644# ifdef NESTING
645 IF (master) WRITE (stdout,10) kout, hss(ng)%Rindex, ng
646# else
647 IF (master) WRITE (stdout,10) kout, hss(ng)%Rindex
648# endif
649# endif
650!
651! Write out model time (s).
652!
653 CALL pio_netcdf_put_fvar (ng, iadm, hss(ng)%name, &
654 & trim(vname(1,idtime)), time(ng:), &
655 & (/hss(ng)%Rindex/), (/1/), &
656 & piofile = hss(ng)%pioFile, &
657 & piovar = hss(ng)%pioVar(idtime)%vd)
658 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
659!
660! Write out free-surface.
661!
662 scale=1.0_dp
663 IF (hss(ng)%pioVar(idfsur)%dkind.eq.pio_double) THEN
664 iodesc => iodesc_dp_r2dvar(ng)
665 ELSE
666 iodesc => iodesc_sp_r2dvar(ng)
667 END IF
668!
669 status=nf_fwrite2d(ng, iadm, hss(ng)%pioFile, idfsur, &
670 & hss(ng)%pioVar(idfsur), &
671 & hss(ng)%Rindex, &
672 & iodesc, &
673 & lbi, ubi, lbj, ubj, scale, &
674# ifdef MASKING
675 & grid(ng) % rmask, &
676# endif
677 & ocean(ng)% ad_zeta(:,:,kout))
678 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
679 IF (master) THEN
680 WRITE (stdout,20) trim(vname(1,idfsur)), hss(ng)%Rindex
681 END IF
682 exit_flag=3
683 ioerror=status
684 RETURN
685 END IF
686
687# ifdef ADJUST_BOUNDARY
688!
689! Write out free-surface open boundaries.
690!
691 IF (any(lobc(:,isfsur,ng))) THEN
692 scale=1.0_dp
693 IF (hss(ng)%pioVar(idsbry(isfsur))%dkind.eq.pio_double) THEN
694 iodesc => iodesc_dp_r2dobc(ng)
695 ELSE
696 iodesc => iodesc_sp_r2dobc(ng)
697 END IF
698!
699 status=nf_fwrite2d_bry(ng, iadm, hss(ng)%name, &
700 & hss(ng)%pioFile, &
701 & vname(1,idsbry(isfsur)), &
702 & hss(ng)%pioVar(idsbry(isfsur)), &
703 & hss(ng)%Rindex, &
704 & iodesc, &
705 & lbij, ubij, nbrec(ng), scale, &
706 & boundary(ng) % ad_zeta_obc(lbij:,:,:, &
707 & kout))
708 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
709 IF (master) THEN
710 WRITE (stdout,20) trim(vname(1,idsbry(isfsur))), &
711 & hss(ng)%Rindex
712 END IF
713 exit_flag=3
714 ioerror=status
715 RETURN
716 END IF
717 END IF
718# endif
719!
720! Write out 2D U-momentum component.
721!
722 scale=1.0_dp
723 IF (hss(ng)%pioVar(idubar)%dkind.eq.pio_double) THEN
724 iodesc => iodesc_dp_u2dvar(ng)
725 ELSE
726 iodesc => iodesc_sp_u2dvar(ng)
727 END IF
728!
729 status=nf_fwrite2d(ng, iadm, hss(ng)%pioFile, idubar, &
730 & hss(ng)%pioVar(idubar), &
731 & hss(ng)%Rindex, &
732 & iodesc, &
733 & lbi, ubi, lbj, ubj, scale, &
734# ifdef MASKING
735 & grid(ng) % umask_full, &
736# endif
737 & ocean(ng) % ad_ubar(:,:,kout))
738 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
739 IF (master) THEN
740 WRITE (stdout,20) trim(vname(1,idubar)), hss(ng)%Rindex
741 END IF
742 exit_flag=3
743 ioerror=status
744 RETURN
745 END IF
746
747# ifdef ADJUST_BOUNDARY
748!
749! Write out 2D U-momentum component open boundaries.
750!
751 IF (any(lobc(:,isubar,ng))) THEN
752 scale=1.0_dp
753 IF (hss(ng)%pioVar(idsbry(isubar))%dkind.eq.pio_double) THEN
754 iodesc => iodesc_dp_u2dobc(ng)
755 ELSE
756 iodesc => iodesc_sp_u2dobc(ng)
757 END IF
758!
759 status=nf_fwrite2d_bry(ng, iadm, hss(ng)%name, &
760 & hss(ng)%pioFile, &
761 & vname(1,idsbry(isubar)), &
762 & hss(ng)%pioVar(idsbry(isubar)), &
763 & hss(ng)%Rindex, &
764 & iodesc, &
765 & lbij, ubij, nbrec(ng), scale, &
766 & boundary(ng) % ad_ubar_obc(lbij:,:,:, &
767 & kout))
768 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
769 IF (master) THEN
770 WRITE (stdout,20) trim(vname(1,idsbry(isubar))), &
771 & hss(ng)%Rindex
772 END IF
773 exit_flag=3
774 ioerror=status
775 RETURN
776 END IF
777 END IF
778# endif
779!
780! Write out 2D V-momentum component.
781!
782 scale=1.0_dp
783 IF (hss(ng)%pioVar(idvbar)%dkind.eq.pio_double) THEN
784 iodesc => iodesc_dp_v2dvar(ng)
785 ELSE
786 iodesc => iodesc_sp_v2dvar(ng)
787 END IF
788!
789 status=nf_fwrite2d(ng, iadm, hss(ng)%pioFile, idvbar, &
790 & hss(ng)%pioVar(idvbar), &
791 & hss(ng)%Rindex, &
792 & iodesc, &
793 & lbi, ubi, lbj, ubj, scale, &
794# ifdef MASKING
795 & grid(ng) % vmask_full, &
796# endif
797 & ocean(ng) % ad_vbar(:,:,kout))
798 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
799 IF (master) THEN
800 WRITE (stdout,20) trim(vname(1,idvbar)), hss(ng)%Rindex
801 END IF
802 exit_flag=3
803 ioerror=status
804 RETURN
805 END IF
806
807# ifdef ADJUST_BOUNDARY
808!
809! Write out 2D V-momentum component open boundaries.
810!
811 IF (any(lobc(:,isvbar,ng))) THEN
812 scale=1.0_dp
813 IF (hss(ng)%pioVar(idsbry(isvbar))%dkind.eq.pio_double) THEN
814 iodesc => iodesc_dp_v2dobc(ng)
815 ELSE
816 iodesc => iodesc_sp_v2dobc(ng)
817 END IF
818!
819 status=nf_fwrite2d_bry(ng, iadm, hss(ng)%name, &
820 & hss(ng)%pioFile, &
821 & vname(1,idsbry(isvbar)), &
822 & hss(ng)%pioVar(idsbry(isvbar)), &
823 & hss(ng)%Rindex, &
824 & iodesc, &
825 & lbij, ubij, nbrec(ng), scale, &
826 & boundary(ng) % ad_vbar_obc(lbij:,:,:, &
827 & kout))
828 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
829 IF (master) THEN
830 WRITE (stdout,20) trim(vname(1,idsbry(isvbar))), &
831 & hss(ng)%Rindex
832 END IF
833 exit_flag=3
834 ioerror=status
835 RETURN
836 END IF
837 END IF
838# endif
839
840# ifdef SOLVE3D
841!
842! Write out 3D U-momentum component.
843!
844 scale=1.0_dp
845 IF (hss(ng)%pioVar(iduvel)%dkind.eq.pio_double) THEN
846 iodesc => iodesc_dp_u3dvar(ng)
847 ELSE
848 iodesc => iodesc_sp_u3dvar(ng)
849 END IF
850!
851 status=nf_fwrite3d(ng, iadm, hss(ng)%pioFile, iduvel, &
852 & hss(ng)%pioVar(iduvel), &
853 & hss(ng)%Rindex, &
854 & iodesc, &
855 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
856# ifdef MASKING
857 & grid(ng) % umask_full, &
858# endif
859 & ocean(ng) % ad_u(:,:,:,nout))
860 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
861 IF (master) THEN
862 WRITE (stdout,20) trim(vname(1,iduvel)), hss(ng)%Rindex
863 END IF
864 exit_flag=3
865 ioerror=status
866 RETURN
867 END IF
868
869# ifdef ADJUST_BOUNDARY
870!
871! Write out 3D U-momentum component open boundaries.
872!
873 IF (any(lobc(:,isuvel,ng))) THEN
874 scale=1.0_dp
875 IF (hss(ng)%pioVar(idsbry(isuvel))%dkind.eq.pio_double) THEN
876 iodesc => iodesc_dp_u3dobc(ng)
877 ELSE
878 iodesc => iodesc_sp_u3dobc(ng)
879 END IF
880!
881 status=nf_fwrite3d_bry(ng, iadm, hss(ng)%name, &
882 & hss(ng)%pioFile, &
883 & vname(1,idsbry(isuvel)), &
884 & hss(ng)%pioVar(idsbry(isuvel)), &
885 & hss(ng)%Rindex, &
886 & iodesc, &
887 & lbij, ubij, 1, n(ng), nbrec(ng), scale, &
888 & boundary(ng) % ad_u_obc(lbij:,:,:,:, &
889 & nout))
890 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
891 IF (master) THEN
892 WRITE (stdout,20) trim(vname(1,idsbry(isuvel))), &
893 & hss(ng)%Rindex
894 END IF
895 exit_flag=3
896 ioerror=status
897 RETURN
898 END IF
899 END IF
900# endif
901!
902! Write out 3D V-momentum component.
903!
904 scale=1.0_dp
905 IF (hss(ng)%pioVar(idvvel)%dkind.eq.pio_double) THEN
906 iodesc => iodesc_dp_v3dvar(ng)
907 ELSE
908 iodesc => iodesc_sp_v3dvar(ng)
909 END IF
910!
911 status=nf_fwrite3d(ng, iadm, hss(ng)%pioFile, idvvel, &
912 & hss(ng)%pioVar(idvvel), &
913 & hss(ng)%Rindex, &
914 & iodesc, &
915 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
916# ifdef MASKING
917 & grid(ng) % vmask_full, &
918# endif
919 & ocean(ng) % ad_v(:,:,:,nout))
920 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
921 IF (master) THEN
922 WRITE (stdout,20) trim(vname(1,idvvel)), hss(ng)%Rindex
923 END IF
924 exit_flag=3
925 ioerror=status
926 RETURN
927 END IF
928
929# ifdef ADJUST_BOUNDARY
930!
931! Write out 3D V-momentum component open boundaries.
932!
933 IF (any(lobc(:,isvvel,ng))) THEN
934 scale=1.0_dp
935 IF (hss(ng)%pioVar(idsbry(isvvel))%dkind.eq.pio_double) THEN
936 iodesc => iodesc_dp_v3dobc(ng)
937 ELSE
938 iodesc => iodesc_sp_v3dobc(ng)
939 END IF
940!
941 status=nf_fwrite3d_bry(ng, iadm, hss(ng)%name, &
942 & hss(ng)%pioFile, &
943 & vname(1,idsbry(isvvel)), &
944 & hss(ng)%pioVar(idsbry(isvvel)), &
945 & hss(ng)%Rindex, &
946 & iodesc, &
947 & lbij, ubij, 1, n(ng), nbrec(ng), scale, &
948 & boundary(ng) % ad_v_obc(lbij:,:,:,:, &
949 & nout))
950 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
951 IF (master) THEN
952 WRITE (stdout,20) trim(vname(1,idsbry(isvvel))), &
953 & hss(ng)%Rindex
954 END IF
955 exit_flag=3
956 ioerror=status
957 RETURN
958 END IF
959 END IF
960# endif
961!
962! Write out tracer type variables.
963!
964 DO itrc=1,nt(ng)
965 scale=1.0_dp
966 IF (hss(ng)%pioTrc(itrc)%dkind.eq.pio_double) THEN
967 iodesc => iodesc_dp_r3dvar(ng)
968 ELSE
969 iodesc => iodesc_sp_r3dvar(ng)
970 END IF
971!
972 status=nf_fwrite3d(ng, iadm, hss(ng)%pioFile, idtvar(itrc), &
973 & hss(ng)%pioTrc(itrc), &
974 & hss(ng)%Rindex, &
975 & iodesc, &
976 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
977# ifdef MASKING
978 & grid(ng) % rmask, &
979# endif
980 & ocean(ng) % ad_t(:,:,:,nout,itrc))
981 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
982 IF (master) THEN
983 WRITE (stdout,20) trim(vname(1,idtvar(itrc))), &
984 & hss(ng)%Rindex
985 END IF
986 exit_flag=3
987 ioerror=status
988 RETURN
989 END IF
990 END DO
991
992# ifdef ADJUST_BOUNDARY
993!
994! Write out tracers open boundaries.
995!
996 DO itrc=1,nt(ng)
997 IF (any(lobc(:,istvar(itrc),ng))) THEN
998 scale=1.0_dp
999 ifield=idsbry(istvar(itrc))
1000 IF (hss(ng)%pioVar(ifield)%dkind.eq.pio_double) THEN
1001 iodesc => iodesc_dp_r3dobc(ng)
1002 ELSE
1003 iodesc => iodesc_sp_r3dobc(ng)
1004 END IF
1005!
1006 status=nf_fwrite3d_bry(ng, iadm, hss(ng)%name, &
1007 & hss(ng)%pioFile, &
1008 & vname(1,ifield), &
1009 & hss(ng)%pioVar(ifield), &
1010 & hss(ng)%Rindex, &
1011 & iodesc, &
1012 & lbij, ubij, 1, n(ng), nbrec(ng), &
1013 & scale, &
1014 & boundary(ng) % ad_t_obc(lbij:,:,:,:, &
1015 & nout,itrc))
1016 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1017 IF (master) THEN
1018 WRITE (stdout,20) trim(vname(1,ifield)), hss(ng)%Rindex
1019 END IF
1020 exit_flag=3
1021 ioerror=status
1022 RETURN
1023 END IF
1024 END IF
1025 END DO
1026# endif
1027
1028# ifdef ADJUST_STFLUX
1029!
1030! Write out surface net tracers fluxes. Notice that fluxes have their
1031! own fixed time-dimension (of size Nfrec) to allow 4DVAR adjustments
1032! at other times in addition to initialization time.
1033!
1034 DO itrc=1,nt(ng)
1035 IF (lstflux(itrc,ng)) THEN
1036 scale=1.0_dp
1037 IF (hss(ng)%pioVar(idtsur(itrc))%dkind.eq.pio_double) THEN
1038 iodesc => iodesc_dp_r2dfrc(ng)
1039 ELSE
1040 iodesc => iodesc_sp_r2dfrc(ng)
1041 END IF
1042!
1043 status=nf_fwrite3d(ng, iadm, hss(ng)%pioFile, idtsur(itrc), &
1044 & hss(ng)%pioVar(idtsur(itrc)), &
1045 & hss(ng)%Rindex, &
1046 & iodesc, &
1047 & lbi, ubi, lbj, ubj, 1, nfrec(ng), scale, &
1048# ifdef MASKING
1049 & grid(ng) % rmask, &
1050# endif
1051 & forces(ng) % ad_tflux(:,:,:,kout,itrc))
1052 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1053 IF (master) THEN
1054 WRITE (stdout,20) trim(vname(1,idtsur(itrc))), &
1055 & hss(ng)%Rindex
1056 END IF
1057 exit_flag=3
1058 ioerror=status
1059 RETURN
1060 END IF
1061 END IF
1062 END DO
1063# endif
1064# endif
1065# ifdef ADJUST_WSTRESS
1066!
1067! Write out surface U-momentum stress. Notice that the stress has its
1068! own fixed time-dimension (of size Nfrec) to allow 4DVAR adjustments
1069! at other times in addition to initialization time.
1070!
1071 scale=1.0_dp
1072 IF (hss(ng)%pioVar(idusms)%dkind.eq.pio_double) THEN
1073 iodesc => iodesc_dp_u2dfrc(ng)
1074 ELSE
1075 iodesc => iodesc_sp_u2dfrc(ng)
1076 END IF
1077!
1078 status=nf_fwrite3d(ng, iadm, hss(ng)%pioFile, idusms, &
1079 & hss(ng)%pioVar(idusms), &
1080 & hss(ng)%Rindex, &
1081 & iodesc, &
1082 & lbi, ubi, lbj, ubj, 1, nfrec(ng), scale, &
1083# ifdef MASKING
1084 & grid(ng) % umask, &
1085# endif
1086 & forces(ng) % ad_ustr(:,:,:,kout))
1087 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1088 IF (master) THEN
1089 WRITE (stdout,20) trim(vname(1,idusms)), hss(ng)%Rindex
1090 END IF
1091 exit_flag=3
1092 ioerror=status
1093 RETURN
1094 END IF
1095!
1096! Write out surface V-momentum stress.
1097!
1098 scale=1.0_dp
1099 IF (hss(ng)%pioVar(idvsms)%dkind.eq.pio_double) THEN
1100 iodesc => iodesc_dp_v2dfrc(ng)
1101 ELSE
1102 iodesc => iodesc_sp_v2dfrc(ng)
1103 END IF
1104!
1105 status=nf_fwrite3d(ng, iadm, hss(ng)%pioFile, idvsms, &
1106 & hss(ng)%pioVar(idvsms), &
1107 & hss(ng)%Rindex, &
1108 & iodesc, &
1109 & lbi, ubi, lbj, ubj, 1, nfrec(ng), scale, &
1110# ifdef MASKING
1111 & grid(ng) % vmask, &
1112# endif
1113 & forces(ng) % ad_vstr(:,:,:,kout))
1114 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1115 IF (master) THEN
1116 WRITE (stdout,20) trim(vname(1,idvsms)), hss(ng)%Rindex
1117 END IF
1118 exit_flag=3
1119 ioerror=status
1120 RETURN
1121 END IF
1122# endif
1123!
1124!-----------------------------------------------------------------------
1125! Synchronize Hessian eigenvectors NetCDF file to disk to allow other
1126! processes to access data immediately after it is written.
1127!-----------------------------------------------------------------------
1128!
1129 CALL pio_netcdf_sync (ng, iadm, hss(ng)%name, hss(ng)%pioFile)
1130 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1131!
1132 10 FORMAT (2x,'WRT_HESSIAN_PIO - writing Hessian', t42, &
1133# ifdef SOLVE3D
1134# ifdef NESTING
1135 & 'fields (Index=',i1,',',i1,') in record = ',i0,t92,i2.2)
1136# else
1137 & 'fields (Index=',i1,',',i1,') in record = ',i0)
1138# endif
1139# else
1140# ifdef NESTING
1141 & 'fields (Index=',i1,') in record = ',i0,t92,i2.2)
1142# else
1143 & 'fields (Index=',i1,') in record = ',i0)
1144# endif
1145# endif
1146 20 FORMAT (/,' WRT_HESSIAN_PIO - error while writing variable: ',a, &
1147 & /,19x,'into Hessian NetCDF file for time record: ',i0)
1148!
1149 RETURN
type(io_desc_t), dimension(:), pointer iodesc_sp_v2dobc
type(io_desc_t), dimension(:), pointer iodesc_dp_r3dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_u3dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_v2dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_u2dfrc
type(io_desc_t), dimension(:), pointer iodesc_dp_v3dobc
subroutine, public pio_netcdf_sync(ng, model, ncname, piofile)
type(io_desc_t), dimension(:), pointer iodesc_sp_u2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_u3dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_r3dobc
type(io_desc_t), dimension(:), pointer iodesc_dp_u3dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_u2dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_u3dobc
type(io_desc_t), dimension(:), pointer iodesc_dp_u2dfrc
type(io_desc_t), dimension(:), pointer iodesc_dp_r2dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_u2dobc
type(io_desc_t), dimension(:), pointer iodesc_dp_u2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_v2dfrc
type(io_desc_t), dimension(:), pointer iodesc_dp_v3dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_v2dfrc
type(io_desc_t), dimension(:), pointer iodesc_dp_r2dfrc
type(io_desc_t), dimension(:), pointer iodesc_dp_v2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_v3dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_r3dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_r2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_r3dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_v3dobc
type(io_desc_t), dimension(:), pointer iodesc_sp_v2dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dfrc

References mod_boundary::boundary, mod_scalars::exit_flag, mod_forces::forces, strings_mod::founderror(), mod_grid::grid, mod_iounits::hss, mod_param::iadm, mod_ncparam::idfsur, mod_ncparam::idsbry, mod_ncparam::idtime, mod_ncparam::idtsur, mod_ncparam::idtvar, mod_ncparam::idubar, mod_ncparam::idusms, mod_ncparam::iduvel, mod_ncparam::idvbar, mod_ncparam::idvsms, mod_ncparam::idvvel, mod_pio_netcdf::iodesc_dp_r2dfrc, mod_pio_netcdf::iodesc_dp_r2dobc, mod_pio_netcdf::iodesc_dp_r2dvar, mod_pio_netcdf::iodesc_dp_r3dobc, mod_pio_netcdf::iodesc_dp_r3dvar, mod_pio_netcdf::iodesc_dp_u2dfrc, mod_pio_netcdf::iodesc_dp_u2dobc, mod_pio_netcdf::iodesc_dp_u2dvar, mod_pio_netcdf::iodesc_dp_u3dobc, mod_pio_netcdf::iodesc_dp_u3dvar, mod_pio_netcdf::iodesc_dp_v2dfrc, mod_pio_netcdf::iodesc_dp_v2dobc, mod_pio_netcdf::iodesc_dp_v2dvar, mod_pio_netcdf::iodesc_dp_v3dobc, mod_pio_netcdf::iodesc_dp_v3dvar, mod_pio_netcdf::iodesc_sp_r2dfrc, mod_pio_netcdf::iodesc_sp_r2dobc, mod_pio_netcdf::iodesc_sp_r2dvar, mod_pio_netcdf::iodesc_sp_r3dobc, mod_pio_netcdf::iodesc_sp_r3dvar, mod_pio_netcdf::iodesc_sp_u2dfrc, mod_pio_netcdf::iodesc_sp_u2dobc, mod_pio_netcdf::iodesc_sp_u2dvar, mod_pio_netcdf::iodesc_sp_u3dobc, mod_pio_netcdf::iodesc_sp_u3dvar, mod_pio_netcdf::iodesc_sp_v2dfrc, mod_pio_netcdf::iodesc_sp_v2dobc, mod_pio_netcdf::iodesc_sp_v2dvar, mod_pio_netcdf::iodesc_sp_v3dobc, mod_pio_netcdf::iodesc_sp_v3dvar, mod_iounits::ioerror, mod_ncparam::isfsur, mod_ncparam::istvar, mod_ncparam::isubar, mod_ncparam::isuvel, mod_ncparam::isvbar, mod_ncparam::isvvel, mod_scalars::lobc, mod_scalars::lstflux, mod_parallel::master, mod_param::n, mod_scalars::nbrec, mod_scalars::nfrec, mod_scalars::noerror, mod_param::nt, mod_ocean::ocean, mod_pio_netcdf::pio_netcdf_sync(), mod_iounits::sourcefile, mod_iounits::stdout, mod_scalars::time, and mod_ncparam::vname.

Referenced by wrt_hessian().

Here is the call graph for this function:
Here is the caller graph for this function: