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

Functions/Subroutines

subroutine, public wrt_evolved (ng, tile, kout, nout)
 
subroutine, private wrt_evolved_nf90 (ng, tile, kout, nout, lbij, ubij, lbi, ubi, lbj, ubj)
 
subroutine, private wrt_evolved_pio (ng, tile, kout, nout, lbij, ubij, lbi, ubi, lbj, ubj)
 

Function/Subroutine Documentation

◆ wrt_evolved()

subroutine, public wrt_evolved_mod::wrt_evolved ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) kout,
integer, intent(in) nout )

Definition at line 54 of file wrt_evolved.F.

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

References mod_param::bounds, mod_scalars::exit_flag, strings_mod::founderror(), mod_ncparam::io_nf90, mod_ncparam::io_pio, mod_iounits::lze, mod_parallel::master, mod_scalars::noerror, mod_iounits::stdout, wrt_evolved_nf90(), and wrt_evolved_pio().

Referenced by i4dvar_mod::increment().

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

◆ wrt_evolved_nf90()

subroutine, private wrt_evolved_mod::wrt_evolved_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 113 of file wrt_evolved.F.

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

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

◆ wrt_evolved_pio()

subroutine, private wrt_evolved_mod::wrt_evolved_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 583 of file wrt_evolved.F.

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

Referenced by wrt_evolved().

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