MODULE wrt_ini_mod ! !git $Id$ !svn $Id: wrt_ini.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 module writes state variables initial conditions into initial ! ! file using either the standard NetCDF library or the Parallel-IO ! ! (PIO) library. ! ! ! ! On Input: ! ! ! ! ng Nested grid number. ! ! Tindex State variables time index to write. ! ! OutRec NetCDF file unlimited dimension record to write. ! ! ! ! Notice that only momentum is affected by the full time-averaged ! ! masks. If applicable, these mask contains information about ! ! river runoff and time-dependent wetting and drying variations. ! ! ! !======================================================================= ! USE mod_param USE mod_parallel USE mod_fourdvar USE mod_grid USE mod_iounits USE mod_mixing USE mod_ncparam USE mod_ocean USE mod_scalars USE mod_stepping ! USE dateclock_mod, ONLY : time_string USE nf_fwrite2d_mod, ONLY : nf_fwrite2d USE nf_fwrite3d_mod, ONLY : nf_fwrite3d USE strings_mod, ONLY : FoundError ! implicit none ! PUBLIC :: wrt_ini PRIVATE :: wrt_ini_nf90 ! CONTAINS ! !*********************************************************************** SUBROUTINE wrt_ini (ng, tile, Tindex, OutRec) !*********************************************************************** ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, Tindex integer, intent(in), optional :: OutRec ! ! Local variable declarations. ! integer :: LBi, UBi, LBj, UBj ! character (len=*), parameter :: MyFile = & & "ROMS/Utility/wrt_ini.F" ! !----------------------------------------------------------------------- ! Write out nonlinear initial conditions. !----------------------------------------------------------------------- ! LBi=BOUNDS(ng)%LBi(tile) UBi=BOUNDS(ng)%UBi(tile) LBj=BOUNDS(ng)%LBj(tile) UBj=BOUNDS(ng)%UBj(tile) ! SELECT CASE (DIA(ng)%IOtype) CASE (io_nf90) CALL wrt_ini_nf90 (ng, tile, Tindex, & & LBi, UBi, LBj, UBj, & & OutRec) CASE DEFAULT IF (Master) WRITE (stdout,10) DAI(ng)%IOtype exit_flag=3 END SELECT IF (FoundError(exit_flag, NoError, 129, MyFile)) RETURN ! 10 FORMAT (' WRT_INI - Illegal output file type, io_type = ',i0, & & /,11x,'Check KeyWord ''OUT_LIB'' in ''roms.in''.') ! RETURN END SUBROUTINE wrt_ini ! !*********************************************************************** SUBROUTINE wrt_ini_nf90 (ng, tile, Tindex, & & LBi, UBi, LBj, UBj, & & OutRec) !*********************************************************************** ! USE mod_netcdf ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, Tindex integer, intent(in) :: LBi, UBi, LBj, UBj ! integer, intent(in), optional :: OutRec ! ! Local variable declarations. ! logical :: SetFillVal ! integer :: Fcount, gfactor, gtype, i, itrc, status, varid ! real(r8) :: Fmin, Fmax real(dp) :: my_time, scale ! character (len=15) :: Tstring character (len=22) :: t_code character (len=*), parameter :: MyFile = & & "ROMS/Utility/wrt_ini.F"//", wrt_ini_nf90" ! SourceFile=MyFile ! !----------------------------------------------------------------------- ! Set SetFillVal to FALSE. This is essential to ensure identical ! outer-loop solutions when increments are zero. NOT SURE WHY THIS ! IS NECESSARY - PERHAPS SOMETHING WRONG WITH THE LOGIC OF ! nf_fwrite routines. !----------------------------------------------------------------------- ! SetFillVal=.FALSE. ! !----------------------------------------------------------------------- ! Write out initial conditions. !----------------------------------------------------------------------- ! IF (FoundError(exit_flag, NoError, 182, MyFile)) RETURN ! ! Set grid type factor to write full (gfactor=1) fields or water ! points (gfactor=-1) fields only. ! gfactor=1 ! ! Set time record index. ! IF (PRESENT(OutRec)) THEN INI(ng)%Rindex=OutRec ELSE INI(ng)%Rindex=INI(ng)%Rindex+1 END IF Fcount=INI(ng)%Fcount INI(ng)%Nrec(Fcount)=INI(ng)%Nrec(Fcount)+1 ! ! Write out model time (s). Use the "tdays" variable here because of ! the management of the "time" variable due to nesting. ! my_time=tdays(ng)*day2sec CALL netcdf_put_fvar (ng, iNLM, INI(ng)%name, & & TRIM(Vname(1,idtime)), my_time, & & (/INI(ng)%Rindex/), (/1/), & & ncid = INI(ng)%ncid, & & varid = INI(ng)%Vid(idtime)) IF (FoundError(exit_flag, NoError, 213, MyFile)) THEN IF (Master) THEN WRITE (stdout,10) TRIM(Vname(1,idFsur)), TRIM(INI(ng)%name), & & INI(ng)%Rindex END IF RETURN ELSE IF (Master) THEN CALL time_string (my_time, t_code) WRITE (Tstring,'(f15.4)') my_time*sec2day WRITE (stdout,20) t_code, ng, TRIM(ADJUSTL(Tstring)), & & TRIM(INI(ng)%name), INI(ng)%Rindex, Tindex END IF END IF ! ! Write out free-surface (m) ! scale=1.0_dp gtype=gfactor*r2dvar status=nf_fwrite2d(ng, iNLM, INI(ng)%ncid, idFsur, & & INI(ng)%Vid(idFsur), & & INI(ng)%Rindex, gtype, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % rmask, & & OCEAN(ng) % zeta(:,:,Tindex), & & SetFillVal = SetFillVal, & & MinValue = Fmin, & & MaxValue = Fmax) IF (FoundError(status, nf90_noerr, 247, MyFile)) THEN IF (Master) THEN WRITE (stdout,10) TRIM(Vname(1,idFsur)), TRIM(INI(ng)%name), & & INI(ng)%Rindex END IF exit_flag=3 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,30) TRIM(Vname(2,idFsur)), Fmin, Fmax END IF END IF ! ! Write out 2D momentum component (m/s) in the XI-direction. ! scale=1.0_dp gtype=gfactor*u2dvar status=nf_fwrite2d(ng, iNLM, INI(ng)%ncid, idUbar, & & INI(ng)%Vid(idUbar), & & INI(ng)%Rindex, gtype, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % umask_full, & & OCEAN(ng) % ubar(:,:,Tindex), & & SetFillVal = SetFillVal, & & MinValue = Fmin, & & MaxValue = Fmax) IF (FoundError(status, nf90_noerr, 276, MyFile)) THEN IF (Master) THEN WRITE (stdout,10) TRIM(Vname(1,idUbar)), TRIM(INI(ng)%name), & & INI(ng)%Rindex END IF exit_flag=3 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,30) TRIM(Vname(2,idUbar)), Fmin, Fmax END IF END IF ! ! Write out 2D momentum component (m/s) in the ETA-direction. ! scale=1.0_dp gtype=gfactor*v2dvar status=nf_fwrite2d(ng, iNLM, INI(ng)%ncid, idVbar, & & INI(ng)%Vid(idVbar), & & INI(ng)%Rindex, gtype, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % vmask_full, & & OCEAN(ng) % vbar(:,:,Tindex), & & SetFillVal = SetFillVal, & & MinValue = Fmin, & & MaxValue = Fmax) IF (FoundError(status, nf90_noerr, 305, MyFile)) THEN IF (Master) THEN WRITE (stdout,10) TRIM(Vname(1,idVbar)), TRIM(INI(ng)%name), & & INI(ng)%Rindex END IF exit_flag=3 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,30) TRIM(Vname(2,idVbar)), Fmin, Fmax END IF END IF ! ! Write out 3D momentum component (m/s) in the XI-direction. ! scale=1.0_dp gtype=gfactor*u3dvar status=nf_fwrite3d(ng, iNLM, INI(ng)%ncid, idUvel, & & INI(ng)%Vid(idUvel), & & INI(ng)%Rindex, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % umask_full, & & OCEAN(ng) % u(:,:,:,Tindex), & & SetFillVal = SetFillVal, & & MinValue = Fmin, & & MaxValue = Fmax) IF (FoundError(status, nf90_noerr, 336, MyFile)) THEN IF (Master) THEN WRITE (stdout,10) TRIM(Vname(1,idUvel)), TRIM(INI(ng)%name), & & INI(ng)%Rindex END IF exit_flag=3 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,30) TRIM(Vname(2,idUvel)), Fmin, Fmax END IF END IF ! ! Write out 3D momentum component (m/s) in the ETA-direction. ! scale=1.0_dp gtype=gfactor*v3dvar status=nf_fwrite3d(ng, iNLM, INI(ng)%ncid, idVvel, & & INI(ng)%Vid(idVvel), & & INI(ng)%Rindex, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % vmask_full, & & OCEAN(ng) % v(:,:,:,Tindex), & & SetFillVal = SetFillVal, & & MinValue = Fmin, & & MaxValue = Fmax) IF (FoundError(status, nf90_noerr, 365, MyFile)) THEN IF (Master) THEN WRITE (stdout,10) TRIM(Vname(1,idVvel)), TRIM(INI(ng)%name), & & INI(ng)%Rindex END IF exit_flag=3 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,30) TRIM(Vname(2,idVvel)), Fmin, Fmax END IF END IF ! ! Write out tracer type variables. ! DO itrc=1,NT(ng) scale=1.0_dp gtype=gfactor*r3dvar status=nf_fwrite3d(ng, iNLM, INI(ng)%ncid, idTvar(itrc), & & INI(ng)%Tid(itrc), & & INI(ng)%Rindex, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % rmask, & & OCEAN(ng) % t(:,:,:,Tindex,itrc), & & SetFillVal = SetFillVal, & & MinValue = Fmin, & & MaxValue = Fmax) IF (FoundError(status, nf90_noerr, 395, MyFile)) THEN IF (Master) THEN WRITE (stdout,10) TRIM(Vname(1,idTvar(itrc))), & & TRIM(INI(ng)%name), INI(ng)%Rindex END IF exit_flag=3 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,30) TRIM(Vname(2,idTvar(itrc))), Fmin, Fmax END IF END IF END DO ! ! If defined, write out vertical viscosity coefficient. ! IF (INI(ng)%Vid(idVvis).le.0) THEN CALL netcdf_inq_varid (ng, iNLM, INI(ng)%name, Vname(1,idVvis), & & INI(ng)%ncid, INI(ng)%Vid(idVvis)) IF (FoundError(status, nf90_noerr, 418, MyFile)) RETURN END IF scale=1.0_dp gtype=gfactor*w3dvar status=nf_fwrite3d(ng, iNLM, INI(ng)%ncid, idVvis, & & INI(ng)%Vid(idVvis), & & INI(ng)%Rindex, gtype, & & LBi, UBi, LBj, UBj, 0, N(ng), scale, & & GRID(ng) % rmask, & & MIXING(ng) % Akv, & & SetFillVal = .FALSE., & & MinValue = Fmin, & & MaxValue = Fmax) IF (FoundError(status, nf90_noerr, 433, MyFile)) THEN IF (Master) THEN WRITE (stdout,10) TRIM(Vname(1,idVvis)), & & TRIM(INI(ng)%name), INI(ng)%Rindex END IF exit_flag=3 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,30) TRIM(Vname(2,idVvis)), Fmin, Fmax END IF END IF ! ! If defined, write out vertical diffusion coefficient for potential ! temperature. ! IF (INI(ng)%Vid(idTdif).le.0) THEN CALL netcdf_inq_varid (ng, iNLM, INI(ng)%name, Vname(1,idTdif), & & INI(ng)%ncid, INI(ng)%Vid(idTdif)) IF (FoundError(status, nf90_noerr, 453, MyFile)) RETURN END IF scale=1.0_dp gtype=gfactor*w3dvar status=nf_fwrite3d(ng, iNLM, INI(ng)%ncid, idTdif, & & INI(ng)%Vid(idTdif), & & INI(ng)%Rindex, gtype, & & LBi, UBi, LBj, UBj, 0, N(ng), scale, & & GRID(ng) % rmask, & & MIXING(ng) % Akt(:,:,:,itemp), & & SetFillVal = .FALSE., & & MinValue = Fmin, & & MaxValue = Fmax) IF (FoundError(status, nf90_noerr, 468, MyFile)) THEN IF (Master) THEN WRITE (stdout,10) TRIM(Vname(1,idTdif)), & & TRIM(INI(ng)%name), INI(ng)%Rindex END IF exit_flag=3 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,30) TRIM(Vname(2,idTdif)), Fmin, Fmax END IF END IF ! ! If defined, write out vertical diffusion coefficient for salinity. ! IF (INI(ng)%Vid(idSdif).le.0) THEN CALL netcdf_inq_varid (ng, iNLM, INI(ng)%name, Vname(1,idSdif), & & INI(ng)%ncid, INI(ng)%Vid(idSdif)) IF (FoundError(status, nf90_noerr, 489, MyFile)) RETURN END IF scale=1.0_dp gtype=gfactor*w3dvar status=nf_fwrite3d(ng, iNLM, INI(ng)%ncid, idSdif, & & INI(ng)%Vid(idSdif), & & INI(ng)%Rindex, gtype, & & LBi, UBi, LBj, UBj, 0, N(ng), scale, & & GRID(ng) % rmask, & & MIXING(ng) % Akt(:,:,:,isalt), & & SetFillVal = .FALSE., & & MinValue = Fmin, & & MaxValue = Fmax) IF (FoundError(status, nf90_noerr, 504, MyFile)) THEN IF (Master) THEN WRITE (stdout,10) TRIM(Vname(1,idSdif)), & & TRIM(INI(ng)%name), INI(ng)%Rindex END IF exit_flag=3 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,30) TRIM(Vname(2,idSdif)), Fmin, Fmax END IF END IF ! !----------------------------------------------------------------------- ! Synchronize initial NetCDF file to disk to allow other processes ! to access data immediately after it is written. !----------------------------------------------------------------------- ! CALL netcdf_sync (ng, iNLM, INI(ng)%name, INI(ng)%ncid) IF (FoundError(exit_flag, NoError, 527, MyFile)) RETURN ! 10 FORMAT (/,' WRT_INI_NF90 - error while writing variable: ',a, & & /,16x,'into initial NetCDF file for time record: ',i0) 20 FORMAT (2x,'WRT_INI_NF90 - NLM: Writing initial state', & & ' fields,',t75,a, & & /,22x,'(Grid ',i2.2,', t = ',a,', File: ',a, & & ', Rec=',i4.4,', Index=',i1,')') 30 FORMAT (19x,'- ',a,/,22x,'(Min = ',1p,e15.8, & & ' Max = ',1p,e15.8,')') ! RETURN END SUBROUTINE wrt_ini_nf90 END MODULE wrt_ini_mod