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

Functions/Subroutines

subroutine, public wrt_dai (ng, tile)
 
subroutine, private wrt_dai_nf90 (ng, tile, lbi, ubi, lbj, ubj)
 
subroutine, private wrt_dai_pio (ng, tile, lbi, ubi, lbj, ubj)
 

Function/Subroutine Documentation

◆ wrt_dai()

subroutine, public wrt_dai_mod::wrt_dai ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 45 of file wrt_dai.F.

46!***********************************************************************
47!
48! Imported variable declarations.
49!
50 integer, intent(in) :: ng, tile
51!
52! Local variable declarations.
53!
54 integer :: LBi, UBi, LBj, UBj
55!
56 character (len=*), parameter :: MyFile = &
57 & __FILE__
58!
59!-----------------------------------------------------------------------
60! Write out time-averaged fields according to IO type.
61!-----------------------------------------------------------------------
62!
63 lbi=bounds(ng)%LBi(tile)
64 ubi=bounds(ng)%UBi(tile)
65 lbj=bounds(ng)%LBj(tile)
66 ubj=bounds(ng)%UBj(tile)
67!
68 SELECT CASE (dia(ng)%IOtype)
69 CASE (io_nf90)
70 CALL wrt_dai_nf90 (ng, tile, &
71 & lbi, ubi, lbj, ubj)
72
73# if defined PIO_LIB && defined DISTRIBUTE
74 CASE (io_pio)
75 CALL wrt_dai_pio (ng, tile, &
76 & lbi, ubi, lbj, ubj)
77# endif
78 CASE DEFAULT
79 IF (master) WRITE (stdout,10) dai(ng)%IOtype
80 exit_flag=3
81 END SELECT
82 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
83!
84 10 FORMAT (' WRT_DAI - Illegal output file type, io_type = ',i0, &
85 & /,11x,'Check KeyWord ''OUT_LIB'' in ''roms.in''.')
86!
87 RETURN

References mod_param::bounds, mod_iounits::dai, mod_iounits::dia, mod_scalars::exit_flag, strings_mod::founderror(), mod_ncparam::io_nf90, mod_ncparam::io_pio, mod_parallel::master, mod_scalars::noerror, mod_iounits::stdout, wrt_dai_nf90(), and wrt_dai_pio().

Here is the call graph for this function:

◆ wrt_dai_nf90()

subroutine, private wrt_dai_mod::wrt_dai_nf90 ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )
private

Definition at line 91 of file wrt_dai.F.

