MODULE wrt_impulse_mod ! !git $Id$ !svn $Id: wrt_impulse.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 ! !======================================================================= ! ! ! These subroutine read in requested adjoint solution from input ! ! NetCDF (saved at nADJ time-step intervals) and then writes to ! ! output impulse forcing NetCDF file in ascending time order ! ! since it is processed by the TL and RP models. ! ! ! ! On Input: ! ! ! ! ng Nested grid number. ! ! tile Domain partition. ! ! model Calling model identifier. ! ! INPncname Input adjoint solution NetCDF file name. ! ! ! ! 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_grid USE mod_iounits USE mod_ncparam USE mod_ocean USE mod_scalars ! USE nf_fread2d_mod, ONLY : nf_fread2d USE nf_fread3d_mod, ONLY : nf_fread3d 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_impulse PRIVATE :: wrt_impulse_nf90 ! CONTAINS ! !*********************************************************************** SUBROUTINE wrt_impulse (ng, tile, model, INPncname) !*********************************************************************** ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model ! character (len=*), intent(in) :: INPncname ! ! Local variable declarations. ! integer :: LBi, UBi, LBj, UBj ! character (len=*), parameter :: MyFile = & & "ROMS/Utility/wrt_impulse.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 (LZE(ng)%IOtype) CASE (io_nf90) CALL wrt_impulse_nf90 (ng, tile, model, INPncname, & & LBi, UBi, LBj, UBj) CASE DEFAULT IF (Master) WRITE (stdout,10) LZE(ng)%IOtype exit_flag=3 END SELECT IF (FoundError(exit_flag, NoError, 98, MyFile)) RETURN ! 10 FORMAT (' WRT_IMPULSE - Illegal output file type, io_type = ',i0, & & /,15x,'Check KeyWord ''OUT_LIB'' in ''roms.in''.') ! RETURN END SUBROUTINE wrt_impulse ! !*********************************************************************** SUBROUTINE wrt_impulse_nf90 (ng, tile, model, INPncname, & & LBi, UBi, LBj, UBj) !*********************************************************************** ! USE mod_netcdf ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model ! character (len=*), intent(in) :: INPncname ! ! Local variable declarations. ! integer :: LBi, UBi, LBj, UBj integer :: Iinp, Iout, Irec, MyType, Nrec integer :: INPncid, INPvid integer :: i, gtype, status, varid integer :: ibuffer(2), Vsize(4) ! real(r8) :: Fmin, Fmax real(dp) :: scale real(dp) :: inp_time(1) ! character (len=*), parameter :: MyFile = & & "ROMS/Utility/wrt_impulse.F"//", wrt_impulse_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. !----------------------------------------------------------------------- ! ! Inquire about the dimensions and check for consistency. ! CALL netcdf_check_dim (ng, model, INPncname) IF (FoundError(exit_flag, NoError, 148, MyFile)) RETURN Nrec=rec_size ! ! Inquire about the variables. ! CALL netcdf_inq_var (ng, model, INPncname) IF (FoundError(exit_flag, NoError, 154, MyFile)) RETURN ! ! Set Vsize to zero to deactivate interpolation of input data to model ! grid in "nf_fread2d". ! DO i=1,4 Vsize(i)=0 END DO ! IF (Master) WRITE (stdout,10) Nrec-1, TRIM(TLF(ng)%name) ! !----------------------------------------------------------------------- ! Read adjoint solution and convert to impulse forcing. Then, write ! impulse forcing into output NetCDF file. !----------------------------------------------------------------------- ! ! Open input NetCDF file. ! CALL netcdf_open (ng, model, INPncname, 0, INPncid) IF (FoundError(exit_flag, NoError, 177, MyFile)) THEN WRITE (stdout,20) TRIM(INPncname) RETURN END IF ! ! Process each record in input adjoint NetCDF except last. Note that ! the adjoint records are processed backwards (Nrec-1:1) and written ! in ascending time order (Iout initialized to 0) since the weak ! constraint forcing will be read by the TL and RP models. Record ! Nrec is not processed since it is not needed. ! Iinp=1 Iout=0 scale=1.0_dp ! DO Irec=Nrec-1,1,-1 Iout=Iout+1 ! ! Process time. ! IF (find_string(var_name, n_var, Vname(1,idtime), INPvid)) THEN CALL netcdf_get_time (ng, model, INPncname, Vname(1,idtime), & & Rclock%DateNumber, inp_time, & & ncid = INPncid, & & start = (/Irec/), total = (/1/)) IF (FoundError(exit_flag, NoError, 211, MyFile)) RETURN ! CALL netcdf_put_fvar (ng, model, TLF(ng)%name, & & Vname(1,idtime), inp_time, & & (/Iout/), (/1/), & & ncid = TLF(ng)%ncid) IF (FoundError(exit_flag, NoError, 217, MyFile)) RETURN ELSE IF (Master) WRITE (stdout,30) TRIM(Vname(1,idtime)), & & TRIM(INPncname) exit_flag=2 END IF ! ! Process free-surface weak-constraint impulse forcing. ! IF (find_string(var_name, n_var, Vname(1,idZtlf), INPvid)) THEN gtype=var_flag(INPvid)*r2dvar MyType=gtype status=nf_fread2d(ng, model, INPncname, INPncid, & & Vname(1,idZtlf), INPvid, & & Irec, MyType, Vsize, & & LBi, UBi, LBj, UBj, & & scale, Fmin, Fmax, & & GRID(ng) % rmask, & & OCEAN(ng) % ad_zeta(:,:,Iinp)) IF (FoundError(status, nf90_noerr, 243, MyFile)) THEN IF (Master) THEN WRITE (stdout,40) TRIM(Vname(1,idZtlf)), Irec, & & TRIM(INPncname) END IF exit_flag=2 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,50) TRIM(Vname(1,idZtlf)), Irec, & & Fmin, Fmax END IF END IF ! MyType=gtype status=nf_fwrite2d(ng, model, TLF(ng)%ncid, idZtlf, & & TLF(ng)%Vid(idZtlf), & & Iout, MyType, & & LBi, UBi, LBj, UBj, scale, & & GRID(ng) % rmask, & & OCEAN(ng) % ad_zeta(:,:,Iinp)) IF (FoundError(status, nf90_noerr, 271, MyFile)) THEN IF (Master) THEN WRITE (stdout,60) TRIM(Vname(1,idZtlf)), Irec, & & TRIM(TLF(ng)%name) END IF exit_flag=3 ioerror=status RETURN END IF ELSE IF (Master) WRITE (stdout,30) TRIM(Vname(1,idZtlf)), & & TRIM(INPncname) exit_flag=2 RETURN END IF ! ! Process 3D U-momentum weak-constraint impulse forcing. ! IF (find_string(var_name, n_var, Vname(1,idUtlf), INPvid)) THEN gtype=var_flag(INPvid)*u3dvar MyType=gtype status=nf_fread3d(ng, model, INPncname, INPncid, & & Vname(1,idUtlf), INPvid, & & Irec, MyType, Vsize, & & LBi, UBi, LBj, UBj, 1, N(ng), & & scale, Fmin, Fmax, & & GRID(ng) % umask_full, & & OCEAN(ng) % ad_u(:,:,:,Iinp)) IF (FoundError(status, nf90_noerr, 436, MyFile)) THEN IF (Master) THEN WRITE (stdout,40) TRIM(Vname(1,idUtlf)), Irec, & & TRIM(INPncname) END IF exit_flag=2 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,50) TRIM(Vname(1,idUtlf)), Irec, & & Fmin, Fmax END IF END IF ! status=nf_fwrite3d(ng, model, TLF(ng)%ncid, idUtlf, & & TLF(ng)%Vid(idUtlf), & & Iout, MyType, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % umask_full, & & OCEAN(ng) % ad_u(:,:,:,Iinp)) IF (FoundError(status, nf90_noerr, 463, MyFile)) THEN IF (Master) THEN WRITE (stdout,60) TRIM(Vname(1,idUtlf)), Irec, & & TRIM(TLF(ng)%name) END IF exit_flag=3 ioerror=status RETURN END IF ELSE IF (Master) WRITE (stdout,30) TRIM(Vname(1,idUtlf)), & & TRIM(INPncname) exit_flag=2 RETURN END IF ! ! Process 3D V-momentum weak-constraint impulse forcing. ! IF (find_string(var_name, n_var, Vname(1,idVtlf), INPvid)) THEN gtype=var_flag(INPvid)*v3dvar MyType=gtype status=nf_fread3d(ng, model, INPncname, INPncid, & & Vname(1,idVtlf), INPvid, & & Irec, MyType, Vsize, & & LBi, UBi, LBj, UBj, 1, N(ng), & & scale, Fmin, Fmax, & & GRID(ng) % vmask_full, & & OCEAN(ng) % ad_v(:,:,:,Iinp)) IF (FoundError(status, nf90_noerr, 498, MyFile)) THEN IF (Master) THEN WRITE (stdout,40) TRIM(Vname(1,idVtlf)), Irec, & & TRIM(INPncname) END IF exit_flag=2 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,50) TRIM(Vname(1,idVtlf)), Irec, & & Fmin, Fmax END IF END IF ! status=nf_fwrite3d(ng, model, TLF(ng)%ncid, idVtlf, & & TLF(ng)%Vid(idVtlf), & & Iout, MyType, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % vmask_full, & & OCEAN(ng) % ad_v(:,:,:,Iinp)) IF (FoundError(status, nf90_noerr, 525, MyFile)) THEN IF (Master) THEN WRITE (stdout,60) TRIM(Vname(1,idVtlf)), Irec, & & TRIM(TLF(ng)%name) END IF exit_flag=3 ioerror=status RETURN END IF ELSE IF (Master) WRITE (stdout,30) TRIM(Vname(1,idVtlf)), & & TRIM(INPncname) exit_flag=2 RETURN END IF ! ! Process tracer type variables impulses. ! DO i=1,NT(ng) IF (find_string(var_name, n_var, Vname(1,idTtlf(i)), & & INPvid)) THEN gtype=var_flag(INPvid)*r3dvar MyType=gtype status=nf_fread3d(ng, model, INPncname, INPncid, & & Vname(1,idTtlf(i)), INPvid, & & Irec, MyType, Vsize, & & LBi, UBi, LBj, UBj, 1, N(ng), & & scale, Fmin, Fmax, & & GRID(ng) % rmask, & & OCEAN(ng) % ad_t(:,:,:,Iinp,i)) IF (FoundError(status, nf90_noerr, 562, MyFile)) THEN IF (Master) THEN WRITE (stdout,40) TRIM(Vname(1,idTtlf(i))), Irec, & & TRIM(INPncname) END IF exit_flag=2 ioerror=status RETURN ELSE IF (Master) THEN WRITE (stdout,50) TRIM(Vname(1,idTtlf(i))), Irec, & & Fmin, Fmax END IF END IF ! status=nf_fwrite3d(ng, model, TLF(ng)%ncid, idTtlf(i), & & TLF(ng)%Tid(i), & & Iout, MyType, & & LBi, UBi, LBj, UBj, 1, N(ng), scale, & & GRID(ng) % rmask, & & OCEAN(ng) % ad_t(:,:,:,Iinp,i)) IF (FoundError(status, nf90_noerr, 589, MyFile)) THEN IF (Master) THEN WRITE (stdout,60) TRIM(Vname(1,idTtlf(i))), Irec, & & TRIM(TLF(ng)%name) END IF exit_flag=3 ioerror=status RETURN END IF ELSE IF (Master) WRITE (stdout,30) TRIM(Vname(1,idTtlf(i))), & & TRIM(INPncname) exit_flag=2 RETURN END IF END DO END DO ! ! Close input NetCDF file. ! CALL netcdf_close (ng, model, INPncid, INPncname, .FALSE.) IF (FoundError(exit_flag, NoError, 611, MyFile)) THEN WRITE (stdout,70) TRIM(INPncname) RETURN END IF ! !----------------------------------------------------------------------- ! Synchronize impulse NetCDF file to disk to allow other processes ! to access data immediately after it is written. !----------------------------------------------------------------------- ! CALL netcdf_sync (ng, model, INPncname, TLF(ng)%ncid) IF (FoundError(exit_flag, NoError, 622, MyFile)) RETURN ! 10 FORMAT (2x,'WRT_IMPULSE_NF90 - processing convolved adjoint', & & ' impulses, records: 1 to ',i0,/,21x,'file: ',a) 20 FORMAT (/,' WRT_IMPULSE_NF90 - unable to open input NetCDF', & & ' file: ',a) 30 FORMAT (/,' WRT_IMPULSE_NF90 - cannot find state variable: ',a, & & /,20x,'in input NetCDF file: ',a) 40 FORMAT (/,' WRT_IMPULSE_NF90 - error while reading variable: ',a, & & 2x,'at time record = ',i0, & & /,20x,'in input NetCDF file: ',a) 50 FORMAT (19x,'- ',a,/,22x,'(Rec = ',i0,' Min = ',1p,e15.8, & & ' Max = ',1p,e15.8,')') 60 FORMAT (/,' WRT_IMPULSE_NF90 - error while writing variable: ',a, & & 2x,'at time record = ',i0,/,20x,'into NetCDF file: ',a) 70 FORMAT (/,' WRT_IMPULSE_NF90 - unable to close input NetCDF', & & ' file: ',a) ! RETURN END SUBROUTINE wrt_impulse_nf90 END MODULE wrt_impulse_mod