MODULE wrt_aug_imp_mod ! !git $Id$ !svn $Id: wrt_aug_imp.F 1189 2023-08-15 21:26:58Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! This routine writes augmented impulse state variables using the ! ! standard NetCDF library or the Parallel-IO (PIO) library. ! ! ! !======================================================================= ! USE mod_param USE mod_parallel USE mod_grid USE mod_iounits USE mod_ncparam USE mod_ocean USE mod_scalars ! USE nf_fwrite2d_mod, ONLY : nf_fwrite2d USE nf_fwrite3d_mod, ONLY : nf_fwrite3d USE strings_mod, ONLY : FoundError, find_string ! implicit none ! PUBLIC :: wrt_aug_imp PRIVATE :: wrt_aug_imp_nf90 ! CONTAINS ! !*********************************************************************** SUBROUTINE wrt_aug_imp (ng, tile, model, Iinp, Iout) !*********************************************************************** ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model, Iinp, Iout ! ! Local variable declarations. ! integer :: LBi, UBi, LBj, UBj ! character (len=*), parameter :: MyFile = & & "ROMS/Utility/wrt_aug_imp.F" ! !----------------------------------------------------------------------- ! Write out time-averaged 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 (TLF(ng)%IOtype) CASE (io_nf90) CALL wrt_aug_imp_nf90 (ng, tile, model, Iinp, Iout, & & LBi, UBi, LBj, UBj) CASE DEFAULT IF (Master) WRITE (stdout,10) TLF(ng)%IOtype exit_flag=3 END SELECT IF (FoundError(exit_flag, NoError, 80, MyFile)) RETURN ! 10 FORMAT (' WRT_AUG_IMP - Illegal output tile type, io_type = ',i0, & & /,15x,'Check KeyWord ''OUT_LIB'' in ''roms.in''.' ) ! RETURN END SUBROUTINE wrt_aug_imp ! !*********************************************************************** SUBROUTINE wrt_aug_imp_nf90 (ng, tile, model, Iinp, Iout, & & LBi, UBi, LBj, UBj) !*********************************************************************** ! USE mod_netcdf ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model, Iinp, Iout integer, intent(in) :: LBi, UBi, LBj, UBj ! ! Local variable declarations. ! integer :: itrc, gfactor, gtype, status ! real(dp) :: scale ! character (len=*), parameter :: MyFile = & & "ROMS/Utility/wrt_aug_imp.F"//", wrt_aug_imp_nf90" ! !----------------------------------------------------------------------- ! Set lower and upper tile bounds and staggered variables bounds for ! this horizontal domain partition. Notice that if tile=-1, it will ! set the values for the global grid. !----------------------------------------------------------------------- ! integer :: Istr, IstrB, IstrP, IstrR, IstrT, IstrM, IstrU integer :: Iend, IendB, IendP, IendR, IendT integer :: Jstr, JstrB, JstrP, JstrR, JstrT, JstrM, JstrV integer :: Jend, JendB, JendP, JendR, JendT integer :: Istrm3, Istrm2, Istrm1, IstrUm2, IstrUm1 integer :: Iendp1, Iendp2, Iendp2i, Iendp3 integer :: Jstrm3, Jstrm2, Jstrm1, JstrVm2, JstrVm1 integer :: Jendp1, Jendp2, Jendp2i, Jendp3 ! Istr =BOUNDS(ng) % Istr (tile) IstrB =BOUNDS(ng) % IstrB (tile) IstrM =BOUNDS(ng) % IstrM (tile) IstrP =BOUNDS(ng) % IstrP (tile) IstrR =BOUNDS(ng) % IstrR (tile) IstrT =BOUNDS(ng) % IstrT (tile) IstrU =BOUNDS(ng) % IstrU (tile) Iend =BOUNDS(ng) % Iend (tile) IendB =BOUNDS(ng) % IendB (tile) IendP =BOUNDS(ng) % IendP (tile) IendR =BOUNDS(ng) % IendR (tile) IendT =BOUNDS(ng) % IendT (tile) Jstr =BOUNDS(ng) % Jstr (tile) JstrB =BOUNDS(ng) % JstrB (tile) JstrM =BOUNDS(ng) % JstrM (tile) JstrP =BOUNDS(ng) % JstrP (tile) JstrR =BOUNDS(ng) % JstrR (tile) JstrT =BOUNDS(ng) % JstrT (tile) JstrV =BOUNDS(ng) % JstrV (tile) Jend =BOUNDS(ng) % Jend (tile) JendB =BOUNDS(ng) % JendB (tile) JendP =BOUNDS(ng) % JendP (tile) JendR =BOUNDS(ng) % JendR (tile) JendT =BOUNDS(ng) % JendT (tile) ! Istrm3 =BOUNDS(ng) % Istrm3 (tile) ! Istr-3 Istrm2 =BOUNDS(ng) % Istrm2 (tile) ! Istr-2 Istrm1 =BOUNDS(ng) % Istrm1 (tile) ! Istr-1 IstrUm2=BOUNDS(ng) % IstrUm2(tile) ! IstrU-2 IstrUm1=BOUNDS(ng) % IstrUm1(tile) ! IstrU-1 Iendp1 =BOUNDS(ng) % Iendp1 (tile) ! Iend+1 Iendp2 =BOUNDS(ng) % Iendp2 (tile) ! Iend+2 Iendp2i=BOUNDS(ng) % Iendp2i(tile) ! Iend+2 interior Iendp3 =BOUNDS(ng) % Iendp3 (tile) ! Iend+3 Jstrm3 =BOUNDS(ng) % Jstrm3 (tile) ! Jstr-3 Jstrm2 =BOUNDS(ng) % Jstrm2 (tile) ! Jstr-2 Jstrm1 =BOUNDS(ng) % Jstrm1 (tile) ! Jstr-1 JstrVm2=BOUNDS(ng) % JstrVm2(tile) ! JstrV-2 JstrVm1=BOUNDS(ng) % JstrVm1(tile) ! JstrV-1 Jendp1 =BOUNDS(ng) % Jendp1 (tile) ! Jend+1 Jendp2 =BOUNDS(ng) % Jendp2 (tile) ! Jend+2 Jendp2i=BOUNDS(ng) % Jendp2i(tile) ! Jend+2 interior Jendp3 =BOUNDS(ng) % Jendp3 (tile) ! Jend+3 ! SourceFile=MyFile ! !----------------------------------------------------------------------- ! Determine variables to read and their availability. !----------------------------------------------------------------------- ! IF (FoundError(exit_flag, NoError, 117, MyFile)) RETURN ! ! Report. ! IF (Master) WRITE (stdout,10) Iout, TRIM(TLF(ng)%name) ! ! Set grid type factor to write full (gfactor=1) fields or water ! points (gfactor=-1) fields only. ! gfactor=1 ! ! Process free-surface weak-constraint impulse forcing. ! scale=1.0_dp gtype=gfactor*r2dvar status=nf_fwrite2d(ng, model, TLF(ng)%ncid, idZtlf, & & TLF(ng)%Vid(idZtlf), & & Iout, gtype, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % rmask, & & OCEAN(ng) % tl_zeta(:,:,Iinp)) IF (FoundError(status, nf90_noerr, 144, MyFile)) THEN IF (Master) THEN WRITE (stdout,20) TRIM(Vname(1,idZtlf)), Iout, & & TRIM(TLF(ng)%name) END IF exit_flag=3 ioerror=status RETURN END IF ! ! Process 3D U-momentum weak-constraint impulse forcing. ! scale=1.0_dp gtype=gfactor*u3dvar status=nf_fwrite3d(ng, model, TLF(ng)%ncid, idUtlf, & & TLF(ng)%Vid(idUtlf), & & Iout, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % umask_full, & & OCEAN(ng) % tl_u(:,:,:,Iinp)) IF (FoundError(status, nf90_noerr, 214, MyFile)) THEN IF (Master) THEN WRITE (stdout,20) TRIM(Vname(1,idUtlf)), Iout, & & TRIM(TLF(ng)%name) END IF exit_flag=3 ioerror=status RETURN END IF ! ! Process 3D V-momentum weak-constraint impulse forcing. ! scale=1.0_dp gtype=gfactor*v3dvar status=nf_fwrite3d(ng, model, TLF(ng)%ncid, idVtlf, & & TLF(ng)%Vid(idVtlf), & & Iout, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % vmask_full, & & OCEAN(ng) % tl_v(:,:,:,Iinp)) IF (FoundError(status, nf90_noerr, 236, MyFile)) THEN IF (Master) THEN WRITE (stdout,20) TRIM(Vname(1,idVtlf)), Iout, & & TRIM(TLF(ng)%name) END IF exit_flag=3 ioerror=status RETURN END IF ! ! Process tracer type variables impulses. ! scale=1.0_dp gtype=gfactor*r3dvar DO itrc=1,NT(ng) status=nf_fwrite3d(ng, model, TLF(ng)%ncid, idTtlf(itrc), & & TLF(ng)%Tid(itrc), & & Iout, gtype, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % rmask_full, & & OCEAN(ng) % tl_t(:,:,:,Iinp,itrc)) IF (FoundError(status, nf90_noerr, 259, MyFile)) THEN IF (Master) THEN WRITE (stdout,20) TRIM(Vname(1,idTtlf(itrc))), Iout, & & TRIM(TLF(ng)%name) END IF exit_flag=3 ioerror=status RETURN END IF END DO ! !----------------------------------------------------------------------- ! Synchronize impulse NetCDF file to disk to allow other processes ! to access data immediately after it is written. !----------------------------------------------------------------------- ! CALL netcdf_sync (ng, model, TLF(ng)%name, TLF(ng)%ncid) IF (FoundError(exit_flag, NoError, 277, MyFile)) RETURN ! 10 FORMAT (2x,'WRT_AUG_IMP_NF90 - writing augmented adjoint', & & ' impulses, record: ',i0,/,22x,'file: ',a) 20 FORMAT (/,' WRT_AUG_IMP_NF90 - error while writing variable: ',a, & & 2x,'at time record = ',i0, & & /,20x,'into NetCDF file: ',a) ! RETURN END SUBROUTINE wrt_aug_imp_nf90 END MODULE wrt_aug_imp_mod