93!***********************************************************************
94!
95 USE mod_netcdf
96!
97! Imported variable declarations.
98!
99 integer, intent(in) :: ng, tile
100 integer, intent(in) :: LBi, UBi, LBj, UBj
101!
102! Local variable declarations.
103!
104 integer :: i, j, k, itrc
105 integer :: Fcount, gfactor, gtype, status, varid
106!
107 real(dp) :: scale
108!
109 character (len=*), parameter :: MyFile = &
110 & __FILE__//", wrt_dai_nf90"
111
112# include "set_bounds.h"
113!
114 sourcefile=myfile
115!
116!-----------------------------------------------------------------------
117! Write out Data Assimilation initial/restart fields.
118!-----------------------------------------------------------------------
119!
120 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
121!
122! Set grid type factor to write full (gfactor=1) fields.
123!
124 gfactor=1
125!
126! Set time record index.
127!
128 dai(ng)%Rindex=dai(ng)%Rindex+1
129 fcount=dai(ng)%Fcount
130 dai(ng)%Nrec(fcount)=dai(ng)%Nrec(fcount)+1
131!
132! If requested, set time index to recycle time records in restart
133! file.
134!
135 dai(ng)%Rindex=mod(dai(ng)%Rindex-1,2)+1
136!
137! Report.
138!
139# ifdef SOLVE3D
140# ifdef NESTING
141 IF (master) WRITE (stdout,10) kout, nout, dai(ng)%Rindex, ng
142# else
143 IF (master) WRITE (stdout,10) kout, nout, dai(ng)%Rindex
144# endif
145# else
146# ifdef NESTING
147 IF (master) WRITE (stdout,10) kout, dai(ng)%Rindex, ng
148# else
149 IF (master) WRITE (stdout,10) kout, dai(ng)%Rindex
150# endif
151# endif
152
153# ifdef SOLVE3D
154!
155! Write out time independent (unperturb) depths of RHO-points.
156!
157 scale=1.0_dp
158 gtype=gfactor*r3dvar
159 status=nf_fwrite3d(ng, inlm, dai(ng)%ncid, idpthr, &
160 & dai(ng)%Vid(idpthr), &
161 & 0, gtype, &
162 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
163# ifdef MASKING
164 & grid(ng) % rmask, &
165# endif
166 & grid(ng) % z0_r, &
167 & setfillval = .false.)
168 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
169 IF (master) THEN
170 WRITE (stdout,20) trim(vname(1,idpthr))
171 END IF
172 exit_flag=3
173 ioerror=status
174 RETURN
175 END IF
176!
177! Write out time independent (unperturb) depths of U-points.
178!
179 scale=1.0_dp
180 gtype=gfactor*u3dvar
181 DO k=1,n(ng)
182 DO j=jstr-1,jend+1
183 DO i=istru-1,iend+1
184 grid(ng)%z_v(i,j,k)=0.5_r8*(grid(ng)%z0_r(i-1,j,k)+ &
185 & grid(ng)%z0_r(i ,j,k))
186 END DO
187 END DO
188 END DO
189 status=nf_fwrite3d(ng, inlm, dai(ng)%ncid, idpthu, &
190 & dai(ng)%Vid(idpthu), &
191 & 0, gtype, &
192 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
193# ifdef MASKING
194 & grid(ng) % umask, &
195# endif
196 & grid(ng) % z_v, &
197 & setfillval = .false.)
198 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
199 IF (master) THEN
200 WRITE (stdout,20) trim(vname(1,idpthu))
201 END IF
202 exit_flag=3
203 ioerror=status
204 RETURN
205 END IF
206!
207! Write out time independent (unperturb) depths of V-points.
208!
209 scale=1.0_dp
210 gtype=gfactor*v3dvar
211 DO k=1,n(ng)
212 DO j=jstrv-1,jend+1
213 DO i=istr-1,iend+1
214 grid(ng)%z_v(i,j,k)=0.5_r8*(grid(ng)%z0_r(i,j-1,k)+ &
215 & grid(ng)%z0_r(i,j ,k))
216 END DO
217 END DO
218 END DO
219 status=nf_fwrite3d(ng, inlm, dai(ng)%ncid, idpthv, &
220 & dai(ng)%Vid(idpthv), &
221 & 0, gtype, &
222 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
223# ifdef MASKING
224 & grid(ng) % vmask, &
225# endif
226 & grid(ng) % z_v, &
227 & setfillval = .false.)
228 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
229 IF (master) THEN
230 WRITE (stdout,20) trim(vname(1,idpthv))
231 END IF
232 exit_flag=3
233 ioerror=status
234 RETURN
235 END IF
236!
237! Write out time independent (unperturb) depths of W-points.
238!
239 scale=1.0_dp
240 gtype=gfactor*w3dvar
241 status=nf_fwrite3d(ng, inlm, dai(ng)%ncid, idpthw, &
242 & dai(ng)%Vid(idpthw), &
243 & 0, gtype, &
244 & lbi, ubi, lbj, ubj, 0, n(ng), scale, &
245# ifdef MASKING
246 & grid(ng) % rmask, &
247# endif
248 & grid(ng) % z0_w, &
249 & setfillval = .false.)
250 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
251 IF (master) THEN
252 WRITE (stdout,20) trim(vname(1,idpthw))
253 END IF
254 exit_flag=3
255 ioerror=status
256 RETURN
257 END IF
258# endif
259!
260! Write out model time (s).
261!
262 CALL netcdf_put_fvar (ng, inlm, dai(ng)%name, &
263 & trim(vname(1,idtime)), time(ng:), &
264 & (/dai(ng)%Rindex/), (/1/), &
265 & ncid = dai(ng)%ncid, &
266 & varid = dai(ng)%Vid(idtime))
267 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
268!
269! Write out free-surface (m).
270!
271 scale=1.0_dp
272 gtype=gfactor*r2dvar
273 status=nf_fwrite2d(ng, inlm, dai(ng)%ncid, idfsur, &
274 & dai(ng)%Vid(idfsur), &
275 & dai(ng)%Rindex, gtype, &
276 & lbi, ubi, lbj, ubj, scale, &
277# ifdef MASKING
278 & grid(ng) % rmask, &
279# endif
280# if defined R4DVAR || defined SPLIT_R4DVAR
281# ifdef WET_DRY
282 & ocean(ng) % tl_zeta(:,:,kout), &
283 & setfillval = .false.)
284# else
285 & ocean(ng) % tl_zeta(:,:,kout))
286# endif
287# else
288# ifdef WET_DRY
289 & ocean(ng) % zeta(:,:,kout), &
290 & setfillval = .false.)
291# else
292 & ocean(ng) % zeta(:,:,kout))
293# endif
294# endif
295 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
296 IF (master) THEN
297 WRITE (stdout,30) trim(vname(1,idfsur)), dai(ng)%Rindex
298 END IF
299 exit_flag=3
300 ioerror=status
301 RETURN
302 END IF
303!
304! Write out 2D momentum component (m/s) in the XI-direction.
305!
306 scale=1.0_dp
307 gtype=gfactor*u2dvar
308 status=nf_fwrite2d(ng, inlm, dai(ng)%ncid, idubar, &
309 & dai(ng)%Vid(idubar), &
310 & dai(ng)%Rindex, gtype, &
311 & lbi, ubi, lbj, ubj, scale, &
312# ifdef MASKING
313 & grid(ng) % umask_full, &
314# endif
315# if defined R4DVAR || defined SPLIT_R4DVAR
316 & ocean(ng) % tl_ubar(:,:,kout))
317# else
318 & ocean(ng) % ubar(:,:,kout))
319# endif
320 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
321 IF (master) THEN
322 WRITE (stdout,30) trim(vname(1,idubar)), dai(ng)%Rindex
323 END IF
324 exit_flag=3
325 ioerror=status
326 RETURN
327 END IF
328!
329! Write out 2D momentum component (m/s) in the ETA-direction.
330!
331 scale=1.0_dp
332 gtype=gfactor*v2dvar
333 status=nf_fwrite2d(ng, inlm, dai(ng)%ncid, idvbar, &
334 & dai(ng)%Vid(idvbar), &
335 & dai(ng)%Rindex, gtype, &
336 & lbi, ubi, lbj, ubj, scale, &
337# ifdef MASKING
338 & grid(ng) % vmask_full, &
339# endif
340# if defined R4DVAR || defined SPLIT_R4DVAR
341 & ocean(ng) % tl_vbar(:,:,kout))
342# else
343 & ocean(ng) % vbar(:,:,kout))
344# endif
345 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
346 IF (master) THEN
347 WRITE (stdout,30) trim(vname(1,idvbar)), dai(ng)%Rindex
348 END IF
349 exit_flag=3
350 ioerror=status
351 RETURN
352 END IF
353
354# ifdef SOLVE3D
355!
356! Write out 3D momentum component (m/s) in the XI-direction.
357!
358 scale=1.0_dp
359 gtype=gfactor*u3dvar
360 status=nf_fwrite3d(ng, inlm, dai(ng)%ncid, iduvel, &
361 & dai(ng)%Vid(iduvel), &
362 & dai(ng)%Rindex, gtype, &
363 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
364# ifdef MASKING
365 & grid(ng) % umask_full, &
366# endif
367# if defined R4DVAR || defined SPLIT_R4DVAR
368 & ocean(ng) % tl_u(:,:,:,nout))
369# else
370 & ocean(ng) % u(:,:,:,nout))
371# endif
372 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
373 IF (master) THEN
374 WRITE (stdout,30) trim(vname(1,iduvel)), dai(ng)%Rindex
375 END IF
376 exit_flag=3
377 ioerror=status
378 RETURN
379 END IF
380!
381! Write out momentum component (m/s) in the ETA-direction.
382!
383 scale=1.0_dp
384 gtype=gfactor*v3dvar
385 status=nf_fwrite3d(ng, inlm, dai(ng)%ncid, idvvel, &
386 & dai(ng)%Vid(idvvel), &
387 & dai(ng)%Rindex, gtype, &
388 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
389# ifdef MASKING
390 & grid(ng) % vmask_full, &
391# endif
392# if defined R4DVAR || defined SPLIT_R4DVAR
393 & ocean(ng) % tl_v(:,:,:,nout))
394# else
395 & ocean(ng) % v(:,:,:,nout))
396# endif
397 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
398 IF (master) THEN
399 WRITE (stdout,30) trim(vname(1,idvvel)), dai(ng)%Rindex
400 END IF
401 exit_flag=3
402 ioerror=status
403 RETURN
404 END IF
405!
406! Write out tracer type variables.
407!
408 DO itrc=1,nt(ng)
409 scale=1.0_dp
410 gtype=gfactor*r3dvar
411 status=nf_fwrite3d(ng, inlm, dai(ng)%ncid, idtvar(itrc), &
412 & dai(ng)%Tid(itrc), &
413 & dai(ng)%Rindex, gtype, &
414 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
415# ifdef MASKING
416 & grid(ng) % rmask, &
417# endif
418# if defined R4DVAR || defined SPLIT_R4DVAR
419 & ocean(ng) % tl_t(:,:,:,nout,itrc))
420# else
421 & ocean(ng) % t(:,:,:,nout,itrc))
422# endif
423 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
424 IF (master) THEN
425 WRITE (stdout,30) trim(vname(1,idtvar(itrc))), dai(ng)%Rindex
426 END IF
427 exit_flag=3
428 ioerror=status
429 RETURN
430 END IF
431 END DO
432!
433! Write out vertical viscosity coefficient.
434!
435 scale=1.0_dp
436 gtype=gfactor*w3dvar
437 status=nf_fwrite3d(ng, inlm, dai(ng)%ncid, idvvis, &
438 & dai(ng)%Vid(idvvis), &
439 & dai(ng)%Rindex, gtype, &
440 & lbi, ubi, lbj, ubj, 0, n(ng), scale, &
441# ifdef MASKING
442 & grid(ng) % rmask, &
443# endif
444 & mixing(ng) % Akv, &
445 & setfillval = .false.)
446 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
447 IF (master) THEN
448 WRITE (stdout,30) trim(vname(1,idvvis)), dai(ng)%Rindex
449 END IF
450 exit_flag=3
451 ioerror=status
452 RETURN
453 END IF
454!
455! Write out vertical diffusion coefficient for potential temperature.
456!
457 scale=1.0_dp
458 gtype=gfactor*w3dvar
459 status=nf_fwrite3d(ng, inlm, dai(ng)%ncid, idtdif, &
460 & dai(ng)%Vid(idtdif), &
461 & dai(ng)%Rindex, gtype, &
462 & lbi, ubi, lbj, ubj, 0, n(ng), scale, &
463# ifdef MASKING
464 & grid(ng) % rmask, &
465# endif
466 & mixing(ng) % Akt(:,:,:,itemp), &
467 & setfillval = .false.)
468 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
469 IF (master) THEN
470 WRITE (stdout,30) trim(vname(1,idtdif)), dai(ng)%Rindex
471 END IF
472 exit_flag=3
473 ioerror=status
474 RETURN
475 END IF
476
477# ifdef SALINITY
478!
479! Write out vertical diffusion coefficient for salinity.
480!
481 scale=1.0_dp
482 gtype=gfactor*w3dvar
483 status=nf_fwrite3d(ng, inlm, dai(ng)%ncid, idsdif, &
484 & dai(ng)%Vid(idsdif), &
485 & dai(ng)%Rindex, gtype, &
486 & lbi, ubi, lbj, ubj, 0, n(ng), scale, &
487# ifdef MASKING
488 & grid(ng) % rmask, &
489# endif
490 & mixing(ng) % Akt(:,:,:,isalt), &
491 & setfillval = .false.)
492 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
493 IF (master) THEN
494 WRITE (stdout,30) trim(vname(1,idsdif)), dai(ng)%Rindex
495 END IF
496 exit_flag=3
497 ioerror=status
498 RETURN
499 END IF
500# endif
501# endif
502!
503!-----------------------------------------------------------------------
504! Synchronize restart NetCDF file to disk.
505!-----------------------------------------------------------------------
506!
507 CALL netcdf_sync (ng, inlm, dai(ng)%name, dai(ng)%ncid)
508 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
509!
510 10 FORMAT (2x,'WRT_DAI_NF90 - writing DA INI/RST', t42, &
511# ifdef SOLVE3D
512# ifdef NESTING
513 & 'fields (Index=',i1,',',i1,') in record = ',i0,t92,i2.2)
514# else
515 & 'fields (Index=',i1,',',i1,') in record = ',i0)
516# endif
517# else
518# ifdef NESTING
519 & 'fields (Index=',i1,') in record = ',i0,t92,i2.2)
520# else
521 & 'fields (Index=',i1,') in record = ',i0)
522# endif
523# endif
524 20 FORMAT (/,' WRT_DAI_NF90 - error while writing variable: ',a, &
525 & /,11x,'into DA initial/restart NetCDF file.')
526 30 FORMAT (/,' WRT_DAI_NF90 - error while writing variable: ',a, &
527 & /,11x,'into DA initial/rstart NetCDF file for time ', &
528 & 'record: ',i0)
529!
530 RETURN
subroutine, public netcdf_sync(ng, model, ncname, ncid)

