MODULE nf_fwrite3d_mod ! !git $Id$ !svn $Id: nf_fwrite3d.F 1190 2023-08-18 19:51:09Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! This function writes out a generic floating point 3D array into an ! ! output NetCDF file. ! ! ! ! On Input: ! ! ! ! ng Nested grid number ! ! model Calling model identifier ! ! ncid NetCDF file ID ! ! ifield Field metadata index (integer) ! ! ncvarid NetCDF variable ID ! ! tindex NetCDF time record index to write ! ! gtype Grid type. If negative, only write water points ! ! LBi I-dimension Lower bound ! ! UBi I-dimension Upper bound ! ! LBj J-dimension Lower bound ! ! UBj J-dimension Upper bound ! ! LBk K-dimension Lower bound ! ! UBk K-dimension Upper bound ! ! Amask land/Sea mask, if any (real) ! ! Ascl Factor to scale field before writing (real) ! ! Adat Field to write out (real) ! ! SetFillVal Logical switch to set fill value in land areas ! ! (OPTIONAL) ! ! ! ! On Output: ! ! ! ! status Error flag (integer) ! ! MinValue Minimum value (real, OPTIONAL) ! ! MaxValue Maximum value (real, OPTIONAL) ! ! ! !======================================================================= ! USE mod_param USE mod_parallel USE mod_ncparam USE mod_scalars ! implicit none ! INTERFACE nf_fwrite3d MODULE PROCEDURE nf90_fwrite3d END INTERFACE nf_fwrite3d ! CONTAINS ! !*********************************************************************** FUNCTION nf90_fwrite3d (ng, model, ncid, ifield, & & ncvarid, tindex, gtype, & & LBi, UBi, LBj, UBj, LBk, UBk, Ascl, & & Amask, & & Adat, SetFillVal, & & MinValue, MaxValue) RESULT (status) !*********************************************************************** ! USE mod_netcdf ! USE distribute_mod, ONLY : mp_bcasti USE distribute_mod, ONLY : mp_gather3d ! ! Imported variable declarations. ! logical, intent(in), optional :: SetFillVal ! integer, intent(in) :: ng, model, ncid, ncvarid, tindex, gtype integer, intent(in) :: ifield integer, intent(in) :: LBi, UBi, LBj, UBj, LBk, UBk ! real(dp), intent(in) :: Ascl ! real(r8), intent(in) :: Amask(LBi:,LBj:) real(r8), intent(in) :: Adat(LBi:,LBj:,LBk:) real(r8), intent(out), optional :: MinValue real(r8), intent(out), optional :: MaxValue ! ! Local variable declarations. ! logical :: LandFill integer :: i, j, k, ic, Npts, tile integer :: Imin, Imax, Jmin, Jmax, Koff integer :: Ilen, Jlen, Klen, IJlen, MyType integer :: status integer, dimension(4) :: start, total ! real(r8), dimension((Lm(ng)+2)*(Mm(ng)+2)*(UBk-LBk+1)) :: Awrk ! !----------------------------------------------------------------------- ! Set starting and ending indices to process. !----------------------------------------------------------------------- ! status=nf90_noerr ! ! Set first and last grid point according to staggered C-grid ! classification. Set loops offsets. ! MyType=gtype ! SELECT CASE (ABS(MyType)) CASE (p3dvar) Imin=IOBOUNDS(ng)%ILB_psi Imax=IOBOUNDS(ng)%IUB_psi Jmin=IOBOUNDS(ng)%JLB_psi Jmax=IOBOUNDS(ng)%JUB_psi CASE (r3dvar) Imin=IOBOUNDS(ng)%ILB_rho Imax=IOBOUNDS(ng)%IUB_rho Jmin=IOBOUNDS(ng)%JLB_rho Jmax=IOBOUNDS(ng)%JUB_rho CASE (u3dvar) Imin=IOBOUNDS(ng)%ILB_u Imax=IOBOUNDS(ng)%IUB_u Jmin=IOBOUNDS(ng)%JLB_u Jmax=IOBOUNDS(ng)%JUB_u CASE (v3dvar) Imin=IOBOUNDS(ng)%ILB_v Imax=IOBOUNDS(ng)%IUB_v Jmin=IOBOUNDS(ng)%JLB_v Jmax=IOBOUNDS(ng)%JUB_v CASE DEFAULT Imin=IOBOUNDS(ng)%ILB_rho Imax=IOBOUNDS(ng)%IUB_rho Jmin=IOBOUNDS(ng)%JLB_rho Jmax=IOBOUNDS(ng)%JUB_rho END SELECT ! Ilen=Imax-Imin+1 Jlen=Jmax-Jmin+1 Klen=UBk-LBk+1 IJlen=Ilen*Jlen ! IF (LBk.eq.0) THEN Koff=0 ELSE Koff=1 END IF ! ! Set switch to replace land areas with fill value, spval. ! IF (PRESENT(SetFillVal)) THEN LandFill=SetFillVal ELSE LandFill=tindex.gt.0 END IF ! ! If appropriate, initialize minimum and maximum processed values. ! IF (PRESENT(MinValue)) THEN MinValue=spval MaxValue=-spval END IF ! ! Initialize local array to avoid denormalized numbers. This ! facilitates processing and debugging. ! Awrk=0.0_r8 ! !----------------------------------------------------------------------- ! If distributed-memory set-up, collect tile data from all spawned ! nodes and store it into a global scratch 1D array, packed in column- ! major order. ! Overwrite masked points with special value. !----------------------------------------------------------------------- ! CALL mp_gather3d (ng, model, LBi, UBi, LBj, UBj, LBk, UBk, & & tindex, gtype, Ascl, & & Amask, & & Adat, Npts, Awrk, SetFillVal) ! !----------------------------------------------------------------------- ! If applicable, compute output field minimum and maximum values. !----------------------------------------------------------------------- ! IF (PRESENT(MinValue)) THEN IF (OutThread) THEN DO i=1,Npts IF (ABS(Awrk(i)).lt.spval) THEN MinValue=MIN(MinValue,Awrk(i)) MaxValue=MAX(MaxValue,Awrk(i)) END IF END DO END IF END IF ! !----------------------------------------------------------------------- ! Write output buffer into NetCDF file. !----------------------------------------------------------------------- ! IF (OutThread) THEN IF (gtype.gt.0) THEN start(1)=1 total(1)=Ilen start(2)=1 total(2)=Jlen start(3)=1 total(3)=Klen start(4)=tindex total(4)=1 ELSE start(1)=1 total(1)=Npts start(2)=tindex total(2)=1 END IF status=nf90_put_var(ncid, ncvarid, Awrk, start, total) END IF ! !----------------------------------------------------------------------- ! Broadcast IO error flag to all nodes. !----------------------------------------------------------------------- ! CALL mp_bcasti (ng, model, status) ! RETURN END FUNCTION nf90_fwrite3d ! END MODULE nf_fwrite3d_mod