Skip to content

Commit

Permalink
modfielddump: correct handling of klow and ncoarse
Browse files Browse the repository at this point in the history
klow: make the zt and zm coordinates match the selected output range
klow and ncoarse: fix some of the optional 3D output fields
  • Loading branch information
fjansson committed Dec 21, 2023
1 parent 282ba98 commit d35a03e
Show file tree
Hide file tree
Showing 2 changed files with 23 additions and 19 deletions.
17 changes: 7 additions & 10 deletions src/modfielddump.f90
Original file line number Diff line number Diff line change
Expand Up @@ -42,7 +42,7 @@ module modfielddump

real :: dtav, tmin, tmax
integer(kind=longint) :: idtav,tnext,itmax,itmin
integer :: klow,khigh,ncoarse=-1
integer :: klow,khigh,ncoarse=1
logical :: lfielddump= .false. !< switch to enable the fielddump (on/off)
logical :: ldiracc = .false. !< switch for doing direct access writing (on/off)
logical :: lbinary = .false. !< switch for doing direct access writing (on/off)
Expand Down Expand Up @@ -132,9 +132,6 @@ subroutine initfielddump
call D_MPI_BCAST(ltntrl ,1,0,comm3d,ierr)
call D_MPI_BCAST(lsv ,100,0,comm3d,ierr)

if (ncoarse==-1) then
ncoarse = 1
end if
idtav = dtav/tres
itmin = tmin/tres
itmax = tmax/tres
Expand Down Expand Up @@ -254,7 +251,7 @@ subroutine initfielddump

if (nrec==0) then
call define_nc( ncid, 1, tncname)
call writestat_dims_nc(ncid, ncoarse)
call writestat_dims_nc(ncid, ncoarse, klow)
end if
call define_nc( ncid, NVar, ncname(1:NVar,:))
! must slice ncname here because define_nc expects first dimension to be NVar
Expand Down Expand Up @@ -411,7 +408,7 @@ subroutine fielddump
if (lnetcdf .and. lbuoy) then
vars(:,:,:,ind_buoy) = thv0h(2:i1:ncoarse,2:j1:ncoarse,klow:khigh)
do k=klow,khigh
vars(:,:,k,ind_buoy) = vars(:,:,k,ind_buoy) - thvh(k)
vars(:,:,k-klow+1,ind_buoy) = vars(:,:,k-klow+1,ind_buoy) - thvh(k)
end do
end if

Expand Down Expand Up @@ -456,7 +453,7 @@ subroutine fielddump
do k = klow, khigh
do j = 2, j1, ncoarse
do i = 2, i1, ncoarse
vars(i,j,k,ind_hur) = 100 * (qt0(i,j,k) - ql0(i,j,k)) / &
vars(i/ncoarse,j/ncoarse,k-klow+1,ind_hur) = 100 * (qt0(i,j,k) - ql0(i,j,k)) / &
qsat_tab(tmp0(i,j,k), presf(k))
end do
end do
Expand All @@ -465,7 +462,7 @@ subroutine fielddump

