MODULE wrt_state_mod ! !================================================== 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 requested tangent linear or adjoint state ! ! fields using the standard NetCDF library or the Parallel-IO (PIO) ! ! library into an output file described by IO structure S. The NetCDF ! ! is opened and closed, and file ID or descriptor is not updated. ! ! ! ! On Input: ! ! ! ! ng Nested grid number (integer) ! ! model Calling model identifier (integer) ! ! label Identification label (string) ! ! OutRec NetCDF file time record (integer) ! ! i2d 2D state variables time level index (integer) ! ! i3d 3D state variables time level index (integer) ! ! stime Time variable value (real; seconds) ! ! S File derived type structure (TYPE T_IO) ! ! CloseFile Switch to close NetCDF file (logical, OPTIONAL) ! ! (defaut = .TRUE.) ! ! On Output: ! ! ! ! ad_arrays ADM state arrays if model = iADM ! ! or, ! ! tl_arrays TLM state arrays if model = iTLM or iRPM ! ! ! !======================================================================= ! USE mod_param USE mod_parallel USE mod_coupling 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_state PRIVATE :: wrt_state_nf90 ! CONTAINS ! !*********************************************************************** SUBROUTINE wrt_state (ng, tile, model, label, & & OutRec, i2d, i3d, stime, S) !*********************************************************************** ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model, OutRec, i2d, i3d ! real(dp), intent(in) :: stime ! character (len=*), intent(in) :: label ! TYPE(T_IO), intent(inout) :: S(Ngrids) ! ! Local variable declarations. ! integer :: LBi, UBi, LBj, UBj ! character (len=*), parameter :: MyFile = & & "ROMS/Utility/wrt_state.F" ! !----------------------------------------------------------------------- ! Write out history fields according to IO type. !----------------------------------------------------------------------- ! LBi=BOUNDS(ng)%LBi(tile) UBi=BOUNDS(ng)%UBi(tile) LBj=BOUNDS(ng)%LBj(tile) UBj=BOUNDS(ng)%UBj(tile) ! SELECT CASE (HIS(ng)%IOtype) CASE (io_nf90) CALL wrt_state_nf90 (ng, tile, model, label, & & OutRec, i2d, i3d, stime, & & S) CASE DEFAULT IF (Master) THEN WRITE (stdout,10) S(ng)%IOtype END IF exit_flag=3 END SELECT IF (FoundError(exit_flag, NoError, 144, MyFile)) RETURN ! 10 FORMAT (' WRT_STATE - Illegal output file type, io_type = ',i0, & & /,13x,'Check KeyWord ''OUT_LIB'' in ''roms.in''.') ! RETURN END SUBROUTINE wrt_state ! !*********************************************************************** SUBROUTINE wrt_state_nf90 (ng, tile, model, label, & & OutRec, i2d, i3d, stime, & & S) !*********************************************************************** ! USE mod_netcdf ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model, OutRec, i2d, i3d integer :: LBi, UBi, LBj, UBj ! real(dp), intent(in) :: stime ! character (len=*), intent(in) :: label ! TYPE(T_IO), intent(inout) :: S(Ngrids) ! ! Local variable declarations. ! integer :: Sstr, Send integer :: Fcount, gfactor, gtype, ncid, omode, status integer :: i, itrc, j, k ! real(r8) :: Fmin, Fmax real(dp) :: scale real(dp) :: Tval(1) ! character (len=15) :: Tstring character (len=22) :: t_code character (len=50) :: string character (len=*), parameter :: MyFile = & & "ROMS/Utility/wrt_state.F"//", wrt_state_nf90" ! SourceFile=MyFile ! !----------------------------------------------------------------------- ! Write out tangent linear or adjoint state fields. !----------------------------------------------------------------------- ! ! Turn on time wall clock. ! CALL wclock_on (ng, model, 81, 206, MyFile) ! ! Open output NetCDF file. ! omode=1 ! read and write access CALL netcdf_open (ng, model, S(ng)%name, omode, ncid) IF (FoundError(exit_flag, NoError, 213, 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. ! S(ng)%Rindex=S(ng)%Rindex+1 Fcount=S(ng)%Fcount S(ng)%Nrec(Fcount)=S(ng)%Nrec(Fcount)+1 ! ! Report information. ! IF (Master) THEN CALL time_string (stime, t_code) IF (model.eq.iNLM) THEN string='writing NLM state fields,' ELSE IF (model.eq.iTLM) THEN string='writing TLM state fields,' ELSE IF (model.eq.iADM) THEN string='writing ADM state fields,' ELSE IF (model.eq.iRPM) THEN string='writing RPM state fields,' END IF Sstr=SCAN(CalledFrom,'/',BACK=.TRUE.)+1 Send=LEN_TRIM(CalledFrom) WRITE (Tstring,'(f15.4)') stime*sec2day WRITE (stdout,10) TRIM(label), TRIM(string), t_code, & & ng, TRIM(ADJUSTL(Tstring)), TRIM(S(ng)%name), & & i3d, OutRec, CalledFrom(Sstr:Send) END IF ! ! Write out model time (s). ! IF (LwrtPER(ng)) THEN Tval(1)=REAL(OutRec,dp)*day2sec ELSE Tval(1)=stime END IF CALL netcdf_put_fvar (ng, model, S(ng)%name, & & TRIM(Vname(1,idtime)), tval, & & (/OutRec/), (/1/), & & ncid = ncid, & & varid = S(ng)%Vid(idtime)) IF (FoundError(exit_flag, NoError, 263, MyFile)) RETURN ! ! Write out free-surface. ! scale=1.0_dp gtype=gfactor*r2dvar IF (model.eq.iTLM) THEN status=nf_fwrite2d(ng, model, ncid, idFsur, & & S(ng)%Vid(idFsur), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % rmask, & & OCEAN(ng) % tl_zeta(:,:,i2d), & & MinValue = Fmin, & & MaxValue = Fmax) ELSE IF (model.eq.iADM) THEN status=nf_fwrite2d(ng, model, ncid, idFsur, & & S(ng)%Vid(idFsur), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % rmask, & & OCEAN(ng) % ad_zeta(:,:,i2d), & & MinValue = Fmin, & & MaxValue = Fmax) END IF IF (status.ne.nf90_noerr) THEN IF (Master) THEN WRITE (stdout,20) TRIM(Vname(1,idFsur)), TRIM(label), & & TRIM(S(ng)%name), OutRec 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 U-momentum component (m/s). ! scale=1.0_dp gtype=gfactor*u2dvar IF (model.eq.iTLM) THEN status=nf_fwrite2d(ng, model, ncid, idUbar, & & S(ng)%Vid(idUbar), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % umask, & & OCEAN(ng) % tl_ubar(:,:,i2d), & & MinValue = Fmin, & & MaxValue = Fmax) ELSE IF (model.eq.iADM) THEN status=nf_fwrite2d(ng, model, ncid, idUbar, & & S(ng)%Vid(idUbar), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % umask, & & OCEAN(ng) % ad_ubar(:,:,i2d), & & MinValue = Fmin, & & MaxValue = Fmax) END IF IF (status.ne.nf90_noerr) THEN IF (Master) THEN WRITE (stdout,20) TRIM(Vname(1,idUbar)), TRIM(label), & & TRIM(S(ng)%name), OutRec 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 V-momentum component. ! scale=1.0_dp gtype=gfactor*v2dvar IF (model.eq.iTLM) THEN status=nf_fwrite2d(ng, model, ncid, idVbar, & & S(ng)%Vid(idVbar), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % vmask, & & OCEAN(ng) % tl_vbar(:,:,i2d), & & MinValue = Fmin, & & MaxValue = Fmax) ELSE IF (model.eq.iADM) THEN status=nf_fwrite2d(ng, model, ncid, idVbar, & & S(ng)%Vid(idVbar), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % vmask, & & OCEAN(ng) % ad_vbar(:,:,i2d), & & MinValue = Fmin, & & MaxValue = Fmax) END IF IF (status.ne.nf90_noerr) THEN IF (Master) THEN WRITE (stdout,20) TRIM(Vname(1,idVbar)), TRIM(label), & & TRIM(S(ng)%name), OutRec 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 U-momentum component (m/s). ! scale=1.0_dp gtype=gfactor*u3dvar IF (model.eq.iTLM) THEN status=nf_fwrite3d(ng, model, ncid, idUvel, & & S(ng)%Vid(idUvel), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % umask, & & OCEAN(ng) % tl_u(:,:,:,i3d), & & MinValue = Fmin, & & MaxValue = Fmax) ELSE IF (model.eq.iADM) THEN status=nf_fwrite3d(ng, model, ncid, idUvel, & & S(ng)%Vid(idUvel), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % umask, & & OCEAN(ng) % ad_u(:,:,:,i3d), & & MinValue = Fmin, & & MaxValue = Fmax) END IF IF (status.ne.nf90_noerr) THEN IF (Master) THEN WRITE (stdout,20) TRIM(Vname(1,idUvel)), TRIM(label), & & TRIM(S(ng)%name), OutRec 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 V-momentum component (m/s). ! scale=1.0_dp gtype=gfactor*v3dvar IF (model.eq.iTLM) THEN status=nf_fwrite3d(ng, model, ncid, idVvel, & & S(ng)%Vid(idVvel), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % vmask, & & OCEAN(ng) % tl_v(:,:,:,i3d), & & MinValue = Fmin, & & MaxValue = Fmax) ELSE IF (model.eq.iADM) THEN status=nf_fwrite3d(ng, model, ncid, idVvel, & & S(ng)%Vid(idVvel), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % vmask, & & OCEAN(ng) % ad_v(:,:,:,i3d), & & MinValue = Fmin, & & MaxValue = Fmax) END IF IF (status.ne.nf90_noerr) THEN IF (Master) THEN WRITE (stdout,20) TRIM(Vname(1,idVvel)), TRIM(label), & & TRIM(S(ng)%name), OutRec 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 IF (model.eq.iTLM) THEN status=nf_fwrite3d(ng, model, ncid, idTvar(itrc), & & S(ng)%Tid(itrc), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % rmask, & & OCEAN(ng) % tl_t(:,:,:,i3d,itrc), & & MinValue = Fmin, & & MaxValue = Fmax) ELSE IF (model.eq.iADM) THEN status=nf_fwrite3d(ng, model, ncid, idTvar(itrc), & & S(ng)%Tid(itrc), & & OutRec, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % rmask, & & OCEAN(ng) % ad_t(:,:,:,i3d,itrc), & & MinValue = Fmin, & & MaxValue = Fmax) END IF IF (status.ne.nf90_noerr) THEN IF (Master) THEN WRITE (stdout,20) TRIM(Vname(1,idTvar(itrc))), & & TRIM(label), TRIM(S(ng)%name), OutRec 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 ! !----------------------------------------------------------------------- ! Synchronize tangent NetCDF file to disk to allow other processes ! to access data immediately after it is written. Then, close NetCDF ! file. !----------------------------------------------------------------------- ! CALL netcdf_sync (ng, iTLM, S(ng)%name, ncid) IF (FoundError(exit_flag, NoError, 950, MyFile)) RETURN ! CALL netcdf_close (ng, iTLM, ncid, S(ng)%name, .FALSE.) ! ! Turn off time wall clock. ! CALL wclock_off (ng, model, 81, 958, MyFile) ! 10 FORMAT (1x,a,': WRT_STATE_NF90 - ',a,t75,a, & & /,24x,'(Grid ',i2.2,', t = ',a,', File: ',a, & & ', Index=',i1,', Rec=',i0,')', & & /,24x,'Called from ''',a,'''') 20 FORMAT (/,' WRT_STATE_NF90 - error while writing variable: ',a, & & /,18x,'into ',a,' NetCDF file: ',a, & & /,18x,'for time record: ',i0) 30 FORMAT (21x,'- ',a,/,24x,'(Min = ',1p,e15.8, & & ' Max = ',1p,e15.8,')') ! RETURN END SUBROUTINE wrt_state_nf90 END MODULE wrt_state_mod