References mod_iounits::dai, mod_scalars::exit_flag, strings_mod::founderror(), mod_grid::grid, mod_ncparam::idfsur, mod_ncparam::idpthr, mod_ncparam::idpthu, mod_ncparam::idpthv, mod_ncparam::idpthw, mod_ncparam::idsdif, mod_ncparam::idtdif, mod_ncparam::idtime, mod_ncparam::idtvar, mod_ncparam::idubar, mod_ncparam::iduvel, mod_ncparam::idvbar, mod_ncparam::idvvel, mod_ncparam::idvvis, mod_param::inlm, mod_iounits::ioerror, mod_scalars::isalt, mod_scalars::itemp, mod_parallel::master, mod_mixing::mixing, mod_param::n, mod_netcdf::netcdf_sync(), 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, mod_ncparam::vname, and mod_param::w3dvar.

Referenced by wrt_dai().

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

◆ wrt_dai_pio()

subroutine, private wrt_dai_mod::wrt_dai_pio ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj )
private

Definition at line 536 of file wrt_dai.F.

538!***********************************************************************
539!
541!
542! Imported variable declarations.
543!
544 integer, intent(in) :: ng, tile
545 integer, intent(in) :: LBi, UBi, LBj, UBj
546!
547! Local variable declarations.
548!
549 integer :: i, j, k, itrc
550 integer :: Fcount, status
551!
552 real(dp) :: scale
553!
554 character (len=*), parameter :: MyFile = &
555 & __FILE__//", wrt_dai_pio"
556!
557 TYPE (IO_desc_t), pointer :: ioDesc
558
559# include "set_bounds.h"
560!
561 sourcefile=myfile
562!
563!-----------------------------------------------------------------------
564! Write out Data Assimilation initial/restart fields.
565!-----------------------------------------------------------------------
566!
567 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
568!
569! Set time record index.
570!
571 dai(ng)%Rindex=dai(ng)%Rindex+1
572 fcount=dai(ng)%Fcount
573 dai(ng)%Nrec(fcount)=dai(ng)%Nrec(fcount)+1
574!
575! If requested, set time index to recycle time records in restart
576! file.
577!
578 dai(ng)%Rindex=mod(dai(ng)%Rindex-1,2)+1
579!
580! Report
581!
582# ifdef SOLVE3D
583# ifdef NESTING
584 IF (master) WRITE (stdout,10) kout, nout, dai(ng)%Rindex, ng
585# else
586 IF (master) WRITE (stdout,10) kout, nout, dai(ng)%Rindex
587# endif
588# else
589# ifdef NESTING
590 IF (master) WRITE (stdout,10) kout, dai(ng)%Rindex, ng
591# else
592 IF (master) WRITE (stdout,10) kout, dai(ng)%Rindex
593# endif
594# endif
595
596# ifdef SOLVE3D
597!
598! Write out time independent (unperturb) depths of RHO-points.
599!
600 scale=1.0_dp
601 IF (dai(ng)%pioVar(idpthr)%dkind.eq.pio_double) THEN
602 iodesc => iodesc_dp_r3dvar(ng)
603 ELSE
604 iodesc => iodesc_sp_r3dvar(ng)
605 END IF
606!
607 status=nf_fwrite3d(ng, inlm, dai(ng)%pioFile, idpthr, &
608 & dai(ng)%pioVar(idpthr), 0, &
609 & iodesc, &
610 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
611# ifdef MASKING
612 & grid(ng) % rmask, &
613# endif
614 & grid(ng) % z0_r, &
615 & setfillval = .false.)
616 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
617 IF (master) THEN
618 WRITE (stdout,20) trim(vname(1,idpthr))
619 END IF
620 exit_flag=3
621 ioerror=status
622 RETURN
623 END IF
624!
625! Write out time independent (unperturb) depths of U-points.
626!
627 DO k=1,n(ng)
628 DO j=jstr-1,jend+1
629 DO i=istru-1,iend+1
630 grid(ng)%z_v(i,j,k)=0.5_r8*(grid(ng)%z0_r(i-1,j,k)+ &
631 & grid(ng)%z0_r(i ,j,k))
632 END DO
633 END DO
634 END DO
635!
636 scale=1.0_dp
637 IF (dai(ng)%pioVar(idpthu)%dkind.eq.pio_double) THEN
638 iodesc => iodesc_dp_u3dvar(ng)
639 ELSE
640 iodesc => iodesc_sp_u3dvar(ng)
641 END IF
642!
643 status=nf_fwrite3d(ng, inlm, dai(ng)%pioFile, idpthu, &
644 & dai(ng)%pioVar(idpthu), 0, &
645 & iodesc, &
646 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
647# ifdef MASKING
648 & grid(ng) % umask, &
649# endif
650 & grid(ng) % z_v, &
651 & setfillval = .false.)
652 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
653 IF (master) THEN
654 WRITE (stdout,20) trim(vname(1,idpthu))
655 END IF
656 exit_flag=3
657 ioerror=status
658 RETURN
659 END IF
660!
661! Write out time independent (unperturb) depths of V-points.
662!
663 DO k=1,n(ng)
664 DO j=jstrv-1,jend+1
665 DO i=istr-1,iend+1
666 grid(ng)%z_v(i,j,k)=0.5_r8*(grid(ng)%z0_r(i,j-1,k)+ &
667 & grid(ng)%z0_r(i,j ,k))
668 END DO
669 END DO
670 END DO
671!
672 scale=1.0_dp
673 IF (dai(ng)%pioVar(idpthv)%dkind.eq.pio_double) THEN
674 iodesc => iodesc_dp_v3dvar(ng)
675 ELSE
676 iodesc => iodesc_sp_v3dvar(ng)
677 END IF
678 status=nf_fwrite3d(ng, inlm, dai(ng)%pioFile, idpthv, &
679 & dai(ng)%pioVar(idpthv), 0, &
680 & iodesc, &
681 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
682# ifdef MASKING
683 & grid(ng) % vmask, &
684# endif
685 & grid(ng) % z_v, &
686 & setfillval = .false.)
687 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
688 IF (master) THEN
689 WRITE (stdout,20) trim(vname(1,idpthv))
690 END IF
691 exit_flag=3
692 ioerror=status
693 RETURN
694 END IF
695!
696! Write out time independent (unperturb) depths of W-points.
697!
698 scale=1.0_dp
699 IF (dai(ng)%pioVar(idpthw)%dkind.eq.pio_double) THEN
700 iodesc => iodesc_dp_w3dvar(ng)
701 ELSE
702 iodesc => iodesc_sp_w3dvar(ng)
703 END IF
704!
705 status=nf_fwrite3d(ng, inlm, dai(ng)%pioFile, idpthw, &
706 & dai(ng)%pioVar(idpthw), 0, &
707 & iodesc, &
708 & lbi, ubi, lbj, ubj, 0, n(ng), scale, &
709# ifdef MASKING
710 & grid(ng) % rmask, &
711# endif
712 & grid(ng) % z0_w, &
713 & setfillval = .false.)
714 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
715 IF (master) THEN
716 WRITE (stdout,20) trim(vname(1,idpthw))
717 END IF
718 exit_flag=3
719 ioerror=status
720 RETURN
721 END IF
722# endif
723!
724! Write out model time (s).
725!
726 CALL pio_netcdf_put_fvar (ng, inlm, dai(ng)%name, &
727 & trim(vname(1,idtime)), time(ng:), &
728 & (/dai(ng)%Rindex/), (/1/), &
729 & piofile = dai(ng)%pioFile, &
730 & piovar = dai(ng)%pioVar(idtime)%vd)
731 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
732!
733! Write out free-surface (m).
734!
735 scale=1.0_dp
736 IF (dai(ng)%pioVar(idfsur)%dkind.eq.pio_double) THEN
737 iodesc => iodesc_dp_r2dvar(ng)
738 ELSE
739 iodesc => iodesc_sp_r2dvar(ng)
740 END IF
741!
742 status=nf_fwrite2d(ng, inlm, dai(ng)%pioFile, idfsur, &
743 & dai(ng)%pioVar(idfsur), &
744 & dai(ng)%Rindex, &
745 & iodesc, &
746 & lbi, ubi, lbj, ubj, scale, &
747# ifdef MASKING
748 & grid(ng) % rmask, &
749# endif
750# if defined R4DVAR || defined SPLIT_R4DVAR
751# ifdef WET_DRY
752 & ocean(ng) % tl_zeta(:,:,kout), &
753 & setfillval = .false.)
754# else
755 & ocean(ng) % tl_zeta(:,:,kout))
756# endif
757# else
758# ifdef WET_DRY
759 & ocean(ng) % zeta(:,:,kout), &
760 & setfillval = .false.)
761# else
762 & ocean(ng) % zeta(:,:,kout))
763# endif
764# endif
765 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
766 IF (master) THEN
767 WRITE (stdout,30) trim(vname(1,idfsur)), dai(ng)%Rindex
768 END IF
769 exit_flag=3
770 ioerror=status
771 RETURN
772 END IF
773!
774! Write out 2D momentum component (m/s) in the XI-direction.
775!
776 scale=1.0_dp
777 IF (dai(ng)%pioVar(idubar)%dkind.eq.pio_double) THEN
778 iodesc => iodesc_dp_u2dvar(ng)
779 ELSE
780 iodesc => iodesc_sp_u2dvar(ng)
781 END IF
782!
783 status=nf_fwrite2d(ng, inlm, dai(ng)%pioFile, idubar, &
784 & dai(ng)%pioVar(idubar), &
785 & dai(ng)%Rindex, &
786 & iodesc, &
787 & lbi, ubi, lbj, ubj, scale, &
788# ifdef MASKING
789 & grid(ng) % umask_full, &
790# endif
791# if defined R4DVAR || defined SPLIT_R4DVAR
792 & ocean(ng) % tl_ubar(:,:,kout))
793# else
794 & ocean(ng) % ubar(:,:,kout))
795# endif
796 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
797 IF (master) THEN
798 WRITE (stdout,30) trim(vname(1,idubar)), dai(ng)%Rindex
799 END IF
800 exit_flag=3
801 ioerror=status
802 RETURN
803 END IF
804!
805! Write out 2D momentum component (m/s) in the ETA-direction.
806!
807 scale=1.0_dp
808 IF (dai(ng)%pioVar(idvbar)%dkind.eq.pio_double) THEN
809 iodesc => iodesc_dp_v2dvar(ng)
810 ELSE
811 iodesc => iodesc_sp_v2dvar(ng)
812 END IF
813!
814 status=nf_fwrite2d(ng, inlm, dai(ng)%pioFile, idvbar, &
815 & dai(ng)%pioVar(idvbar), &
816 & dai(ng)%Rindex, &
817 & iodesc, &
818 & lbi, ubi, lbj, ubj, scale, &
819# ifdef MASKING
820 & grid(ng) % vmask_full, &
821# endif
822# if defined R4DVAR || defined SPLIT_R4DVAR
823 & ocean(ng) % tl_vbar(:,:,kout))
824# else
825 & ocean(ng) % vbar(:,:,kout))
826# endif
827 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
828 IF (master) THEN
829 WRITE (stdout,30) trim(vname(1,idvbar)), dai(ng)%Rindex
830 END IF
831 exit_flag=3
832 ioerror=status
833 RETURN
834 END IF
835
836# ifdef SOLVE3D
837!
838! Write out 3D momentum component (m/s) in the XI-direction.
839!
840 scale=1.0_dp
841 IF (dai(ng)%pioVar(iduvel)%dkind.eq.pio_double) THEN
842 iodesc => iodesc_dp_u3dvar(ng)
843 ELSE
844 iodesc => iodesc_sp_u3dvar(ng)
845 END IF
846!
847 status=nf_fwrite3d(ng, inlm, dai(ng)%pioFile, iduvel, &
848 & dai(ng)%pioVar(iduvel), &
849 & dai(ng)%Rindex, &
850 & iodesc, &
851 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
852# ifdef MASKING
853 & grid(ng) % umask_full, &
854# endif
855# if defined R4DVAR || defined SPLIT_R4DVAR
856 & ocean(ng) % tl_u(:,:,:,nout))
857# else
858 & ocean(ng) % u(:,:,:,nout))
859# endif
860 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
861 IF (master) THEN
862 WRITE (stdout,30) trim(vname(1,iduvel)), dai(ng)%Rindex
863 END IF
864 exit_flag=3
865 ioerror=status
866 RETURN
867 END IF
868!
869! Write out momentum component (m/s) in the ETA-direction.
870!
871 scale=1.0_dp
872 IF (dai(ng)%pioVar(idfsur)%dkind.eq.pio_double) THEN
873 iodesc => iodesc_dp_v3dvar(ng)
874 ELSE
875 iodesc => iodesc_sp_v3dvar(ng)
876 END IF
877!
878 status=nf_fwrite3d(ng, inlm, dai(ng)%pioFile, idvvel, &
879 & dai(ng)%pioVar(idvvel), &
880 & dai(ng)%Rindex, &
881 & iodesc, &
882 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
883# ifdef MASKING
884 & grid(ng) % vmask_full, &
885# endif
886# if defined R4DVAR || defined SPLIT_R4DVAR
887 & ocean(ng) % tl_v(:,:,:,nout))
888# else
889 & ocean(ng) % v(:,:,:,nout))
890# endif
891 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
892 IF (master) THEN
893 WRITE (stdout,30) trim(vname(1,idvvel)), dai(ng)%Rindex
894 END IF
895 exit_flag=3
896 ioerror=status
897 RETURN
898 END IF
899!
900! Write out tracer type variables.
901!
902 DO itrc=1,nt(ng)
903 scale=1.0_dp
904 IF (dai(ng)%pioTrc(itrc)%dkind.eq.pio_double) THEN
905 iodesc => iodesc_dp_r3dvar(ng)
906 ELSE
907 iodesc => iodesc_sp_r3dvar(ng)
908 END IF
909!
910 status=nf_fwrite3d(ng, inlm, dai(ng)%pioFile, idtvar(itrc), &
911 & dai(ng)%pioTrc(itrc), &
912 & dai(ng)%Rindex, &
913 & iodesc, &
914 & lbi, ubi, lbj, ubj, 1, n(ng), scale, &
915# ifdef MASKING
916 & grid(ng) % rmask, &
917# endif
918# if defined R4DVAR || defined SPLIT_R4DVAR
919 & ocean(ng) % tl_t(:,:,:,nout,itrc))
920# else
921 & ocean(ng) % t(:,:,:,nout,itrc))
922# endif
923 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
924 IF (master) THEN
925 WRITE (stdout,30) trim(vname(1,idtvar(itrc))), dai(ng)%Rindex
926 END IF
927 exit_flag=3
928 ioerror=status
929 RETURN
930 END IF
931 END DO
932!
933! Write out vertical viscosity coefficient.
934!
935 scale=1.0_dp
936 IF (dai(ng)%pioVar(idvvis)%dkind.eq.pio_double) THEN
937 iodesc => iodesc_dp_w3dvar(ng)
938 ELSE
939 iodesc => iodesc_sp_w3dvar(ng)
940 END IF
941!
942 status=nf_fwrite3d(ng, inlm, dai(ng)%pioFile, idvvis, &
943 & dai(ng)%pioVar(idvvis), &
944 & dai(ng)%Rindex, &
945 & iodesc, &
946 & lbi, ubi, lbj, ubj, 0, n(ng), scale, &
947# ifdef MASKING
948 & grid(ng) % rmask, &
949# endif
950 & mixing(ng) % Akv, &
951 & setfillval = .false.)
952 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
953 IF (master) THEN
954 WRITE (stdout,30) trim(vname(1,idvvis)), dai(ng)%Rindex
955 END IF
956 exit_flag=3
957 ioerror=status
958 RETURN
959 END IF
960!
961! Write out vertical diffusion coefficient for potential temperature.
962!
963 scale=1.0_dp
964 IF (dai(ng)%pioVar(idtdif)%dkind.eq.pio_double) THEN
965 iodesc => iodesc_dp_w3dvar(ng)
966 ELSE
967 iodesc => iodesc_sp_w3dvar(ng)
968 END IF
969!
970 status=nf_fwrite3d(ng, inlm, dai(ng)%pioFile, idtdif, &
971 & dai(ng)%pioVar(idtdif), &
972 & dai(ng)%Rindex, &
973 & iodesc, &
974 & lbi, ubi, lbj, ubj, 0, n(ng), scale, &
975# ifdef MASKING
976 & grid(ng) % rmask, &
977# endif
978 & mixing(ng) % Akt(:,:,:,itemp), &
979 & setfillval = .false.)
980 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
981 IF (master) THEN
982 WRITE (stdout,30) trim(vname(1,idtdif)), dai(ng)%Rindex
983 END IF
984 exit_flag=3
985 ioerror=status
986 RETURN
987 END IF
988
989# ifdef SALINITY
990!
991! Write out vertical diffusion coefficient for salinity.
992!
993 scale=1.0_dp
994 IF (dai(ng)%pioVar(idsdif)%dkind.eq.pio_double) THEN
995 iodesc => iodesc_dp_w3dvar(ng)
996 ELSE
997 iodesc => iodesc_sp_w3dvar(ng)
998 END IF
999!
1000 status=nf_fwrite3d(ng, inlm, dai(ng)%pioFile, idsdif, &
1001 & dai(ng)%pioVar(idsdif), &
1002 & dai(ng)%Rindex, &
1003 & iodesc, &
1004 & lbi, ubi, lbj, ubj, 0, n(ng), scale, &
1005# ifdef MASKING
1006 & grid(ng) % rmask, &
1007# endif
1008 & mixing(ng) % Akt(:,:,:,isalt), &
1009 & setfillval = .false.)
1010 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1011 IF (master) THEN
1012 WRITE (stdout,30) trim(vname(1,idsdif)), dai(ng)%Rindex
1013 END IF
1014 exit_flag=3
1015 ioerror=status
1016 RETURN
1017 END IF
1018# endif
1019# endif
1020!
1021!-----------------------------------------------------------------------
1022! Synchronize restart NetCDF file to disk.
1023!-----------------------------------------------------------------------
1024!
1025 CALL pio_netcdf_sync (ng, inlm, dai(ng)%name, dai(ng)%pioFile)
1026 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
1027!
1028 10 FORMAT (2x,'WRT_DAI_PIO - writing DA INI/RST', t42, &
1029# ifdef SOLVE3D
1030# ifdef NESTING
1031 & 'fields (Index=',i1,',',i1,') in record = ',i0,t92,i2.2)
1032# else
1033 & 'fields (Index=',i1,',',i1,') in record = ',i0)
1034# endif
1035# else
1036# ifdef NESTING
1037 & 'fields (Index=',i1,') in record = ',i0,t92,i2.2)
1038# else
1039 & 'fields (Index=',i1,') in record = ',i0)
1040# endif
1041# endif
1042 20 FORMAT (/,' WRT_DAI_PIO - error while writing variable: ',a, &
1043 & /,11x,'into DA initial/restart NetCDF file.')
1044 30 FORMAT (/,' WRT_DAI_PIO - error while writing variable: ',a, &
1045 & /,11x,'into DA initial/rstart NetCDF file for time ', &
1046 & 'record: ',i0)
1047!
1048 RETURN
type(io_desc_t), dimension(:), pointer iodesc_sp_w3dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_r3dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_u3dvar
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_u3dvar
type(io_desc_t), dimension(:), pointer iodesc_sp_r2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_u2dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_v3dvar
type(io_desc_t), dimension(:), pointer iodesc_dp_w3dvar
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_sp_v2dvar

