3#if defined FOUR_DVAR || defined ENKF_RESTART
38# if defined PIO_LIB && defined DISTRIBUTE
50 integer,
intent(in) :: ng, tile
54 integer :: lbi, ubi, lbj, ubj
56 character (len=*),
parameter :: myfile = &
68 SELECT CASE (
dia(ng)%IOtype)
73# if defined PIO_LIB && defined DISTRIBUTE
84 10
FORMAT (
' WRT_DAI - Illegal output file type, io_type = ',i0, &
85 & /,11x,
'Check KeyWord ''OUT_LIB'' in ''roms.in''.')
99 integer,
intent(in) :: ng, tile
100 integer,
intent(in) :: lbi, ubi, lbj, ubj
104 integer :: i, j, k, itrc
105 integer :: fcount, gfactor, gtype, status, varid
109 character (len=*),
parameter :: myfile = &
110 & __FILE__//
", wrt_dai_nf90"
112# include "set_bounds.h"
128 dai(ng)%Rindex=
dai(ng)%Rindex+1
129 fcount=
dai(ng)%Fcount
130 dai(ng)%Nrec(fcount)=
dai(ng)%Nrec(fcount)+1
135 dai(ng)%Rindex=mod(
dai(ng)%Rindex-1,2)+1
162 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
164 &
grid(ng) % rmask, &
167 & setfillval = .false.)
168 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
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))
192 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
194 &
grid(ng) % umask, &
197 & setfillval = .false.)
198 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
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))
222 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
224 &
grid(ng) % vmask, &
227 & setfillval = .false.)
228 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
244 & lbi, ubi, lbj, ubj, 0,
n(ng), scale, &
246 &
grid(ng) % rmask, &
249 & setfillval = .false.)
250 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
264 & (/
dai(ng)%Rindex/), (/1/), &
265 & ncid =
dai(ng)%ncid, &
275 &
dai(ng)%Rindex, gtype, &
276 & lbi, ubi, lbj, ubj, scale, &
278 &
grid(ng) % rmask, &
280# if defined R4DVAR || defined SPLIT_R4DVAR
282 &
ocean(ng) % tl_zeta(:,:,kout), &
283 & setfillval = .false.)
285 &
ocean(ng) % tl_zeta(:,:,kout))
289 &
ocean(ng) % zeta(:,:,kout), &
290 & setfillval = .false.)
292 &
ocean(ng) % zeta(:,:,kout))
295 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
310 &
dai(ng)%Rindex, gtype, &
311 & lbi, ubi, lbj, ubj, scale, &
313 &
grid(ng) % umask_full, &
315# if defined R4DVAR || defined SPLIT_R4DVAR
316 &
ocean(ng) % tl_ubar(:,:,kout))
318 &
ocean(ng) % ubar(:,:,kout))
320 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
335 &
dai(ng)%Rindex, gtype, &
336 & lbi, ubi, lbj, ubj, scale, &
338 &
grid(ng) % vmask_full, &
340# if defined R4DVAR || defined SPLIT_R4DVAR
341 &
ocean(ng) % tl_vbar(:,:,kout))
343 &
ocean(ng) % vbar(:,:,kout))
345 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
362 &
dai(ng)%Rindex, gtype, &
363 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
365 &
grid(ng) % umask_full, &
367# if defined R4DVAR || defined SPLIT_R4DVAR
368 &
ocean(ng) % tl_u(:,:,:,nout))
370 &
ocean(ng) % u(:,:,:,nout))
372 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
387 &
dai(ng)%Rindex, gtype, &
388 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
390 &
grid(ng) % vmask_full, &
392# if defined R4DVAR || defined SPLIT_R4DVAR
393 &
ocean(ng) % tl_v(:,:,:,nout))
395 &
ocean(ng) % v(:,:,:,nout))
397 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
412 &
dai(ng)%Tid(itrc), &
413 &
dai(ng)%Rindex, gtype, &
414 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
416 &
grid(ng) % rmask, &
418# if defined R4DVAR || defined SPLIT_R4DVAR
419 &
ocean(ng) % tl_t(:,:,:,nout,itrc))
421 &
ocean(ng) % t(:,:,:,nout,itrc))
423 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
439 &
dai(ng)%Rindex, gtype, &
440 & lbi, ubi, lbj, ubj, 0,
n(ng), scale, &
442 &
grid(ng) % rmask, &
445 & setfillval = .false.)
446 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
461 &
dai(ng)%Rindex, gtype, &
462 & lbi, ubi, lbj, ubj, 0,
n(ng), scale, &
464 &
grid(ng) % rmask, &
467 & setfillval = .false.)
468 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
485 &
dai(ng)%Rindex, gtype, &
486 & lbi, ubi, lbj, ubj, 0,
n(ng), scale, &
488 &
grid(ng) % rmask, &
491 & setfillval = .false.)
492 IF (
founderror(status, nf90_noerr, __line__, myfile))
THEN
510 10
FORMAT (2x,
'WRT_DAI_NF90 - writing DA INI/RST', t42, &
513 &
'fields (Index=',i1,
',',i1,
') in record = ',i0,t92,i2.2)
515 &
'fields (Index=',i1,
',',i1,
') in record = ',i0)
519 &
'fields (Index=',i1,
') in record = ',i0,t92,i2.2)
521 &
'fields (Index=',i1,
') in record = ',i0)
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 ', &
533# if defined PIO_LIB && defined DISTRIBUTE
537 & LBi, UBi, LBj, UBj)
544 integer,
intent(in) :: ng, tile
545 integer,
intent(in) :: lbi, ubi, lbj, ubj
549 integer :: i, j, k, itrc
550 integer :: fcount, status
554 character (len=*),
parameter :: myfile = &
555 & __FILE__//
", wrt_dai_pio"
557 TYPE (io_desc_t),
pointer :: iodesc
559# include "set_bounds.h"
571 dai(ng)%Rindex=
dai(ng)%Rindex+1
572 fcount=
dai(ng)%Fcount
573 dai(ng)%Nrec(fcount)=
dai(ng)%Nrec(fcount)+1
578 dai(ng)%Rindex=mod(
dai(ng)%Rindex-1,2)+1
601 IF (
dai(ng)%pioVar(
idpthr)%dkind.eq.pio_double)
THEN
610 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
612 &
grid(ng) % rmask, &
615 & setfillval = .false.)
616 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
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))
637 IF (
dai(ng)%pioVar(
idpthu)%dkind.eq.pio_double)
THEN
646 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
648 &
grid(ng) % umask, &
651 & setfillval = .false.)
652 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
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))
673 IF (
dai(ng)%pioVar(
idpthv)%dkind.eq.pio_double)
THEN
681 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
683 &
grid(ng) % vmask, &
686 & setfillval = .false.)
687 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
699 IF (
dai(ng)%pioVar(
idpthw)%dkind.eq.pio_double)
THEN
708 & lbi, ubi, lbj, ubj, 0,
n(ng), scale, &
710 &
grid(ng) % rmask, &
713 & setfillval = .false.)
714 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
728 & (/
dai(ng)%Rindex/), (/1/), &
729 & piofile =
dai(ng)%pioFile, &
736 IF (
dai(ng)%pioVar(
idfsur)%dkind.eq.pio_double)
THEN
746 & lbi, ubi, lbj, ubj, scale, &
748 &
grid(ng) % rmask, &
750# if defined R4DVAR || defined SPLIT_R4DVAR
752 &
ocean(ng) % tl_zeta(:,:,kout), &
753 & setfillval = .false.)
755 &
ocean(ng) % tl_zeta(:,:,kout))
759 &
ocean(ng) % zeta(:,:,kout), &
760 & setfillval = .false.)
762 &
ocean(ng) % zeta(:,:,kout))
765 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
777 IF (
dai(ng)%pioVar(
idubar)%dkind.eq.pio_double)
THEN
787 & lbi, ubi, lbj, ubj, scale, &
789 &
grid(ng) % umask_full, &
791# if defined R4DVAR || defined SPLIT_R4DVAR
792 &
ocean(ng) % tl_ubar(:,:,kout))
794 &
ocean(ng) % ubar(:,:,kout))
796 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
808 IF (
dai(ng)%pioVar(
idvbar)%dkind.eq.pio_double)
THEN
818 & lbi, ubi, lbj, ubj, scale, &
820 &
grid(ng) % vmask_full, &
822# if defined R4DVAR || defined SPLIT_R4DVAR
823 &
ocean(ng) % tl_vbar(:,:,kout))
825 &
ocean(ng) % vbar(:,:,kout))
827 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
841 IF (
dai(ng)%pioVar(
iduvel)%dkind.eq.pio_double)
THEN
851 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
853 &
grid(ng) % umask_full, &
855# if defined R4DVAR || defined SPLIT_R4DVAR
856 &
ocean(ng) % tl_u(:,:,:,nout))
858 &
ocean(ng) % u(:,:,:,nout))
860 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
872 IF (
dai(ng)%pioVar(
idfsur)%dkind.eq.pio_double)
THEN
882 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
884 &
grid(ng) % vmask_full, &
886# if defined R4DVAR || defined SPLIT_R4DVAR
887 &
ocean(ng) % tl_v(:,:,:,nout))
889 &
ocean(ng) % v(:,:,:,nout))
891 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
904 IF (
dai(ng)%pioTrc(itrc)%dkind.eq.pio_double)
THEN
911 &
dai(ng)%pioTrc(itrc), &
914 & lbi, ubi, lbj, ubj, 1,
n(ng), scale, &
916 &
grid(ng) % rmask, &
918# if defined R4DVAR || defined SPLIT_R4DVAR
919 &
ocean(ng) % tl_t(:,:,:,nout,itrc))
921 &
ocean(ng) % t(:,:,:,nout,itrc))
923 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
936 IF (
dai(ng)%pioVar(
idvvis)%dkind.eq.pio_double)
THEN
946 & lbi, ubi, lbj, ubj, 0,
n(ng), scale, &
948 &
grid(ng) % rmask, &
951 & setfillval = .false.)
952 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
964 IF (
dai(ng)%pioVar(
idtdif)%dkind.eq.pio_double)
THEN
974 & lbi, ubi, lbj, ubj, 0,
n(ng), scale, &
976 &
grid(ng) % rmask, &
979 & setfillval = .false.)
980 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
994 IF (
dai(ng)%pioVar(
idsdif)%dkind.eq.pio_double)
THEN
1004 & lbi, ubi, lbj, ubj, 0,
n(ng), scale, &
1006 &
grid(ng) % rmask, &
1009 & setfillval = .false.)
1010 IF (
founderror(status, pio_noerr, __line__, myfile))
THEN
1028 10
FORMAT (2x,
'WRT_DAI_PIO - writing DA INI/RST', t42, &
1031 &
'fields (Index=',i1,
',',i1,
') in record = ',i0,t92,i2.2)
1033 &
'fields (Index=',i1,
',',i1,
') in record = ',i0)
1037 &
'fields (Index=',i1,
') in record = ',i0,t92,i2.2)
1039 &
'fields (Index=',i1,
') in record = ',i0)
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 ', &
type(t_grid), dimension(:), allocatable grid
type(t_io), dimension(:), allocatable dai
character(len=256) sourcefile
type(t_io), dimension(:), allocatable dia
type(t_mixing), dimension(:), allocatable mixing
integer, parameter io_nf90
integer, parameter io_pio
integer, dimension(:), allocatable idtvar
character(len=maxlen), dimension(6, 0:nv) vname
subroutine, public netcdf_sync(ng, model, ncname, ncid)
type(t_ocean), dimension(:), allocatable ocean
integer, dimension(:), allocatable n
type(t_bounds), dimension(:), allocatable bounds
integer, parameter r3dvar
integer, parameter u3dvar
integer, parameter u2dvar
integer, parameter w3dvar
integer, dimension(:), allocatable nt
integer, parameter r2dvar
integer, parameter v2dvar
integer, parameter v3dvar
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
real(dp), dimension(:), allocatable time
logical function, public founderror(flag, noerr, line, routine)
subroutine, private wrt_dai_pio(ng, tile, lbi, ubi, lbj, ubj)
subroutine, private wrt_dai_nf90(ng, tile, lbi, ubi, lbj, ubj)
subroutine, public wrt_dai(ng, tile)