if (lnetcdf .and. ltntr) then
do k = klow,khigh
vars(:,:,k,ind_tntr) = ( & !! need to loop over k here
vars(:,:,k-klow+1,ind_tntr) = ( & !! need to loop over k here
-swd(2:i1:ncoarse,2:j1:ncoarse,k+1) &
-swu(2:i1:ncoarse,2:j1:ncoarse,k+1) &
+swd(2:i1:ncoarse,2:j1:ncoarse,k) &
Expand All @@ -480,7 +477,7 @@ subroutine fielddump

if (lnetcdf .and. ltntrs) then
do k = klow,khigh
vars(:,:,k,ind_tntrs) = ( & !! need to loop over k here
vars(:,:,k-klow+1,ind_tntrs) = ( & !! need to loop over k here
-swd(2:i1:ncoarse,2:j1:ncoarse,k+1) &
-swu(2:i1:ncoarse,2:j1:ncoarse,k+1) &
+swd(2:i1:ncoarse,2:j1:ncoarse,k) &
Expand All @@ -491,7 +488,7 @@ subroutine fielddump

if (lnetcdf .and. ltntrl) then
do k = klow,khigh
vars(:,:,k,ind_tntrl) = ( & !! need to loop over k here
vars(:,:,k-klow+1,ind_tntrl) = ( & !! need to loop over k here
-lwd(2:i1:ncoarse,2:j1:ncoarse,k+1) &
-lwu(2:i1:ncoarse,2:j1:ncoarse,k+1) &
+lwd(2:i1:ncoarse,2:j1:ncoarse,k) &
Expand Down
25 changes: 16 additions & 9 deletions src/modstat_nc.f90
Original file line number Diff line number Diff line change
Expand Up @@ -34,7 +34,7 @@ module modstat_nc
logical :: lsync = .false. ! Sync NetCDF file after each writestat_*_nc
logical :: lclassic = .false. ! Create netCDF in CLASSIC format (less RAM usage, compression not supported)
integer :: deflate = 2 ! Deflate level for netCDF files (only for NETCDF4 format)

integer, save :: timeID=0, ztID=0, zmID=0, xtID=0, xmID=0, ytID=0, ymID=0,ztsID=0, zqID=0
real(kind=4) :: nc_fillvalue = -999.
!> The only interface necessary to write data to netcdf, regardless of the dimensions.
Expand Down Expand Up @@ -70,7 +70,7 @@ subroutine initstat_nc
call D_MPI_BCAST(lsync ,1, 0,comm3d,mpierr)
call D_MPI_BCAST(lclassic ,1, 0,comm3d,mpierr)
call D_MPI_BCAST(deflate ,1, 0,comm3d,mpierr)

end subroutine initstat_nc
!
! ----------------------------------------------------------------------
Expand Down Expand Up @@ -192,7 +192,7 @@ subroutine open_nc (fname, ncid,nrec,n1, n2, n3, ns,nq)
iret = nf90_enddef(ncID)
! Fails with NetCDF: Operation not allowed in data mode
! when the netCDF files already exist

end subroutine open_nc

!
Expand Down Expand Up @@ -339,20 +339,27 @@ subroutine exitstat_nc(ncid)
call nchandle_error(status)
end subroutine exitstat_nc

subroutine writestat_dims_nc(ncid, ncoarse)
use modglobal, only : dx,dy,zf,zh,jmax,imax
subroutine writestat_dims_nc(ncid, ncoarse, klow)
! optional arguments ncoarse (coarsegraining in the horizontal directions)
! klow (lower bound for z. Upper bound is taken from the size of the dimension)
use modglobal, only : dx,dy,zf,zh,jmax,imax,kmax
use modsurfdata, only : zsoilc,isurf
use modmpi, only : myidx,myidy
implicit none
integer, intent(in) :: ncid
integer, optional, intent(in) :: ncoarse
integer, optional, intent(in) :: ncoarse, klow
integer :: i=0,iret,length,varid, nc

integer :: kl
if (present(ncoarse)) then
nc = ncoarse
else
nc = 1
end if
if (present(klow)) then
kl = klow
else
kl = 1
end if

iret = nf90_inq_varid(ncid, 'xt', VarID)
if (iret==0) iret=nf90_inquire_dimension(ncid, xtID, len=length)
Expand All @@ -370,10 +377,10 @@ subroutine writestat_dims_nc(ncid, ncoarse)

iret = nf90_inq_varid(ncid, 'zt', VarID)
if (iret==0) iret=nf90_inquire_dimension(ncid,ztID, len=length)
if (iret==0) iret = nf90_put_var(ncid, varID, zf(1:length),(/1/))
if (iret==0) iret = nf90_put_var(ncid, varID, zf(kl:kl+length-1),(/1/))
iret = nf90_inq_varid(ncid, 'zm', VarID)
if (iret==0) iret=nf90_inquire_dimension(ncid, zmID, len=length)
if (iret==0) iret = nf90_put_var(ncid, varID, zh(1:length),(/1/))
if (iret==0) iret = nf90_put_var(ncid, varID, zh(kl:kl+length-1),(/1/))
if (isurf==1) then
iret = nf90_inq_varid(ncid, 'zts', VarID)
if (iret==0) iret = nf90_inquire_dimension(ncid, ztsID, len=length)
Expand Down

0 comments on commit d35a03e

Please sign in to comment.