References mod_iounits::dai, mod_scalars::exit_flag, strings_mod::founderror(), mod_grid::grid, mod_ncparam::idfsur, mod_ncparam::idpthr, mod_ncparam::idpthu, mod_ncparam::idpthv, mod_ncparam::idpthw, mod_ncparam::idsdif, mod_ncparam::idtdif, mod_ncparam::idtime, mod_ncparam::idtvar, mod_ncparam::idubar, mod_ncparam::iduvel, mod_ncparam::idvbar, mod_ncparam::idvvel, mod_ncparam::idvvis, mod_param::inlm, mod_pio_netcdf::iodesc_dp_r2dvar, mod_pio_netcdf::iodesc_dp_r3dvar, mod_pio_netcdf::iodesc_dp_u2dvar, mod_pio_netcdf::iodesc_dp_u3dvar, mod_pio_netcdf::iodesc_dp_v2dvar, mod_pio_netcdf::iodesc_dp_v3dvar, mod_pio_netcdf::iodesc_dp_w3dvar, mod_pio_netcdf::iodesc_sp_r2dvar, mod_pio_netcdf::iodesc_sp_r3dvar, mod_pio_netcdf::iodesc_sp_u2dvar, mod_pio_netcdf::iodesc_sp_u3dvar, mod_pio_netcdf::iodesc_sp_v2dvar, mod_pio_netcdf::iodesc_sp_v3dvar, mod_pio_netcdf::iodesc_sp_w3dvar, mod_iounits::ioerror, mod_scalars::isalt, mod_scalars::itemp, mod_parallel::master, mod_mixing::mixing, mod_param::n, 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_dai().

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