#include "cppdefs.h" MODULE tl_def_his_mod #if defined TANGENT || defined TL_IOMS ! !git $Id$ !svn $Id: tl_def_his.F 1151 2023-02-09 03:08:53Z 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 creates tangent linear history file using either ! ! the standard NetCDF library or the Parallel-IO (PIO) library. ! ! It defines its dimensions, attributes, and variables. ! ! ! !======================================================================= ! USE mod_param USE mod_parallel # ifdef BIOLOGY USE mod_biology # endif # ifdef FOUR_DVAR USE mod_fourdvar # endif USE mod_iounits USE mod_ncparam USE mod_scalars # if defined SEDIMENT_NOT_YET || defined BBL_MODEL_NOT_YET USE mod_sediment # endif ! USE def_dim_mod, ONLY : def_dim USE def_info_mod, ONLY : def_info USE def_var_mod, ONLY : def_var USE strings_mod, ONLY : FoundError USE wrt_info_mod, ONLY : wrt_info ! implicit none ! PUBLIC :: tl_def_his PRIVATE :: tl_def_his_nf90 #if defined PIO_LIB && defined DISTRIBUTE PRIVATE :: tl_def_his_pio #endif ! CONTAINS ! !*********************************************************************** SUBROUTINE tl_def_his (ng, ldef) !*********************************************************************** ! ! Imported variable declarations. ! logical, intent(in) :: ldef ! integer, intent(in) :: ng ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & __FILE__ ! !----------------------------------------------------------------------- ! Create a new history file according to IO type. !----------------------------------------------------------------------- ! SELECT CASE (TLM(ng)%IOtype) CASE (io_nf90) CALL tl_def_his_nf90 (ng, ldef) # if defined PIO_LIB && defined DISTRIBUTE CASE (io_pio) CALL tl_def_his_pio (ng, ldef) # endif CASE DEFAULT IF (Master) WRITE (stdout,10) TLM(ng)%IOtype exit_flag=3 END SELECT IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! 10 FORMAT (' TL_DEF_HIS - Illegal output file type, io_type = ',i0, & & /,14x,'Check KeyWord ''OUT_LIB'' in ''roms.in''.') ! RETURN END SUBROUTINE tl_def_his ! !*********************************************************************** SUBROUTINE tl_def_his_nf90 (ng, ldef) !*********************************************************************** ! USE mod_netcdf ! ! Imported variable declarations. ! integer, intent(in) :: ng ! logical, intent(in) :: ldef ! ! Local variable declarations. ! logical :: got_var(NV) ! integer, parameter :: Natt = 25 integer :: i, j, ifield, itrc, nvd3, nvd4 integer :: recdim, status, varid # ifdef ADJUST_BOUNDARY integer :: IorJdim, brecdim # endif # if defined ADJUST_STFLUX || defined ADJUST_WSTRESS integer :: frecdim # endif integer :: DimIDs(nDimID) integer :: t2dgrd(3), u2dgrd(3), v2dgrd(3) # ifdef ADJUST_BOUNDARY integer :: t2dobc(4) # endif # ifdef SOLVE3D integer :: t3dgrd(4), u3dgrd(4), v3dgrd(4), w3dgrd(4) # ifdef ADJUST_BOUNDARY integer :: t3dobc(5) # endif # ifdef ADJUST_STFLUX integer :: t3dfrc(4) # endif # endif # ifdef ADJUST_WSTRESS integer :: u3dfrc(4), v3dfrc(4) # endif ! real(r8) :: Aval(6) ! character (len=256) :: ncname character (len=MaxLen) :: Vinfo(Natt) ! character (len=*), parameter :: MyFile = & & __FILE__//", tl_def_his_nf90" ! SourceFile=MyFile ! !----------------------------------------------------------------------- ! Set and report file name. !----------------------------------------------------------------------- ! IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ncname=TLM(ng)%name ! IF (Master) THEN IF (ldef) THEN WRITE (stdout,10) ng, TRIM(ncname) ELSE WRITE (stdout,20) ng, TRIM(ncname) END IF END IF ! !======================================================================= ! Create a new tangent linear history file. !======================================================================= ! DEFINE : IF (ldef) THEN CALL netcdf_create (ng, iTLM, TRIM(ncname), TLM(ng)%ncid) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) THEN IF (Master) WRITE (stdout,30) TRIM(ncname) RETURN END IF ! !----------------------------------------------------------------------- ! Define file dimensions. !----------------------------------------------------------------------- ! DimIDs=0 ! status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xi_rho', & & IOBOUNDS(ng)%xi_rho, DimIDs( 1)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xi_u', & & IOBOUNDS(ng)%xi_u, DimIDs( 2)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xi_v', & & IOBOUNDS(ng)%xi_v, DimIDs( 3)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xi_psi', & & IOBOUNDS(ng)%xi_psi, DimIDs( 4)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'eta_rho', & & IOBOUNDS(ng)%eta_rho, DimIDs( 5)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'eta_u', & & IOBOUNDS(ng)%eta_u, DimIDs( 6)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'eta_v', & & IOBOUNDS(ng)%eta_v, DimIDs( 7)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'eta_psi', & & IOBOUNDS(ng)%eta_psi, DimIDs( 8)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef ADJUST_BOUNDARY status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'IorJ', & & IOBOUNDS(ng)%IorJ, IorJdim) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # if defined WRITE_WATER && defined MASKING status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xy_rho', & & IOBOUNDS(ng)%xy_rho, DimIDs(17)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xy_u', & & IOBOUNDS(ng)%xy_u, DimIDs(18)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xy_v', & & IOBOUNDS(ng)%xy_v, DimIDs(19)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # ifdef SOLVE3D # if defined WRITE_WATER && defined MASKING status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xyz_rho', & & IOBOUNDS(ng)%xy_rho*N(ng), DimIDs(20)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xyz_u', & & IOBOUNDS(ng)%xy_u*N(ng), DimIDs(21)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xyz_v', & & IOBOUNDS(ng)%xy_v*N(ng), DimIDs(22)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xyz_w', & & IOBOUNDS(ng)%xy_rho*(N(ng)+1), DimIDs(23)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'N', & & N(ng), DimIDs( 9)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 's_rho', & & N(ng), DimIDs( 9)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 's_w', & & N(ng)+1, DimIDs(10)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'tracer', & & NT(ng), DimIDs(11)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef SEDIMENT status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'NST', & & NST, DimIDs(32)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'Nbed', & & Nbed, DimIDs(16)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # if defined WRITE_WATER && defined MASKING status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'xybed', & & IOBOUNDS(ng)%xy_rho*Nbed, DimIDs(24)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # endif # ifdef ECOSIM status=def_dim(ng, iNLM, TLM(ng)%ncid, ncname, 'Nbands', & & NBands, DimIDs(33)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'Nphy', & & Nphy, DimIDs(25)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'Nbac', & & Nbac, DimIDs(26)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'Ndom', & & Ndom, DimIDs(27)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'Nfec', & & Nfec, DimIDs(28)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # endif status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'boundary', & & 4, DimIDs(14)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef FOUR_DVAR status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'Nstate', & & NstateVar(ng), DimIDs(29)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # if defined ADJUST_STFLUX || defined ADJUST_WSTRESS status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'frc_adjust', & & Nfrec(ng), DimIDs(30)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # ifdef ADJUST_BOUNDARY status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, 'obc_adjust', & & Nbrec(ng), DimIDs(31)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif status=def_dim(ng, iTLM, TLM(ng)%ncid, ncname, & & TRIM(ADJUSTL(Vname(5,idtime))), & & nf90_unlimited, DimIDs(12)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN recdim=DimIDs(12) # if defined ADJUST_STFLUX || defined ADJUST_WSTRESS frecdim=DimIDs(30) # endif # ifdef ADJUST_BOUNDARY brecdim=DimIDs(31) # endif ! ! Set number of dimensions for output variables. ! # if defined WRITE_WATER && defined MASKING nvd3=2 nvd4=2 # else nvd3=3 nvd4=4 # endif ! ! Define dimension vectors for staggered tracer type variables. ! # if defined WRITE_WATER && defined MASKING t2dgrd(1)=DimIDs(17) t2dgrd(2)=DimIDs(12) # ifdef SOLVE3D t3dgrd(1)=DimIDs(20) t3dgrd(2)=DimIDs(12) # endif # else t2dgrd(1)=DimIDs( 1) t2dgrd(2)=DimIDs( 5) t2dgrd(3)=DimIDs(12) # ifdef SOLVE3D t3dgrd(1)=DimIDs( 1) t3dgrd(2)=DimIDs( 5) t3dgrd(3)=DimIDs( 9) t3dgrd(4)=DimIDs(12) # endif # ifdef ADJUST_STFLUX t3dfrc(1)=DimIDs( 1) t3dfrc(2)=DimIDs( 5) t3dfrc(3)=frecdim t3dfrc(4)=DimIDs(12) # endif # endif # ifdef ADJUST_BOUNDARY t2dobc(1)=IorJdim t2dobc(2)=DimIDs(14) t2dobc(3)=brecdim t2dobc(4)=DimIDs(12) # ifdef SOLVE3D t3dobc(1)=IorJdim t3dobc(2)=DimIDs( 9) t3dobc(3)=DimIDs(14) t3dobc(4)=brecdim t3dobc(5)=DimIDs(12) # endif # endif ! ! Define dimension vectors for staggered u-momentum type variables. ! # if defined WRITE_WATER && defined MASKING u2dgrd(1)=DimIDs(18) u2dgrd(2)=DimIDs(12) # ifdef SOLVE3D u3dgrd(1)=DimIDs(21) u3dgrd(2)=DimIDs(12) # endif # else u2dgrd(1)=DimIDs( 2) u2dgrd(2)=DimIDs( 6) u2dgrd(3)=DimIDs(12) # ifdef SOLVE3D u3dgrd(1)=DimIDs( 2) u3dgrd(2)=DimIDs( 6) u3dgrd(3)=DimIDs( 9) u3dgrd(4)=DimIDs(12) # endif # ifdef ADJUST_WSTRESS u3dfrc(1)=DimIDs( 2) u3dfrc(2)=DimIDs( 6) u3dfrc(3)=frecdim u3dfrc(4)=DimIDs(12) # endif # endif ! ! Define dimension vectors for staggered v-momentum type variables. ! # if defined WRITE_WATER && defined MASKING v2dgrd(1)=DimIDs(19) v2dgrd(2)=DimIDs(12) # ifdef SOLVE3D v3dgrd(1)=DimIDs(22) v3dgrd(2)=DimIDs(12) # endif # else v2dgrd(1)=DimIDs( 3) v2dgrd(2)=DimIDs( 7) v2dgrd(3)=DimIDs(12) # ifdef SOLVE3D v3dgrd(1)=DimIDs( 3) v3dgrd(2)=DimIDs( 7) v3dgrd(3)=DimIDs( 9) v3dgrd(4)=DimIDs(12) # endif # ifdef ADJUST_WSTRESS v3dfrc(1)=DimIDs( 3) v3dfrc(2)=DimIDs( 7) v3dfrc(3)=frecdim v3dfrc(4)=DimIDs(12) # endif # endif # ifdef SOLVE3D ! ! Define dimension vector for staggered w-momentum type variables. ! # if defined WRITE_WATER && defined MASKING w3dgrd(1)=DimIDs(23) w3dgrd(2)=DimIDs(12) # else w3dgrd(1)=DimIDs( 1) w3dgrd(2)=DimIDs( 5) w3dgrd(3)=DimIDs(10) w3dgrd(4)=DimIDs(12) # endif # endif ! ! Initialize unlimited time record dimension. ! TLM(ng)%Rindex=0 ! ! Initialize local information variable arrays. ! DO i=1,Natt DO j=1,LEN(Vinfo(1)) Vinfo(i)(j:j)=' ' END DO END DO DO i=1,6 Aval(i)=0.0_r8 END DO ! !----------------------------------------------------------------------- ! Define time-recordless information variables. !----------------------------------------------------------------------- ! CALL def_info (ng, iTLM, TLM(ng)%ncid, ncname, DimIDs) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! !----------------------------------------------------------------------- ! Define time-varying variables. !----------------------------------------------------------------------- ! ! Define model time. ! Vinfo( 1)=Vname(1,idtime) Vinfo( 2)=Vname(2,idtime) WRITE (Vinfo( 3),'(a,a)') 'seconds since ', TRIM(Rclock%string) Vinfo( 4)=TRIM(Rclock%calendar) Vinfo(14)=Vname(4,idtime) Vinfo(21)=Vname(6,idtime) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idtime), & & NF_TYPE, 1, (/recdim/), Aval, Vinfo, ncname, & & SetParAccess = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef PROPAGATOR ! ! Define Ritz eigenvalues and Ritz eigenvectors Euclidean norm. ! Vinfo( 1)='Ritz_rvalue' Vinfo( 2)='real Ritz eigenvalues' status=def_var(ng, iTLM, TLM(ng)%ncid, varid, NF_TYPE, & & 1, (/recdim/), Aval, Vinfo, ncname, & & SetParAccess = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # if defined FT_EIGENMODES Vinfo( 1)='Ritz_ivalue' Vinfo( 2)='imaginary Ritz eigenvalues' status=def_var(ng, iTLM, TLM(ng)%ncid, varid, NF_TYPE, & & 1, (/recdim/), Aval, Vinfo, ncname, & & SetParAccess = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif Vinfo( 1)='Ritz_norm' Vinfo( 2)='Ritz eigenvectors Euclidean norm' status=def_var(ng, iTLM, TLM(ng)%ncid, varid, NF_TYPE, & & 1, (/recdim/), Aval, Vinfo, ncname, & & SetParAccess = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # ifdef ADJUST_WSTRESS ! ! Define surface U-momentum stress. Notice that the stress has its ! own fixed time-dimension (of size Nfrec) to allow 4DVAR adjustments ! at other times in addition to initialization time. ! Vinfo( 1)=Vname(1,idUsms) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUsms)) Vinfo( 3)='meter2 second-2' Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idUsms,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idUsms), & & NF_FOUT, nvd4, u3dfrc, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! ! Define surface V-momentum stress. ! Vinfo( 1)=Vname(1,idVsms) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVsms)) Vinfo( 3)='meter2 second-2' Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVsms,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idVsms), & & NF_FOUT, nvd4, v3dfrc, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # if defined FORCING_SV || defined STOCHASTIC_OPT || \ defined HESSIAN_SO || defined HESSIAN_FSV ! ! Define surface U-momentum stress. ! IF (Hout(idUsms,ng)) THEN Vinfo( 1)=Vname(1,idUsms) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUsms)) Vinfo( 3)=Vname(3,idUsms) Vinfo(14)=Vname(4,idUsms) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idUsms) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idUsms,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idUsms), & & NF_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF ! ! Define surface V-momentum stress. ! IF (Hout(idVsms,ng)) THEN Vinfo( 1)=Vname(1,idVsms) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVsms)) Vinfo( 3)=Vname(3,idVsms) Vinfo(14)=Vname(4,idVsms) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idVsms) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVsms,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idVsms), & & NF_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF ! ! Define surface tracer fluxes. ! DO itrc=1,NT(ng) IF (Hout(idTsur(itrc),ng)) THEN Vinfo( 1)=Vname(1,idTsur(itrc)) WRITE (Vinfo( 2),40) TRIM(Vname(2,idTsur(itrc))) Vinfo( 3)=Vname(3,idTsur(itrc)) IF (itrc.eq.itemp) THEN Vinfo(11)='upward flux, cooling' Vinfo(12)='downward flux, heating' ELSE IF (itrc.eq.isalt) THEN Vinfo(11)='upward flux, freshening (net precipitation)' Vinfo(12)='downward flux, salting (net evaporation)' END IF Vinfo(14)=Vname(4,idTsur(itrc)) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idTsur(itrc)) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idTsur(itrc),ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, & & TLM(ng)%Vid(idTsur(itrc)), NF_FOUT, & & nvd3, t2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF END DO # endif # if defined ADJUST_STFLUX && defined SOLVE3D ! ! Define surface net heat flux. Notice that different tracer fluxes ! are written at their own fixed time-dimension (of size Nfrec) to ! allow 4DVAR adjustments at other times in addition to initial time. ! DO itrc=1,NT(ng) IF (Lstflux(itrc,ng)) THEN Vinfo( 1)=Vname(1,idTsur(itrc)) WRITE (Vinfo( 2),40) TRIM(Vname(2,idTsur(itrc))) IF (itrc.eq.itemp) THEN Vinfo( 3)='Celsius meter second-1' Vinfo(11)='upward flux, cooling' Vinfo(12)='downward flux, heating' ELSE IF (itrc.eq.isalt) THEN Vinfo( 3)='meter second-1' Vinfo(11)='upward flux, freshening (net precipitation)' Vinfo(12)='downward flux, salting (net evaporation)' END IF Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idTsur(itrc),ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, & & TLM(ng)%Vid(idTsur(itrc)), NF_FOUT, & & nvd4, t3dfrc, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF END DO # endif ! ! Define free-surface. ! IF (Hout(idFsur,ng)) THEN Vinfo( 1)=Vname(1,idFsur) WRITE (Vinfo( 2),40) TRIM(Vname(2,idFsur)) Vinfo( 3)=Vname(3,idFsur) Vinfo(14)=Vname(4,idFsur) Vinfo(16)=Vname(1,idtime) # if !defined WET_DRY && (defined WRITE_WATER && defined MASKING) Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idFsur) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idFsur,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idFsur), & # ifdef WET_DRY & NF_FOUT, nvd3, t2dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) # else & NF_FOUT, nvd3, t2dgrd, Aval, Vinfo, ncname) # endif IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # if defined FORWARD_WRITE && defined FORWARD_RHS Vinfo( 1)=Vname(1,idRzet) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRzet)) Vinfo( 3)=Vname(3,idRzet) Vinfo(14)=Vname(4,idRzet) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idRzet) Vinfo(22)='coordinates' Aval(5)=REAL(r2dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idRzet), & & NF_FOUT, nvd3, t2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif END IF # ifdef ADJUST_BOUNDARY ! ! Define free-surface open boundaries. ! IF (ANY(Lobc(:,isFsur,ng))) THEN ifield=idSbry(isFsur) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(ifield), & & NF_FOUT, 4, t2dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif ! ! Define 2D U-momentum component. ! IF (Hout(idUbar,ng)) THEN Vinfo( 1)=Vname(1,idUbar) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUbar)) Vinfo( 3)=Vname(3,idUbar) Vinfo(14)=Vname(4,idUbar) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idUbar) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idUbar,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idUbar), & & NF_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef FORWARD_WRITE # ifdef FORWARD_RHS Vinfo( 1)=Vname(1,idRu2d) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRu2d)) Vinfo( 3)=Vname(3,idRu2d) Vinfo(14)=Vname(4,idRu2d) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idRu2d) Vinfo(22)='coordinates' Aval(5)=REAL(u2dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idRu2d), & & NF_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # ifdef SOLVE3D # ifdef FORWARD_RHS Vinfo( 1)=Vname(1,idRuct) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRuct)) Vinfo( 3)=Vname(3,idRuct) Vinfo(14)=Vname(4,idRuct) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idRuct) Vinfo(22)='coordinates' Aval(5)=REAL(u2dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idRuct), & & NF_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif Vinfo( 1)=Vname(1,idUfx1) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUfx1)) Vinfo( 3)=Vname(3,idUfx1) Vinfo(14)=Vname(4,idUfx1) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idUfx1) Vinfo(22)='coordinates' Aval(5)=REAL(u2dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idUfx1), & & NF_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN Vinfo( 1)=Vname(1,idUfx2) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUfx2)) Vinfo( 3)=Vname(3,idUfx2) Vinfo(14)=Vname(4,idUfx2) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idUfx2) Vinfo(22)='coordinates' Aval(5)=REAL(u2dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idUfx2), & & NF_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # endif END IF # ifdef ADJUST_BOUNDARY ! ! Define 2D U-momentum component open boundaries. ! IF (ANY(Lobc(:,isUbar,ng))) THEN ifield=idSbry(isUbar) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(ifield), & & NF_FOUT, 4, t2dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif ! ! Define 2D V-momentum component. ! IF (Hout(idVbar,ng)) THEN Vinfo( 1)=Vname(1,idVbar) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVbar)) Vinfo( 3)=Vname(3,idVbar) Vinfo(14)=Vname(4,idVbar) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idVbar) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVbar,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idVbar), & & NF_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef FORWARD_WRITE # ifdef FORWARD_RHS Vinfo( 1)=Vname(1,idRv2d) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRv2d)) Vinfo( 3)=Vname(3,idRv2d) Vinfo(14)=Vname(4,idRv2d) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idRv2d) Vinfo(22)='coordinates' Aval(5)=REAL(v2dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idRv2d), & & NF_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # ifdef SOLVE3D # ifdef FORWARD_RHS Vinfo( 1)=Vname(1,idRvct) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRvct)) Vinfo( 3)=Vname(3,idRvct) Vinfo(14)=Vname(4,idRvct) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idRvct) Vinfo(22)='coordinates' Aval(5)=REAL(v2dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idRvct), & & NF_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif Vinfo( 1)=Vname(1,idVfx1) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVfx1)) Vinfo( 3)=Vname(3,idVfx1) Vinfo(14)=Vname(4,idVfx1) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idVfx1) Vinfo(22)='coordinates' Aval(5)=REAL(v2dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idVfx1), & & NF_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN Vinfo( 1)=Vname(1,idVfx2) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVfx2)) Vinfo( 3)=Vname(3,idVfx2) Vinfo(14)=Vname(4,idVfx2) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idVfx2) Vinfo(22)='coordinates' Aval(5)=REAL(v2dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idVfx2), & & NF_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # endif END IF # ifdef ADJUST_BOUNDARY ! ! Define 2D V-momentum component open boundaries. ! IF (ANY(Lobc(:,isVbar,ng))) THEN ifield=idSbry(isVbar) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(ifield), & & NF_FOUT, 4, t2dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif # ifdef SOLVE3D ! ! Define 3D U-momentum component. ! IF (Hout(idUvel,ng)) THEN Vinfo( 1)=Vname(1,idUvel) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUvel)) Vinfo( 3)=Vname(3,idUvel) Vinfo(14)=Vname(4,idUvel) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idUvel) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idUvel,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idUvel), & & NF_FOUT, nvd4, u3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # if defined FORWARD_WRITE && defined FORWARD_RHS Vinfo( 1)=Vname(1,idRu3d) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRu3d)) Vinfo( 3)=Vname(3,idRu3d) Vinfo(14)=Vname(4,idRu3d) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idRu3d) Vinfo(22)='coordinates' Aval(5)=REAL(u3dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idRu3d), & & NF_FOUT, nvd4, u3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif END IF # ifdef ADJUST_BOUNDARY ! ! Define 3D U-momentum component open boundaries. ! IF (ANY(Lobc(:,isUvel,ng))) THEN ifield=idSbry(isUvel) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(ifield), & & NF_FOUT, 5, t3dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif ! ! Define 3D V-momentum component. ! IF (Hout(idVvel,ng)) THEN Vinfo( 1)=Vname(1,idVvel) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVvel)) Vinfo( 3)=Vname(3,idVvel) Vinfo(14)=Vname(4,idVvel) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idVvel) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVvel,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idVvel), & & NF_FOUT, nvd4, v3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # if defined FORWARD_WRITE && defined FORWARD_RHS Vinfo( 1)=Vname(1,idRv3d) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRv3d)) Vinfo( 3)=Vname(3,idRv3d) Vinfo(14)=Vname(4,idRv3d) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idRv3d) Vinfo(22)='coordinates' Aval(5)=REAL(v3dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idRv3d), & & NF_FOUT, nvd4, v3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif END IF # ifdef ADJUST_BOUNDARY ! ! Define 3D V-momentum component open boundaries. ! IF (ANY(Lobc(:,isVvel,ng))) THEN ifield=idSbry(isVvel) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(ifield), & & NF_FOUT, 5, t3dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif ! ! Define tracer type variables. ! DO itrc=1,NT(ng) IF (Hout(idTvar(itrc),ng)) THEN Vinfo( 1)=Vname(1,idTvar(itrc)) WRITE (Vinfo( 2),40) TRIM(Vname(2,idTvar(itrc))) Vinfo( 3)=Vname(3,idTvar(itrc)) Vinfo(14)=Vname(4,idTvar(itrc)) Vinfo(16)=Vname(1,idtime) # ifdef SEDIMENT_NOT_YET DO i=1,NST IF (itrc.eq.idsed(i)) THEN WRITE (Vinfo(19),50) 1000.0_r8*Sd50(i,ng) END IF END DO # endif # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idTvar(itrc)) Vinfo(22)='coordinates' Aval(5)=REAL(r3dvar,r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Tid(itrc), & & NF_FOUT, nvd4, t3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF END DO # ifdef ADJUST_BOUNDARY ! ! Define tracer type variables open boundaries. ! DO itrc=1,NT(ng) IF (ANY(Lobc(:,isTvar(itrc),ng))) THEN ifield=idSbry(isTvar(itrc)) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) # ifdef SEDIMENT DO i=1,NST IF (itrc.eq.idsed(i)) THEN WRITE (Vinfo(19),50) 1000.0_r8*Sd50(i,ng) END IF END DO # endif Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(ifield), & & NF_FOUT, 5, t3dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF END DO # endif ! ! Define density anomaly. ! IF (Hout(idDano,ng)) THEN Vinfo( 1)=Vname(1,idDano) WRITE (Vinfo( 2),40) TRIM(Vname(2,idDano)) Vinfo( 3)=Vname(3,idDano) Vinfo(14)=Vname(4,idDano) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idDano) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idDano,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idDano), & & NF_FOUT, nvd4, t3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # if defined FORWARD_MIXING && \ (defined BVF_MIXING || defined GLS_MIXING || \ defined LMD_MIXING || defined MY25_MIXING) ! ! Define vertical viscosity coefficient. ! IF (Hout(idVvis,ng)) THEN Vinfo( 1)=Vname(1,idVvis) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVvis)) Vinfo( 3)=Vname(3,idVvis) Vinfo(14)=Vname(4,idVvis) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idVvis) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVvis,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idVvis), & & NF_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF ! ! Define vertical diffusion coefficient for potential temperature. ! IF (Hout(idTdif,ng)) THEN Vinfo( 1)=Vname(1,idTdif) WRITE (Vinfo( 2),40) TRIM(Vname(2,idTdif)) Vinfo( 3)=Vname(3,idTdif) Vinfo(14)=Vname(4,idTdif) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idTdif) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idTdif,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idTdif), & & NF_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # ifdef SALINITY ! ! Define vertical diffusion coefficient for salinity. ! IF (Hout(idSdif,ng)) THEN Vinfo( 1)=Vname(1,idSdif) WRITE (Vinfo( 2),40) TRIM(Vname(2,idSdif)) Vinfo( 3)=Vname(3,idSdif) Vinfo(14)=Vname(4,idSdif) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idSdif) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idSdif,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idSdif), & & NF_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif # if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET ! ! Define turbulent kinetic energy. ! IF (Hout(idMtke,ng)) THEN Vinfo( 1)=Vname(1,idMtke) WRITE (Vinfo( 2),40) TRIM(Vname(2,idMtke)) Vinfo( 3)=Vname(3,idMtke) Vinfo(14)=Vname(4,idMtke) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idMtke) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idMtke,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idMtke), & & NF_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN Vinfo( 1)=Vname(1,idVmKK) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVmKK)) Vinfo( 3)=Vname(3,idVmKK) Vinfo(14)=Vname(4,idVmKK) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idVmKK) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVmKK,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idVmKK), & & NF_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF ! ! Define turbulent kinetic energy time length scale. ! IF (Hout(idMtls,ng)) THEN Vinfo( 1)=Vname(1,idMtls) WRITE (Vinfo( 2),40) TRIM(Vname(2,idMtls)) Vinfo( 3)=Vname(3,idMtls) Vinfo(14)=Vname(4,idMtls) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idMtls) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idMtls,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idMtls), & & NF_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN Vinfo( 1)=Vname(1,idVmLS) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVmLS)) Vinfo( 3)=Vname(3,idVmLS) Vinfo(14)=Vname(4,idVmLS) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idVmLS) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVmLS,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idVmLS), & & NF_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef GLS_MIXING_NOT_YET Vinfo( 1)=Vname(1,idVmKP) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVmKP)) Vinfo( 3)=Vname(3,idVmKP) Vinfo(14)=Vname(4,idVmKP) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idVmKP) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVmKP,ng),r8) status=def_var(ng, iTLM, TLM(ng)%ncid, TLM(ng)%Vid(idVmKP), & & NF_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif END IF # endif # endif # endif ! !----------------------------------------------------------------------- ! Leave definition mode. !----------------------------------------------------------------------- ! CALL netcdf_enddef (ng, iTLM, ncname, TLM(ng)%ncid) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! !----------------------------------------------------------------------- ! Write out time-recordless, information variables. !----------------------------------------------------------------------- ! CALL wrt_info (ng, iTLM, TLM(ng)%ncid, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF DEFINE ! !======================================================================= ! Open an existing tangent file, check its contents, and prepare for ! appending data. !======================================================================= ! QUERY : IF (.not.ldef) THEN ncname=TLM(ng)%name ! ! Open tangent linear history file for read/write. ! CALL netcdf_open (ng, iTLM, ncname, 1, TLM(ng)%ncid) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) THEN WRITE (stdout,60) TRIM(ncname) RETURN END IF ! ! Inquire about the dimensions and check for consistency. ! CALL netcdf_check_dim (ng, iTLM, ncname, & & ncid = TLM(ng)%ncid) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! ! Inquire about the variables. ! CALL netcdf_inq_var (ng, iTLM, ncname, & & ncid = TLM(ng)%ncid) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! ! Initialize logical switches. ! DO i=1,NV got_var(i)=.FALSE. END DO ! ! Scan variable list from input NetCDF and activate switches for ! tangent variables. Get variable IDs. ! DO i=1,n_var IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idtime))) THEN got_var(idtime)=.TRUE. TLM(ng)%Vid(idtime)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idFsur))) THEN got_var(idFsur)=.TRUE. TLM(ng)%Vid(idFsur)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idUbar))) THEN got_var(idUbar)=.TRUE. TLM(ng)%Vid(idUbar)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVbar))) THEN got_var(idVbar)=.TRUE. TLM(ng)%Vid(idVbar)=var_id(i) # ifdef ADJUST_BOUNDARY ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isFsur)))) THEN got_var(idSbry(isFsur))=.TRUE. TLM(ng)%Vid(idSbry(isFsur))=var_id(i) ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isUbar)))) THEN got_var(idSbry(isUbar))=.TRUE. TLM(ng)%Vid(idSbry(isUbar))=var_id(i) ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isVbar)))) THEN got_var(idSbry(isVbar))=.TRUE. TLM(ng)%Vid(idSbry(isVbar))=var_id(i) # endif # ifdef FORWARD_WRITE # ifdef FORWARD_RHS ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRzet))) THEN got_var(idRzet)=.TRUE. TLM(ng)%Vid(idRzet)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRu2d))) THEN got_var(idRu2d)=.TRUE. TLM(ng)%Vid(idRu2d)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRv2d))) THEN got_var(idRv2d)=.TRUE. TLM(ng)%Vid(idRv2d)=var_id(i) # endif # ifdef SOLVE3D # ifdef FORWARD_RHS ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRuct))) THEN got_var(idRuct)=.TRUE. TLM(ng)%Vid(idRuct)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRvct))) THEN got_var(idRvct)=.TRUE. TLM(ng)%Vid(idRvct)=var_id(i) # endif ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idUfx1))) THEN got_var(idUfx1)=.TRUE. TLM(ng)%Vid(idUfx1)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idUfx2))) THEN got_var(idUfx2)=.TRUE. TLM(ng)%Vid(idUfx2)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVfx1))) THEN got_var(idVfx1)=.TRUE. TLM(ng)%Vid(idVfx1)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVfx2))) THEN got_var(idVfx2)=.TRUE. TLM(ng)%Vid(idVfx2)=var_id(i) # ifdef FORWARD_RHS ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRu3d))) THEN got_var(idRu3d)=.TRUE. TLM(ng)%Vid(idRu3d)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRv3d))) THEN got_var(idRv3d)=.TRUE. TLM(ng)%Vid(idRv3d)=var_id(i) # endif # endif # endif # ifdef SOLVE3D ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idUvel))) THEN got_var(idUvel)=.TRUE. TLM(ng)%Vid(idUvel)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVvel))) THEN got_var(idVvel)=.TRUE. TLM(ng)%Vid(idVvel)=var_id(i) # ifdef ADJUST_BOUNDARY ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isUvel)))) THEN got_var(idSbry(isUvel))=.TRUE. TLM(ng)%Vid(idSbry(isUvel))=var_id(i) ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isVvel)))) THEN got_var(idSbry(isVvel))=.TRUE. TLM(ng)%Vid(idSbry(isVvel))=var_id(i) # endif ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idDano))) THEN got_var(idDano)=.TRUE. TLM(ng)%Vid(idDano)=var_id(i) # if defined FORWARD_MIXING && \ (defined BVF_MIXING || defined GLS_MIXING || \ defined LMD_MIXING || defined MY25_MIXING) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVvis))) THEN got_var(idVvis)=.TRUE. TLM(ng)%Vid(idVvis)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idTdif))) THEN got_var(idTdif)=.TRUE. TLM(ng)%Vid(idTdif)=var_id(i) # ifdef SALINITY ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idSdif))) THEN got_var(idSdif)=.TRUE. TLM(ng)%Vid(idSdif)=var_id(i) # endif # if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idMtke))) THEN got_var(idMtke)=.TRUE. TLM(ng)%Vid(idMtke)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVmKK))) THEN got_var(idVmKK)=.TRUE. TLM(ng)%Vid(idVmKK)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idMtls))) THEN got_var(idMtls)=.TRUE. TLM(ng)%Vid(idMtls)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVmLS))) THEN got_var(idVmLS)=.TRUE. TLM(ng)%Vid(idVmLS)=var_id(i) # ifdef GLS_MIXING_NOT_YET ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVmKP))) THEN got_var(idVmKP)=.TRUE. TLM(ng)%Vid(idVmKP)=var_id(i) # endif # endif # endif # endif # ifdef ADJUST_WSTRESS ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idUsms))) THEN got_var(idUsms)=.TRUE. TLM(ng)%Vid(idUsms)=var_id(i) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVsms))) THEN got_var(idVsms)=.TRUE. TLM(ng)%Vid(idVsms)=var_id(i) # endif END IF # ifdef SOLVE3D DO itrc=1,NT(ng) IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idTvar(itrc)))) THEN got_var(idTvar(itrc))=.TRUE. TLM(ng)%Tid(itrc)=var_id(i) # ifdef ADJUST_BOUNDARY ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isTvar(itrc))))) THEN got_var(idSbry(isTvar(itrc)))=.TRUE. TLM(ng)%Vid(idSbry(isTvar(itrc)))=var_id(i) # endif # ifdef ADJUST_STFLUX ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idTsur(itrc)))) THEN got_var(idTsur(itrc))=.TRUE. TLM(ng)%Vid(idTsur(itrc))=var_id(i) # endif END IF END DO # endif END DO ! ! Check if tangent variables are available in input NetCDF file. ! IF (.not.got_var(idtime)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idtime)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idFsur).and.Hout(idFsur,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idFsur)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idUbar).and.Hout(idUbar,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idUbar)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVbar).and.Hout(idVbar,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVbar)), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef ADJUST_BOUNDARY IF (.not.got_var(idSbry(isFsur))) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSbry(isFsur))), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idSbry(isUbar))) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSbry(isUbar))), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idSbry(isVbar))) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSbry(isVbar))), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # ifdef FORWARD_WRITE # ifdef FORWARD_RHS IF (.not.got_var(idRzet)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRzet)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idRu2d)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRu2d)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idRv2d)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRv2d)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # ifdef SOLVE3D # ifdef FORWARD_RHS IF (.not.got_var(idRuct)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRuct)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idRvct)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRvct)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif IF (.not.got_var(idUfx1)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idUfx1)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idUfx2)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idUfx2)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVfx1)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVfx1)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVfx2)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVfx2)), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef FORWARD_RHS IF (.not.got_var(idRu3d)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRu3d)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idRv3d)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRv3d)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # endif # endif # ifdef SOLVE3D IF (.not.got_var(idUvel).and.Hout(idUvel,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idUvel)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVvel).and.Hout(idVvel,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVvel)), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef ADJUST_BOUNDARY IF (.not.got_var(idSbry(isUvel))) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSbry(isUvel))), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idSbry(isVvel))) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSbry(isVvel))), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif IF (.not.got_var(idDano).and.Hout(idDano,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idDano)), & & TRIM(ncname) exit_flag=3 RETURN END IF # if defined FORWARD_MIXING && \ (defined BVF_MIXING || defined GLS_MIXING || \ defined LMD_MIXING || defined MY25_MIXING) IF (.not.got_var(idVvis).and.Hout(idVvis,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVvis)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idTdif).and.Hout(idTdif,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idTdif)), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef SALINITY IF (.not.got_var(idSdif).and.Hout(idSdif,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSdif)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET IF (.not.got_var(idMtke).and.Hout(idMtke,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idMtke)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVmKK).and.Hout(idVmKK,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVmKK)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idMtls).and.Hout(idMtls,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idMtls)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVmLS).and.Hout(idVmLS,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVmLS)), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef GSL_MIXING IF (.not.got_var(idVmKP).and.Hout(idVmKP,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVmKP)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # endif # endif # endif # ifdef ADJUST_WSTRESS IF (.not.got_var(idUsms).and.Hout(idUsms,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idUsms)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVsms).and.Hout(idVsms,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVsms)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # ifdef SOLVE3D DO itrc=1,NT(ng) IF (.not.got_var(idTvar(itrc)).and.Hout(idTvar(itrc),ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idTvar(itrc))), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef ADJUST_BOUNDARY IF (.not.got_var(idSbry(isTvar(itrc)))) THEN IF (Master) WRITE (stdout,70) & & TRIM(Vname(1,idSbry(isTvar(itrc)))), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # ifdef ADJUST_STFLUX IF (.not.got_var(idTsur(itrc)).and.Hout(idTsur(itrc),ng).and. & & Lstflux(itrc,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idTsur(itrc))), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif END DO # endif ! ! Set unlimited time record dimension to the appropriate value. ! IF (ndefTLM(ng).gt.0) THEN TLM(ng)%Rindex=((ntstart(ng)-1)- & & ndefTLM(ng)*((ntstart(ng)-1)/ndefTLM(ng)))/ & & nTLM(ng) ELSE TLM(ng)%Rindex=(ntstart(ng)-1)/nTLM(ng) END IF TLM(ng)%Rindex=MIN(TLM(ng)%Rindex,rec_size) END IF QUERY ! 10 FORMAT (2x,'TL_DEF_HIS_NF90 - creating tangent file,',t56, & & 'Grid ',i2.2,': ',a) 20 FORMAT (2x,'TL_DEF_HIS_NF90 - inquiring tangent file,',t56, & & 'Grid ',i2.2,': ',a) 30 FORMAT (/,' TL_DEF_HIS_NF90 - unable to create tangent NetCDF', & & ' file: ',a) 40 FORMAT ('tangent linear',1x,a) 50 FORMAT (1pe11.4,1x,'millimeter') 60 FORMAT (/,' TL_DEF_HIS_NF90 - unable to open tangent NetCDF', & & ' file: ',a) 70 FORMAT (/,' TL_DEF_HIS_NF90 - unable to find variable: ',a,2x, & & ' in tangent NetCDF file: ',a) ! RETURN END SUBROUTINE tl_def_his_nf90 # if defined PIO_LIB && defined DISTRIBUTE ! !*********************************************************************** SUBROUTINE tl_def_his_pio (ng, ldef) !*********************************************************************** ! USE mod_pio_netcdf ! ! Imported variable declarations. ! integer, intent(in) :: ng logical, intent(in) :: ldef ! ! Local variable declarations. ! logical :: got_var(NV) ! integer, parameter :: Natt = 25 integer :: i, j, ifield, itrc, nvd3, nvd4 integer :: recdim, status, varid # ifdef ADJUST_BOUNDARY integer :: IorJdim, brecdim # endif # if defined ADJUST_STFLUX || defined ADJUST_WSTRESS integer :: frecdim # endif integer :: DimIDs(nDimID) integer :: t2dgrd(3), u2dgrd(3), v2dgrd(3) # ifdef ADJUST_BOUNDARY integer :: t2dobc(4) # endif # ifdef SOLVE3D integer :: t3dgrd(4), u3dgrd(4), v3dgrd(4), w3dgrd(4) # ifdef ADJUST_BOUNDARY integer :: t3dobc(5) # endif # ifdef ADJUST_STFLUX integer :: t3dfrc(4) # endif # endif # ifdef ADJUST_WSTRESS integer :: u3dfrc(4), v3dfrc(4) # endif ! real(r8) :: Aval(6) ! character (len=256) :: ncname character (len=MaxLen) :: Vinfo(Natt) ! character (len=*), parameter :: MyFile = & & __FILE__//", tl_def_his_pio" ! TYPE (Var_desc_t) :: varDesc ! SourceFile=MyFile ! !----------------------------------------------------------------------- ! Set and report file name. !----------------------------------------------------------------------- ! IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ncname=TLM(ng)%name ! IF (Master) THEN IF (ldef) THEN WRITE (stdout,10) ng, TRIM(ncname) ELSE WRITE (stdout,20) ng, TRIM(ncname) END IF END IF ! !======================================================================= ! Create a new tangent linear history file. !======================================================================= ! DEFINE : IF (ldef) THEN CALL pio_netcdf_create (ng, iTLM, TRIM(ncname), TLM(ng)%pioFile) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) THEN IF (Master) WRITE (stdout,30) TRIM(ncname) RETURN END IF ! !----------------------------------------------------------------------- ! Define file dimensions. !----------------------------------------------------------------------- ! DimIDs=0 ! status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xi_rho', & & IOBOUNDS(ng)%xi_rho, DimIDs( 1)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xi_u', & & IOBOUNDS(ng)%xi_u, DimIDs( 2)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xi_v', & & IOBOUNDS(ng)%xi_v, DimIDs( 3)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xi_psi', & & IOBOUNDS(ng)%xi_psi, DimIDs( 4)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'eta_rho', & & IOBOUNDS(ng)%eta_rho, DimIDs( 5)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'eta_u', & & IOBOUNDS(ng)%eta_u, DimIDs( 6)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'eta_v', & & IOBOUNDS(ng)%eta_v, DimIDs( 7)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'eta_psi', & & IOBOUNDS(ng)%eta_psi, DimIDs( 8)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef ADJUST_BOUNDARY status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'IorJ', & & IOBOUNDS(ng)%IorJ, IorJdim) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # if defined WRITE_WATER && defined MASKING status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xy_rho', & & IOBOUNDS(ng)%xy_rho, DimIDs(17)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xy_u', & & IOBOUNDS(ng)%xy_u, DimIDs(18)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xy_v', & & IOBOUNDS(ng)%xy_v, DimIDs(19)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # ifdef SOLVE3D # if defined WRITE_WATER && defined MASKING status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xyz_rho', & & IOBOUNDS(ng)%xy_rho*N(ng), DimIDs(20)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xyz_u', & & IOBOUNDS(ng)%xy_u*N(ng), DimIDs(21)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xyz_v', & & IOBOUNDS(ng)%xy_v*N(ng), DimIDs(22)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xyz_w', & & IOBOUNDS(ng)%xy_rho*(N(ng)+1), DimIDs(23)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'N', & & N(ng), DimIDs( 9)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 's_rho', & & N(ng), DimIDs( 9)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 's_w', & & N(ng)+1, DimIDs(10)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'tracer', & & NT(ng), DimIDs(11)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef SEDIMENT status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'NST', & & NST, DimIDs(32)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'Nbed', & & Nbed, DimIDs(16)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # if defined WRITE_WATER && defined MASKING status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'xybed', & & IOBOUNDS(ng)%xy_rho*Nbed, DimIDs(24)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # endif # ifdef ECOSIM status=def_dim(ng, iNLM, TLM(ng)%pioFile, ncname, 'Nbands', & & NBands, DimIDs(33)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'Nphy', & & Nphy, DimIDs(25)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'Nbac', & & Nbac, DimIDs(26)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'Ndom', & & Ndom, DimIDs(27)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'Nfec', & & Nfec, DimIDs(28)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # endif status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'boundary', & & 4, DimIDs(14)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef FOUR_DVAR status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'Nstate', & & NstateVar(ng), DimIDs(29)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # if defined ADJUST_STFLUX || defined ADJUST_WSTRESS status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'frc_adjust', & & Nfrec(ng), DimIDs(30)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # ifdef ADJUST_BOUNDARY status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, 'obc_adjust', & & Nbrec(ng), DimIDs(31)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif status=def_dim(ng, iTLM, TLM(ng)%pioFile, ncname, & & TRIM(ADJUSTL(Vname(5,idtime))), & & nf90_unlimited, DimIDs(12)) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN recdim=DimIDs(12) # if defined ADJUST_STFLUX || defined ADJUST_WSTRESS frecdim=DimIDs(30) # endif # ifdef ADJUST_BOUNDARY brecdim=DimIDs(31) # endif ! ! Set number of dimensions for output variables. ! # if defined WRITE_WATER && defined MASKING nvd3=2 nvd4=2 # else nvd3=3 nvd4=4 # endif ! ! Define dimension vectors for staggered tracer type variables. ! # if defined WRITE_WATER && defined MASKING t2dgrd(1)=DimIDs(17) t2dgrd(2)=DimIDs(12) # ifdef SOLVE3D t3dgrd(1)=DimIDs(20) t3dgrd(2)=DimIDs(12) # endif # else t2dgrd(1)=DimIDs( 1) t2dgrd(2)=DimIDs( 5) t2dgrd(3)=DimIDs(12) # ifdef SOLVE3D t3dgrd(1)=DimIDs( 1) t3dgrd(2)=DimIDs( 5) t3dgrd(3)=DimIDs( 9) t3dgrd(4)=DimIDs(12) # endif # ifdef ADJUST_STFLUX t3dfrc(1)=DimIDs( 1) t3dfrc(2)=DimIDs( 5) t3dfrc(3)=frecdim t3dfrc(4)=DimIDs(12) # endif # endif # ifdef ADJUST_BOUNDARY t2dobc(1)=IorJdim t2dobc(2)=DimIDs(14) t2dobc(3)=brecdim t2dobc(4)=DimIDs(12) # ifdef SOLVE3D t3dobc(1)=IorJdim t3dobc(2)=DimIDs( 9) t3dobc(3)=DimIDs(14) t3dobc(4)=brecdim t3dobc(5)=DimIDs(12) # endif # endif ! ! Define dimension vectors for staggered u-momentum type variables. ! # if defined WRITE_WATER && defined MASKING u2dgrd(1)=DimIDs(18) u2dgrd(2)=DimIDs(12) # ifdef SOLVE3D u3dgrd(1)=DimIDs(21) u3dgrd(2)=DimIDs(12) # endif # else u2dgrd(1)=DimIDs( 2) u2dgrd(2)=DimIDs( 6) u2dgrd(3)=DimIDs(12) # ifdef SOLVE3D u3dgrd(1)=DimIDs( 2) u3dgrd(2)=DimIDs( 6) u3dgrd(3)=DimIDs( 9) u3dgrd(4)=DimIDs(12) # endif # ifdef ADJUST_WSTRESS u3dfrc(1)=DimIDs( 2) u3dfrc(2)=DimIDs( 6) u3dfrc(3)=frecdim u3dfrc(4)=DimIDs(12) # endif # endif ! ! Define dimension vectors for staggered v-momentum type variables. ! # if defined WRITE_WATER && defined MASKING v2dgrd(1)=DimIDs(19) v2dgrd(2)=DimIDs(12) # ifdef SOLVE3D v3dgrd(1)=DimIDs(22) v3dgrd(2)=DimIDs(12) # endif # else v2dgrd(1)=DimIDs( 3) v2dgrd(2)=DimIDs( 7) v2dgrd(3)=DimIDs(12) # ifdef SOLVE3D v3dgrd(1)=DimIDs( 3) v3dgrd(2)=DimIDs( 7) v3dgrd(3)=DimIDs( 9) v3dgrd(4)=DimIDs(12) # endif # ifdef ADJUST_WSTRESS v3dfrc(1)=DimIDs( 3) v3dfrc(2)=DimIDs( 7) v3dfrc(3)=frecdim v3dfrc(4)=DimIDs(12) # endif # endif # ifdef SOLVE3D ! ! Define dimension vector for staggered w-momentum type variables. ! # if defined WRITE_WATER && defined MASKING w3dgrd(1)=DimIDs(23) w3dgrd(2)=DimIDs(12) # else w3dgrd(1)=DimIDs( 1) w3dgrd(2)=DimIDs( 5) w3dgrd(3)=DimIDs(10) w3dgrd(4)=DimIDs(12) # endif # endif ! ! Initialize unlimited time record dimension. ! TLM(ng)%Rindex=0 ! ! Initialize local information variable arrays. ! DO i=1,Natt DO j=1,LEN(Vinfo(1)) Vinfo(i)(j:j)=' ' END DO END DO DO i=1,6 Aval(i)=0.0_r8 END DO ! !----------------------------------------------------------------------- ! Define time-recordless information variables. !----------------------------------------------------------------------- ! CALL def_info (ng, iTLM, TLM(ng)%pioFile, ncname, DimIDs) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! !----------------------------------------------------------------------- ! Define time-varying variables. !----------------------------------------------------------------------- ! ! Define model time. ! Vinfo( 1)=Vname(1,idtime) Vinfo( 2)=Vname(2,idtime) WRITE (Vinfo( 3),'(a,a)') 'seconds since ', TRIM(Rclock%string) Vinfo( 4)=TRIM(Rclock%calendar) Vinfo(14)=Vname(4,idtime) Vinfo(21)=Vname(6,idtime) TLM(ng)%pioVar(idtime)%dkind=PIO_TOUT TLM(ng)%pioVar(idtime)%gtype=0 ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idtime)%vd, PIO_TOUT, & & 1, (/recdim/), Aval, Vinfo, ncname, & & SetParAccess = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef PROPAGATOR ! ! Define Ritz eigenvalues and Ritz eigenvectors Euclidean norm. ! Vinfo( 1)='Ritz_rvalue' Vinfo( 2)='real Ritz eigenvalues' status=def_var(ng, iTLM, TLM(ng)%pioFile, varDesc, PIO_TYPE, & & 1, (/recdim/), Aval, Vinfo, ncname, & & SetParAccess = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! # if defined FT_EIGENMODES Vinfo( 1)='Ritz_ivalue' Vinfo( 2)='imaginary Ritz eigenvalues' status=def_var(ng, iTLM, TLM(ng)%pioFile, varDesc, PIO_TYPE, & & 1, (/recdim/), Aval, Vinfo, ncname, & & SetParAccess = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! # endif Vinfo( 1)='Ritz_norm' Vinfo( 2)='Ritz eigenvectors Euclidean norm' status=def_var(ng, iTLM, TLM(ng)%pioFile, varDesc, PIO_TYPE, & & 1, (/recdim/), Aval, Vinfo, ncname, & & SetParAccess = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # ifdef ADJUST_WSTRESS ! ! Define surface U-momentum stress. Notice that the stress has its ! own fixed time-dimension (of size Nfrec) to allow 4DVAR adjustments ! at other times in addition to initialization time. ! Vinfo( 1)=Vname(1,idUsms) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUsms)) Vinfo( 3)='meter2 second-2' Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idUsms,ng),r8) TLM(ng)%pioVar(idUsms)%dkind=PIO_FOUT TLM(ng)%pioVar(idUsms)%gtype=u2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idUsms)%vd, & & PIO_FOUT, nvd4, u3dfrc, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! ! Define surface V-momentum stress. ! Vinfo( 1)=Vname(1,idVsms) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVsms)) Vinfo( 3)='meter2 second-2' Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVsms,ng),r8) TLM(ng)%pioVar(idVsms)%dkind=PIO_FOUT TLM(ng)%pioVar(idVsms)%gtype=v2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idVsms)%vd, & & PIO_FOUT, nvd4, v3dfrc, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # if defined FORCING_SV || defined STOCHASTIC_OPT || \ defined HESSIAN_SO || defined HESSIAN_FSV ! ! Define surface U-momentum stress. ! IF (Hout(idUsms,ng)) THEN Vinfo( 1)=Vname(1,idUsms) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUsms)) Vinfo( 3)=Vname(3,idUsms) Vinfo(14)=Vname(4,idUsms) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idUsms) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idUsms,ng),r8) TLM(ng)%pioVar(idUsms)%dkind=PIO_FOUT TLM(ng)%pioVar(idUsms)%gtype=u2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idUsms)%vd, & & PIO_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF ! ! Define surface V-momentum stress. ! IF (Hout(idVsms,ng)) THEN Vinfo( 1)=Vname(1,idVsms) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVsms)) Vinfo( 3)=Vname(3,idVsms) Vinfo(14)=Vname(4,idVsms) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idVsms) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVsms,ng),r8) TLM(ng)%pioVar(idVsms)%dkind=PIO_FOUT TLM(ng)%pioVar(idVsms)%gtype=v2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idVsms)%vd, & & PIO_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF ! ! Define surface tracer fluxes. ! DO itrc=1,NT(ng) IF (Hout(idTsur(itrc),ng)) THEN Vinfo( 1)=Vname(1,idTsur(itrc)) WRITE (Vinfo( 2),40) TRIM(Vname(2,idTsur(itrc))) Vinfo( 3)=Vname(3,idTsur(itrc)) IF (itrc.eq.itemp) THEN Vinfo(11)='upward flux, cooling' Vinfo(12)='downward flux, heating' ELSE IF (itrc.eq.isalt) THEN Vinfo(11)='upward flux, freshening (net precipitation)' Vinfo(12)='downward flux, salting (net evaporation)' END IF Vinfo(14)=Vname(4,idTsur(itrc)) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idTsur(itrc)) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idTsur(itrc),ng),r8) TLM(ng)%pioVar(idTsur(itrc))%dkind=PIO_FOUT TLM(ng)%pioVar(idTsur(itrc))%gtype=r2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idTsur(itrc))%vd, & & PIO_FOUT, nvd3, t2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF END DO # endif # if defined ADJUST_STFLUX && defined SOLVE3D ! ! Define surface net heat flux. Notice that different tracer fluxes ! are written at their own fixed time-dimension (of size Nfrec) to ! allow 4DVAR adjustments at other times in addition to initial time. ! DO itrc=1,NT(ng) IF (Lstflux(itrc,ng)) THEN Vinfo( 1)=Vname(1,idTsur(itrc)) WRITE (Vinfo( 2),40) TRIM(Vname(2,idTsur(itrc))) IF (itrc.eq.itemp) THEN Vinfo( 3)='Celsius meter second-1' Vinfo(11)='upward flux, cooling' Vinfo(12)='downward flux, heating' ELSE IF (itrc.eq.isalt) THEN Vinfo( 3)='meter second-1' Vinfo(11)='upward flux, freshening (net precipitation)' Vinfo(12)='downward flux, salting (net evaporation)' END IF Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idTsur(itrc),ng),r8) TLM(ng)%pioVar(idTsur(itrc))%dkind=PIO_FOUT TLM(ng)%pioVar(idTsur(itrc))%gtype=r2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idTsur(itrc))%vd, & & PIO_FOUT, nvd4, t3dfrc, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF END DO # endif ! ! Define free-surface. ! IF (Hout(idFsur,ng)) THEN Vinfo( 1)=Vname(1,idFsur) WRITE (Vinfo( 2),40) TRIM(Vname(2,idFsur)) Vinfo( 3)=Vname(3,idFsur) Vinfo(14)=Vname(4,idFsur) Vinfo(16)=Vname(1,idtime) # if !defined WET_DRY && (defined WRITE_WATER && defined MASKING) Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idFsur) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idFsur,ng),r8) TLM(ng)%pioVar(idFsur)%dkind=PIO_FOUT TLM(ng)%pioVar(idFsur)%gtype=r2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idFsur)%vd, & # ifdef WET_DRY & PIO_FOUT, nvd3, t2dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) # else & PIO_FOUT, nvd3, t2dgrd, Aval, Vinfo, ncname) # endif IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # if defined FORWARD_WRITE && defined FORWARD_RHS ! Vinfo( 1)=Vname(1,idRzet) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRzet)) Vinfo( 3)=Vname(3,idRzet) Vinfo(14)=Vname(4,idRzet) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idRzet) Vinfo(22)='coordinates' Aval(5)=REAL(r2dvar,r8) TLM(ng)%pioVar(idRzet)%dkind=PIO_FOUT TLM(ng)%pioVar(idRzet)%gtype=r2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idRzet)%vd, & & PIO_FOUT, nvd3, t2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif END IF # ifdef ADJUST_BOUNDARY ! ! Define free-surface open boundaries. ! IF (ANY(Lobc(:,isFsur,ng))) THEN ifield=idSbry(isFsur) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) TLM(ng)%pioVar(ifield)%dkind=PIO_FOUT TLM(ng)%pioVar(ifield)%gtype=r2dobc ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(ifield)%vd, & & PIO_FOUT, 4, t2dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif ! ! Define 2D U-momentum component. ! IF (Hout(idUbar,ng)) THEN Vinfo( 1)=Vname(1,idUbar) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUbar)) Vinfo( 3)=Vname(3,idUbar) Vinfo(14)=Vname(4,idUbar) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idUbar) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idUbar,ng),r8) TLM(ng)%pioVar(idUbar)%dkind=PIO_FOUT TLM(ng)%pioVar(idUbar)%gtype=u2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idUbar)%vd, & & PIO_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef FORWARD_WRITE # ifdef FORWARD_RHS ! Vinfo( 1)=Vname(1,idRu2d) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRu2d)) Vinfo( 3)=Vname(3,idRu2d) Vinfo(14)=Vname(4,idRu2d) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idRu2d) Vinfo(22)='coordinates' Aval(5)=REAL(u2dvar,r8) TLM(ng)%pioVar(idRu2d)%dkind=PIO_FOUT TLM(ng)%pioVar(idRu2d)%gtype=u2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idRu2d)%vd, & & PIO_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # ifdef SOLVE3D # ifdef FORWARD_RHS ! Vinfo( 1)=Vname(1,idRuct) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRuct)) Vinfo( 3)=Vname(3,idRuct) Vinfo(14)=Vname(4,idRuct) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idRuct) Vinfo(22)='coordinates' Aval(5)=REAL(u2dvar,r8) TLM(ng)%pioVar(idRuct)%dkind=PIO_FOUT TLM(ng)%pioVar(idRuct)%gtype=u2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idRuct)%vd, & & PIO_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif ! Vinfo( 1)=Vname(1,idUfx1) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUfx1)) Vinfo( 3)=Vname(3,idUfx1) Vinfo(14)=Vname(4,idUfx1) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idUfx1) Vinfo(22)='coordinates' Aval(5)=REAL(u2dvar,r8) TLM(ng)%pioVar(idUfx1)%dkind=PIO_FOUT TLM(ng)%pioVar(idUfx1)%gtype=u2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idUfx1)%vd, & & PIO_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! Vinfo( 1)=Vname(1,idUfx2) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUfx2)) Vinfo( 3)=Vname(3,idUfx2) Vinfo(14)=Vname(4,idUfx2) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idUfx2) Vinfo(22)='coordinates' Aval(5)=REAL(u2dvar,r8) TLM(ng)%pioVar(idUfx2)%dkind=PIO_FOUT TLM(ng)%pioVar(idUfx2)%gtype=u2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idUfx2)%vd, & & PIO_FOUT, nvd3, u2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # endif END IF # ifdef ADJUST_BOUNDARY ! ! Define 2D U-momentum component open boundaries. ! IF (ANY(Lobc(:,isUbar,ng))) THEN ifield=idSbry(isUbar) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) TLM(ng)%pioVar(ifield)%dkind=PIO_FOUT TLM(ng)%pioVar(ifield)%gtype=u2dobc ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(ifield)%vd, & & PIO_FOUT, 4, t2dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif ! ! Define 2D V-momentum component. ! IF (Hout(idVbar,ng)) THEN Vinfo( 1)=Vname(1,idVbar) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVbar)) Vinfo( 3)=Vname(3,idVbar) Vinfo(14)=Vname(4,idVbar) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idVbar) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVbar,ng),r8) TLM(ng)%pioVar(idVbar)%dkind=PIO_FOUT TLM(ng)%pioVar(idVbar)%gtype=v2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idVbar)%vd, & & PIO_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef FORWARD_WRITE # ifdef FORWARD_RHS ! Vinfo( 1)=Vname(1,idRv2d) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRv2d)) Vinfo( 3)=Vname(3,idRv2d) Vinfo(14)=Vname(4,idRv2d) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idRv2d) Vinfo(22)='coordinates' Aval(5)=REAL(v2dvar,r8) TLM(ng)%pioVar(idRv2d)%dkind=PIO_FOUT TLM(ng)%pioVar(idRv2d)%gtype=v2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idRv2d)%vd, & & PIO_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # ifdef SOLVE3D # ifdef FORWARD_RHS ! Vinfo( 1)=Vname(1,idRvct) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRvct)) Vinfo( 3)=Vname(3,idRvct) Vinfo(14)=Vname(4,idRvct) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idRvct) Vinfo(22)='coordinates' Aval(5)=REAL(v2dvar,r8) TLM(ng)%pioVar(idRvct)%dkind=PIO_FOUT TLM(ng)%pioVar(idRvct)%gtype=v2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idRvct)%vd, & & PIO_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif ! Vinfo( 1)=Vname(1,idVfx1) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVfx1)) Vinfo( 3)=Vname(3,idVfx1) Vinfo(14)=Vname(4,idVfx1) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idVfx1) Vinfo(22)='coordinates' Aval(5)=REAL(v2dvar,r8) TLM(ng)%pioVar(idVfx1)%dkind=PIO_FOUT TLM(ng)%pioVar(idVfx1)%gtype=v2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idVfx1)%vd, & & PIO_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! Vinfo( 1)=Vname(1,idVfx2) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVfx2)) Vinfo( 3)=Vname(3,idVfx2) Vinfo(14)=Vname(4,idVfx2) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idVfx2) Vinfo(22)='coordinates' Aval(5)=REAL(v2dvar,r8) TLM(ng)%pioVar(idVfx2)%dkind=PIO_FOUT TLM(ng)%pioVar(idVfx2)%gtype=v2dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idVfx2)%vd, & & PIO_FOUT, nvd3, v2dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif # endif END IF # ifdef ADJUST_BOUNDARY ! ! Define 2D V-momentum component open boundaries. ! IF (ANY(Lobc(:,isVbar,ng))) THEN ifield=idSbry(isVbar) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) TLM(ng)%pioVar(ifield)%dkind=PIO_FOUT TLM(ng)%pioVar(ifield)%gtype=v2dobc ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(ifield)%vd, & & PIO_FOUT, 4, t2dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif # ifdef SOLVE3D ! ! Define 3D U-momentum component. ! IF (Hout(idUvel,ng)) THEN Vinfo( 1)=Vname(1,idUvel) WRITE (Vinfo( 2),40) TRIM(Vname(2,idUvel)) Vinfo( 3)=Vname(3,idUvel) Vinfo(14)=Vname(4,idUvel) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idUvel) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idUvel,ng),r8) TLM(ng)%pioVar(idUvel)%dkind=PIO_FOUT TLM(ng)%pioVar(idUvel)%gtype=u3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idUvel)%vd, & & PIO_FOUT, nvd4, u3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # if defined FORWARD_WRITE && defined FORWARD_RHS ! Vinfo( 1)=Vname(1,idRu3d) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRu3d)) Vinfo( 3)=Vname(3,idRu3d) Vinfo(14)=Vname(4,idRu3d) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_u' # endif Vinfo(21)=Vname(6,idRu3d) Vinfo(22)='coordinates' Aval(5)=REAL(u3dvar,r8) TLM(ng)%pioVar(idRu3d)%dkind=PIO_FOUT TLM(ng)%pioVar(idRu3d)%gtype=u3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idRu3d)%vd, & & PIO_FOUT, nvd4, u3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif END IF # ifdef ADJUST_BOUNDARY ! ! Define 3D U-momentum component open boundaries. ! IF (ANY(Lobc(:,isUvel,ng))) THEN ifield=idSbry(isUvel) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) TLM(ng)%pioVar(ifield)%dkind=PIO_FOUT TLM(ng)%pioVar(ifield)%gtype=u3dobc ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(ifield)%vd, & & PIO_FOUT, 5, t3dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif ! ! Define 3D V-momentum component. ! IF (Hout(idVvel,ng)) THEN Vinfo( 1)=Vname(1,idVvel) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVvel)) Vinfo( 3)=Vname(3,idVvel) Vinfo(14)=Vname(4,idVvel) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idVvel) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVvel,ng),r8) TLM(ng)%pioVar(idVvel)%dkind=PIO_FOUT TLM(ng)%pioVar(idVvel)%gtype=v3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idVvel)%vd, & & PIO_FOUT, nvd4, v3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # if defined FORWARD_WRITE && defined FORWARD_RHS ! Vinfo( 1)=Vname(1,idRv3d) WRITE (Vinfo( 2),40) TRIM(Vname(2,idRv3d)) Vinfo( 3)=Vname(3,idRv3d) Vinfo(14)=Vname(4,idRv3d) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_v' # endif Vinfo(21)=Vname(6,idRv3d) Vinfo(22)='coordinates' Aval(5)=REAL(v3dvar,r8) TLM(ng)%pioVar(idRv3d)%dkind=PIO_FOUT TLM(ng)%pioVar(idRv3d)%gtype=v3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idRv3d)%vd, & & PIO_FOUT, nvd4, v3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif END IF # ifdef ADJUST_BOUNDARY ! ! Define 3D V-momentum component open boundaries. ! IF (ANY(Lobc(:,isVvel,ng))) THEN ifield=idSbry(isVvel) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) TLM(ng)%pioVar(ifield)%dkind=PIO_FOUT TLM(ng)%pioVar(ifield)%gtype=v3dobc ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(ifield)%vd, & & PIO_FOUT, 5, t3dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif ! ! Define tracer type variables. ! DO itrc=1,NT(ng) IF (Hout(idTvar(itrc),ng)) THEN Vinfo( 1)=Vname(1,idTvar(itrc)) WRITE (Vinfo( 2),40) TRIM(Vname(2,idTvar(itrc))) Vinfo( 3)=Vname(3,idTvar(itrc)) Vinfo(14)=Vname(4,idTvar(itrc)) Vinfo(16)=Vname(1,idtime) # ifdef SEDIMENT_NOT_YET DO i=1,NST IF (itrc.eq.idsed(i)) THEN WRITE (Vinfo(19),50) 1000.0_r8*Sd50(i,ng) END IF END DO # endif # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idTvar(itrc)) Vinfo(22)='coordinates' Aval(5)=REAL(r3dvar,r8) TLM(ng)%pioTrc(itrc)%dkind=PIO_FOUT TLM(ng)%pioTrc(itrc)%gtype=r3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioTrc(itrc)%vd, & & PIO_FOUT, nvd4, t3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF END DO # ifdef ADJUST_BOUNDARY ! ! Define tracer type variables open boundaries. ! DO itrc=1,NT(ng) IF (ANY(Lobc(:,isTvar(itrc),ng))) THEN ifield=idSbry(isTvar(itrc)) Vinfo( 1)=Vname(1,ifield) WRITE (Vinfo( 2),40) TRIM(Vname(2,ifield)) Vinfo( 3)=Vname(3,ifield) Vinfo(14)=Vname(4,ifield) Vinfo(16)=Vname(1,idtime) # ifdef SEDIMENT DO i=1,NST IF (itrc.eq.idsed(i)) THEN WRITE (Vinfo(19),50) 1000.0_r8*Sd50(i,ng) END IF END DO # endif Vinfo(21)=Vname(6,ifield) Aval(5)=REAL(Iinfo(1,ifield,ng),r8) TLM(ng)%pioVar(ifield)%dkind=PIO_FOUT TLM(ng)%pioVar(ifield)%gtype=r3dobc ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(ifield)%vd, & & PIO_FOUT, 5, t3dobc, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF END DO # endif ! ! Define density anomaly. ! IF (Hout(idDano,ng)) THEN Vinfo( 1)=Vname(1,idDano) WRITE (Vinfo( 2),40) TRIM(Vname(2,idDano)) Vinfo( 3)=Vname(3,idDano) Vinfo(14)=Vname(4,idDano) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idDano) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idDano,ng),r8) TLM(ng)%pioVar(idDano)%dkind=PIO_FOUT TLM(ng)%pioVar(idDano)%gtype=r3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idDano)%vd, & & PIO_FOUT, nvd4, t3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # if defined FORWARD_MIXING && \ (defined BVF_MIXING || defined GLS_MIXING || \ defined LMD_MIXING || defined MY25_MIXING) ! ! Define vertical viscosity coefficient. ! IF (Hout(idVvis,ng)) THEN Vinfo( 1)=Vname(1,idVvis) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVvis)) Vinfo( 3)=Vname(3,idVvis) Vinfo(14)=Vname(4,idVvis) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idVvis) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVvis,ng),r8) TLM(ng)%pioVar(idVvis)%dkind=PIO_FOUT TLM(ng)%pioVar(idVvis)%gtype=w3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idVvis)%vd, & & PIO_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF ! ! Define vertical diffusion coefficient for potential temperature. ! IF (Hout(idTdif,ng)) THEN Vinfo( 1)=Vname(1,idTdif) WRITE (Vinfo( 2),40) TRIM(Vname(2,idTdif)) Vinfo( 3)=Vname(3,idTdif) Vinfo(14)=Vname(4,idTdif) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idTdif) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idTdif,ng),r8) TLM(ng)%pioVar(idTdif)%dkind=PIO_FOUT TLM(ng)%pioVar(idTdif)%gtype=w3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idTdif)%vd, & & PIO_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # ifdef SALINITY ! ! Define vertical diffusion coefficient for salinity. ! IF (Hout(idSdif,ng)) THEN Vinfo( 1)=Vname(1,idSdif) WRITE (Vinfo( 2),40) TRIM(Vname(2,idSdif)) Vinfo( 3)=Vname(3,idSdif) Vinfo(14)=Vname(4,idSdif) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idSdif) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idSdif,ng),r8) TLM(ng)%pioVar(idSdif)%dkind=PIO_FOUT TLM(ng)%pioVar(idSdif)%gtype=w3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idSdif)%vd, & & PIO_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF # endif # if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET ! ! Define turbulent kinetic energy. ! IF (Hout(idMtke,ng)) THEN Vinfo( 1)=Vname(1,idMtke) WRITE (Vinfo( 2),40) TRIM(Vname(2,idMtke)) Vinfo( 3)=Vname(3,idMtke) Vinfo(14)=Vname(4,idMtke) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idMtke) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idMtke,ng),r8) TLM(ng)%pioVar(idMtke)%dkind=PIO_FOUT TLM(ng)%pioVar(idMtke)%gtype=w3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idMtke)%vd, & & PIO_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! Vinfo( 1)=Vname(1,idVmKK) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVmKK)) Vinfo( 3)=Vname(3,idVmKK) Vinfo(14)=Vname(4,idVmKK) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idVmKK) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVmKK,ng),r8) TLM(ng)%pioVar(idVmKK)%dkind=PIO_FOUT TLM(ng)%pioVar(idVmKK)%gtype=w3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idVmKK)%vd, & & PIO_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF ! ! Define turbulent kinetic energy time length scale. ! IF (Hout(idMtls,ng)) THEN Vinfo( 1)=Vname(1,idMtls) WRITE (Vinfo( 2),40) TRIM(Vname(2,idMtls)) Vinfo( 3)=Vname(3,idMtls) Vinfo(14)=Vname(4,idMtls) Vinfo(16)=Vname(1,idtime) Vinfo(21)=Vname(6,idMtls) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idMtls,ng),r8) TLM(ng)%pioVar(idMtls)%dkind=PIO_FOUT TLM(ng)%pioVar(idMtls)%gtype=w3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idMtls)%vd, & & PIO_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname, & & SetFillVal = .FALSE.) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! Vinfo( 1)=Vname(1,idVmLS) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVmLS)) Vinfo( 3)=Vname(3,idVmLS) Vinfo(14)=Vname(4,idVmLS) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idVmLS) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVmLS,ng),r8) TLM(ng)%pioVar(idVmLS)%dkind=PIO_FOUT TLM(ng)%pioVar(idVmLS)%gtype=w3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, & & TLM(ng)%pioVar(idVmLS)%vd, & & PIO_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # ifdef GLS_MIXING_NOT_YET ! Vinfo( 1)=Vname(1,idVmKP) WRITE (Vinfo( 2),40) TRIM(Vname(2,idVmKP)) Vinfo( 3)=Vname(3,idVmKP) Vinfo(14)=Vname(4,idVmKP) Vinfo(16)=Vname(1,idtime) # if defined WRITE_WATER && defined MASKING Vinfo(20)='mask_rho' # endif Vinfo(21)=Vname(6,idVmKP) Vinfo(22)='coordinates' Aval(5)=REAL(Iinfo(1,idVmKP,ng),r8) TLM(ng)%pioVar(idVmKP)%dkind=PIO_FOUT TLM(ng)%pioVar(idVmKP)%gtype=w3dvar ! status=def_var(ng, iTLM, TLM(ng)%pioFile, TLM(ng)%pioVar(idVmKP), & & PIO_FOUT, nvd4, w3dgrd, Aval, Vinfo, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN # endif END IF # endif # endif # endif ! !----------------------------------------------------------------------- ! Leave definition mode. !----------------------------------------------------------------------- ! CALL pio_netcdf_enddef (ng, iTLM, ncname, TLM(ng)%pioFile) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! !----------------------------------------------------------------------- ! Write out time-recordless, information variables. !----------------------------------------------------------------------- ! CALL wrt_info (ng, iTLM, TLM(ng)%pioFile, ncname) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN END IF DEFINE ! !======================================================================= ! Open an existing tangent file, check its contents, and prepare for ! appending data. !======================================================================= ! QUERY : IF (.not.ldef) THEN ncname=TLM(ng)%name ! ! Open tangent linear history file for read/write. ! CALL pio_netcdf_open (ng, iTLM, ncname, 1, TLM(ng)%pioFile) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) THEN WRITE (stdout,60) TRIM(ncname) RETURN END IF ! ! Inquire about the dimensions and check for consistency. ! CALL pio_netcdf_check_dim (ng, iTLM, ncname, & & pioFile = TLM(ng)%pioFile) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! ! Inquire about the variables. ! CALL pio_netcdf_inq_var (ng, iTLM, ncname, & & pioFile = TLM(ng)%pioFile) IF (FoundError(exit_flag, NoError, __LINE__, MyFile)) RETURN ! ! Initialize logical switches. ! DO i=1,NV got_var(i)=.FALSE. END DO ! ! Scan variable list from input NetCDF and activate switches for ! tangent variables. Get variable IDs. ! DO i=1,n_var IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idtime))) THEN got_var(idtime)=.TRUE. TLM(ng)%pioVar(idtime)%vd=var_desc(i) TLM(ng)%pioVar(idtime)%dkind=PIO_TOUT TLM(ng)%pioVar(idtime)%gtype=0 ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idFsur))) THEN got_var(idFsur)=.TRUE. TLM(ng)%pioVar(idFsur)%vd=var_desc(i) TLM(ng)%pioVar(idFsur)%dkind=PIO_FOUT TLM(ng)%pioVar(idFsur)%gtype=r2dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idUbar))) THEN got_var(idUbar)=.TRUE. TLM(ng)%pioVar(idUbar)%vd=var_desc(i) TLM(ng)%pioVar(idUbar)%dkind=PIO_FOUT TLM(ng)%pioVar(idUbar)%gtype=u2dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVbar))) THEN got_var(idVbar)=.TRUE. TLM(ng)%pioVar(idVbar)%vd=var_desc(i) TLM(ng)%pioVar(idVbar)%dkind=PIO_FOUT TLM(ng)%pioVar(idVbar)%gtype=v2dvar # ifdef ADJUST_BOUNDARY ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isFsur)))) THEN got_var(idSbry(isFsur))=.TRUE. TLM(ng)%pioVar(idSbry(isFsur))%vd=var_desc(i) TLM(ng)%pioVar(idSbry(isFsur))%dkind=PIO_FOUT TLM(ng)%pioVar(idSbry(isFsur))%gtype=r2dobc ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isUbar)))) THEN got_var(idSbry(isUbar))=.TRUE. TLM(ng)%pioVar(idSbry(isUbar))%vd=var_desc(i) TLM(ng)%pioVar(idSbry(isUbar))%dkind=PIO_FOUT TLM(ng)%pioVar(idSbry(isUbar))%gtype=u2dobc ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isVbar)))) THEN got_var(idSbry(isVbar))=.TRUE. TLM(ng)%pioVar(idSbry(isVbar))%vd=var_desc(i) TLM(ng)%pioVar(idSbry(isVbar))%dkind=PIO_FOUT TLM(ng)%pioVar(idSbry(isVbar))%gtype=v2dobc # endif # ifdef FORWARD_WRITE # ifdef FORWARD_RHS ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRzet))) THEN got_var(idRzet)=.TRUE. TLM(ng)%pioVar(idRzet)%vd=var_desc(i) TLM(ng)%pioVar(idRzet)%dkind=PIO_FOUT TLM(ng)%pioVar(idRzet)%gtype=r2dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRu2d))) THEN got_var(idRu2d)=.TRUE. TLM(ng)%pioVar(idRu2d)%vd=var_desc(i) TLM(ng)%pioVar(idRu2d)%dkind=PIO_FOUT TLM(ng)%pioVar(idRu2d)%gtype=u2dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRv2d))) THEN got_var(idRv2d)=.TRUE. TLM(ng)%pioVar(idRv2d)%vd=var_desc(i) TLM(ng)%pioVar(idRv2d)%dkind=PIO_FOUT TLM(ng)%pioVar(idRv2d)%gtype=v2dvar # endif # ifdef SOLVE3D # ifdef FORWARD_RHS ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRuct))) THEN got_var(idRuct)=.TRUE. TLM(ng)%pioVar(idRuct)%vd=var_desc(i) TLM(ng)%pioVar(idRuct)%dkind=PIO_FOUT TLM(ng)%pioVar(idRuct)%gtype=u2dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRvct))) THEN got_var(idRvct)=.TRUE. TLM(ng)%pioVar(idRvct)%vd=var_desc(i) TLM(ng)%pioVar(idRvct)%dkind=PIO_FOUT TLM(ng)%pioVar(idRvct)%gtype=v2dvar # endif ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idUfx1))) THEN got_var(idUfx1)=.TRUE. TLM(ng)%pioVar(idUfx1)%vd=var_desc(i) TLM(ng)%pioVar(idUfx1)%dkind=PIO_FOUT TLM(ng)%pioVar(idUfx1)%gtype=u2dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idUfx2))) THEN got_var(idUfx2)=.TRUE. TLM(ng)%pioVar(idUfx2)%vd=var_desc(i) TLM(ng)%pioVar(idUfx2)%dkind=PIO_FOUT TLM(ng)%pioVar(idUfx2)%gtype=u2dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVfx1))) THEN got_var(idVfx1)=.TRUE. TLM(ng)%pioVar(idVfx1)%vd=var_desc(i) TLM(ng)%pioVar(idVfx1)%dkind=PIO_FOUT TLM(ng)%pioVar(idVfx1)%gtype=v2dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVfx2))) THEN got_var(idVfx2)=.TRUE. TLM(ng)%pioVar(idVfx2)%vd=var_desc(i) TLM(ng)%pioVar(idVfx2)%dkind=PIO_FOUT TLM(ng)%pioVar(idVfx2)%gtype=v2dvar # ifdef FORWARD_RHS ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRu3d))) THEN got_var(idRu3d)=.TRUE. TLM(ng)%pioVar(idRu3d)%vd=var_desc(i) TLM(ng)%pioVar(idRu3d)%dkind=PIO_FOUT TLM(ng)%pioVar(idRu3d)%gtype=u3dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idRv3d))) THEN got_var(idRv3d)=.TRUE. TLM(ng)%pioVar(idRv3d)%vd=var_desc(i) TLM(ng)%pioVar(idRv3d)%dkind=PIO_FOUT TLM(ng)%pioVar(idRv3d)%gtype=v3dvar # endif # endif # endif # ifdef SOLVE3D ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idUvel))) THEN got_var(idUvel)=.TRUE. TLM(ng)%pioVar(idUvel)%vd=var_desc(i) TLM(ng)%pioVar(idUvel)%dkind=PIO_FOUT TLM(ng)%pioVar(idUvel)%gtype=u3dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVvel))) THEN got_var(idVvel)=.TRUE. TLM(ng)%pioVar(idVvel)%vd=var_desc(i) TLM(ng)%pioVar(idVvel)%dkind=PIO_FOUT TLM(ng)%pioVar(idVvel)%gtype=v3dvar # ifdef ADJUST_BOUNDARY ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isUvel)))) THEN got_var(idSbry(isUvel))=.TRUE. TLM(ng)%pioVar(idSbry(isUvel))%vd=var_desc(i) TLM(ng)%pioVar(idSbry(isUvel))%dkind=PIO_FOUT TLM(ng)%pioVar(idSbry(isUvel))%gtype=u3dobc ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isVvel)))) THEN got_var(idSbry(isVvel))=.TRUE. TLM(ng)%pioVar(idSbry(isVvel))%vd=var_desc(i) TLM(ng)%pioVar(idSbry(isVvel))%dkind=PIO_FOUT TLM(ng)%pioVar(idSbry(isVvel))%gtype=v3dobc # endif ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idDano))) THEN got_var(idDano)=.TRUE. TLM(ng)%pioVar(idDano)%vd=var_desc(i) TLM(ng)%pioVar(idDano)%dkind=PIO_FOUT TLM(ng)%pioVar(idDano)%gtype=r3dvar # if defined FORWARD_MIXING && \ (defined BVF_MIXING || defined GLS_MIXING || \ defined LMD_MIXING || defined MY25_MIXING) ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVvis))) THEN got_var(idVvis)=.TRUE. TLM(ng)%pioVar(idVvis)%vd=var_desc(i) TLM(ng)%pioVar(idVvis)%dkind=PIO_FOUT TLM(ng)%pioVar(idVvis)%gtype=w3dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idTdif))) THEN got_var(idTdif)=.TRUE. TLM(ng)%pioVar(idTdif)%vd=var_desc(i) TLM(ng)%pioVar(idTdif)%dkind=PIO_FOUT TLM(ng)%pioVar(idTdif)%gtype=w3dvar # ifdef SALINITY ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idSdif))) THEN got_var(idSdif)=.TRUE. TLM(ng)%pioVar(idSdif)%vd=var_desc(i) TLM(ng)%pioVar(idSdif)%dkind=PIO_FOUT TLM(ng)%pioVar(idSdif)%gtype=w3dvar # endif # if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idMtke))) THEN got_var(idMtke)=.TRUE. TLM(ng)%pioVar(idMtke)%vd=var_desc(i) TLM(ng)%pioVar(idMtke)%dkind=PIO_FOUT TLM(ng)%pioVar(idMtke)%gtype=w3dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVmKK))) THEN got_var(idVmKK)=.TRUE. TLM(ng)%pioVar(idVmKK)%vd=var_desc(i) TLM(ng)%pioVar(idVmKK)%dkind=PIO_FOUT TLM(ng)%pioVar(idVmKK)%gtype=w3dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idMtls))) THEN got_var(idMtls)=.TRUE. TLM(ng)%pioVar(idMtls)%vd=var_desc(i) TLM(ng)%pioVar(idMtls)%dkind=PIO_FOUT TLM(ng)%pioVar(idMtls)%gtype=w3dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVmLS))) THEN got_var(idVmLS)=.TRUE. TLM(ng)%pioVar(idVmLS)%vd=var_desc(i) TLM(ng)%pioVar(idVmLS)%dkind=PIO_FOUT TLM(ng)%pioVar(idVmLS)%gtype=w3dvar # ifdef GLS_MIXING_NOT_YET ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVmKP))) THEN got_var(idVmKP)=.TRUE. TLM(ng)%pioVar(idVmKP)%vd=var_desc(i) TLM(ng)%pioVar(idVmKP)%dkind=PIO_FOUT TLM(ng)%pioVar(idVmKP)%gtype=w3dvar # endif # endif # endif # endif # ifdef ADJUST_WSTRESS ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idUsms))) THEN got_var(idUsms)=.TRUE. TLM(ng)%pioVar(idUsms)%vd=var_desc(i) TLM(ng)%pioVar(idUsms)%dkind=PIO_FOUT TLM(ng)%pioVar(idUsms)%gtype=u2dvar ELSE IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idVsms))) THEN got_var(idVsms)=.TRUE. TLM(ng)%pioVar(idVsms)%vd=var_desc(i) TLM(ng)%pioVar(idVsms)%dkind=PIO_FOUT TLM(ng)%pioVar(idVsms)%gtype=v2dvar # endif END IF # ifdef SOLVE3D DO itrc=1,NT(ng) IF (TRIM(var_name(i)).eq.TRIM(Vname(1,idTvar(itrc)))) THEN got_var(idTvar(itrc))=.TRUE. TLM(ng)%pioTrc(itrc)%vd=var_desc(i) TLM(ng)%pioTrc(itrc)%dkind=PIO_FOUT TLM(ng)%pioTrc(itrc)%gtype=r3dvar # ifdef ADJUST_BOUNDARY ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idSbry(isTvar(itrc))))) THEN got_var(idSbry(isTvar(itrc)))=.TRUE. TLM(ng)%pioVar(idSbry(isTvar(itrc)))%vd=var_desc(i) TLM(ng)%pioVar(idSbry(isTvar(itrc)))%dkind=PIO_FOUT TLM(ng)%pioVar(idSbry(isTvar(itrc)))%gtype=r3dobc # endif # ifdef ADJUST_STFLUX ELSE IF (TRIM(var_name(i)).eq. & & TRIM(Vname(1,idTsur(itrc)))) THEN got_var(idTsur(itrc))=.TRUE. TLM(ng)%pioVar(idTsur(itrc))%vd=var_desc(i) TLM(ng)%pioVar(idTsur(itrc))%dkind=PIO_FOUT TLM(ng)%pioVar(idTsur(itrc))%gtype=r2dvar # endif END IF END DO # endif END DO ! ! Check if tangent variables are available in input NetCDF file. ! IF (.not.got_var(idtime)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idtime)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idFsur).and.Hout(idFsur,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idFsur)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idUbar).and.Hout(idUbar,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idUbar)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVbar).and.Hout(idVbar,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVbar)), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef ADJUST_BOUNDARY IF (.not.got_var(idSbry(isFsur))) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSbry(isFsur))), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idSbry(isUbar))) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSbry(isUbar))), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idSbry(isVbar))) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSbry(isVbar))), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # ifdef FORWARD_WRITE # ifdef FORWARD_RHS IF (.not.got_var(idRzet)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRzet)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idRu2d)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRu2d)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idRv2d)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRv2d)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # ifdef SOLVE3D # ifdef FORWARD_RHS IF (.not.got_var(idRuct)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRuct)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idRvct)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRvct)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif IF (.not.got_var(idUfx1)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idUfx1)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idUfx2)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idUfx2)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVfx1)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVfx1)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVfx2)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVfx2)), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef FORWARD_RHS IF (.not.got_var(idRu3d)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRu3d)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idRv3d)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idRv3d)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # endif # endif # ifdef SOLVE3D IF (.not.got_var(idUvel).and.Hout(idUvel,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idUvel)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVvel).and.Hout(idVvel,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVvel)), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef ADJUST_BOUNDARY IF (.not.got_var(idSbry(isUvel))) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSbry(isUvel))), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idSbry(isVvel))) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSbry(isVvel))), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif IF (.not.got_var(idDano).and.Hout(idDano,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idDano)), & & TRIM(ncname) exit_flag=3 RETURN END IF # if defined FORWARD_MIXING && \ (defined BVF_MIXING || defined GLS_MIXING || \ defined LMD_MIXING || defined MY25_MIXING) IF (.not.got_var(idVvis).and.Hout(idVvis,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVvis)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idTdif).and.Hout(idTdif,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idTdif)), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef SALINITY IF (.not.got_var(idSdif).and.Hout(idSdif,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idSdif)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # if defined GLS_MIXING_NOT_YET || defined MY25_MIXING_NOT_YET IF (.not.got_var(idMtke).and.Hout(idMtke,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idMtke)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVmKK).and.Hout(idVmKK,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVmKK)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idMtls).and.Hout(idMtls,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idMtls)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVmLS).and.Hout(idVmLS,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVmLS)), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef GSL_MIXING IF (.not.got_var(idVmKP).and.Hout(idVmKP,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVmKP)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # endif # endif # endif # ifdef ADJUST_WSTRESS IF (.not.got_var(idUsms).and.Hout(idUsms,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idUsms)), & & TRIM(ncname) exit_flag=3 RETURN END IF IF (.not.got_var(idVsms).and.Hout(idVsms,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idVsms)), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # ifdef SOLVE3D DO itrc=1,NT(ng) IF (.not.got_var(idTvar(itrc)).and.Hout(idTvar(itrc),ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idTvar(itrc))), & & TRIM(ncname) exit_flag=3 RETURN END IF # ifdef ADJUST_BOUNDARY IF (.not.got_var(idSbry(isTvar(itrc)))) THEN IF (Master) WRITE (stdout,70) & & TRIM(Vname(1,idSbry(isTvar(itrc)))), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif # ifdef ADJUST_STFLUX IF (.not.got_var(idTsur(itrc)).and.Hout(idTsur(itrc),ng).and. & & Lstflux(itrc,ng)) THEN IF (Master) WRITE (stdout,70) TRIM(Vname(1,idTsur(itrc))), & & TRIM(ncname) exit_flag=3 RETURN END IF # endif END DO # endif ! ! Set unlimited time record dimension to the appropriate value. ! IF (ndefTLM(ng).gt.0) THEN TLM(ng)%Rindex=((ntstart(ng)-1)- & & ndefTLM(ng)*((ntstart(ng)-1)/ndefTLM(ng)))/ & & nTLM(ng) ELSE TLM(ng)%Rindex=(ntstart(ng)-1)/nTLM(ng) END IF TLM(ng)%Rindex=MIN(TLM(ng)%Rindex,rec_size) END IF QUERY ! 10 FORMAT (2x,'TL_DEF_HIS_PIO - creating tangent file,',t56, & & 'Grid ',i2.2,': ',a) 20 FORMAT (2x,'TL_DEF_HIS_PIO - inquiring tangent file,',t56, & & 'Grid ',i2.2,': ',a) 30 FORMAT (/,' TL_DEF_HIS_PIO - unable to create tangent NetCDF', & & ' file: ',a) 40 FORMAT ('tangent linear',1x,a) 50 FORMAT (1pe11.4,1x,'millimeter') 60 FORMAT (/,' TL_DEF_HIS_PIO - unable to open tangent NetCDF', & & ' file: ',a) 70 FORMAT (/,' TL_DEF_HIS_PIO - unable to find variable: ',a,2x, & & ' in tangent NetCDF file: ',a) ! RETURN END SUBROUTINE tl_def_his_pio # endif #endif END MODULE tl_def_his_mod