MODULE dtadyn !!====================================================================== !! *** MODULE dtadyn *** !! Off-line : interpolation of the physical fields !!====================================================================== !! History : OPA ! 1992-01 (M. Imbard) Original code !! 8.0 ! 1998-04 (L.Bopp MA Foujols) slopes for isopyc. !! - ! 1998-05 (L. Bopp) read output of coupled run !! 8.2 ! 2001-01 (M. Levy et M. Benjelloul) add netcdf FORMAT !! NEMO 1.0 ! 2005-03 (O. Aumont and A. El Moussaoui) F90 !! - ! 2005-12 (C. Ethe) Adapted for DEGINT !! 3.0 ! 2007-06 (C. Ethe) use of iom module !! - ! 2007-09 (C. Ethe) add swap_dyn_data !! 3.3 ! 2010-11 (C. Ethe) Full reorganization of the off-line: phasing with the on-line !!---------------------------------------------------------------------- !!---------------------------------------------------------------------- !! dta_dyn_init : initialization, namelist read, and parameters control !! dta_dyn : Interpolation of the fields !!---------------------------------------------------------------------- USE oce ! ocean dynamics and tracers variables USE c1d ! 1D configuration: lk_c1d USE dom_oce ! ocean domain: variables USE zdf_oce ! ocean vertical physics: variables USE sbc_oce ! surface module: variables USE phycst ! physical constants USE trabbl ! active tracer: bottom boundary layer USE ldfslp ! lateral diffusion: iso-neutral slopes USE ldfeiv ! eddy induced velocity coef. USE ldftra_oce ! ocean tracer lateral physics USE zdfmxl ! vertical physics: mixed layer depth USE eosbn2 ! equation of state - Brunt Vaisala frequency USE lbclnk ! ocean lateral boundary conditions (or mpp link) USE zpshde ! z-coord. with partial steps: horizontal derivatives USE in_out_manager ! I/O manager USE iom ! I/O library USE lib_mpp ! distributed memory computing library USE prtctl ! print control IMPLICIT NONE PRIVATE PUBLIC dta_dyn_init ! called by opa.F90 PUBLIC dta_dyn ! called by step.F90 LOGICAL, PUBLIC :: lperdyn = .TRUE. !: boolean for periodic fields or not LOGICAL, PUBLIC :: lfirdyn = .TRUE. !: boolean for the first call or not INTEGER, PUBLIC :: ndtadyn = 73 !: Number of dat in one year INTEGER, PUBLIC :: ndtatot = 73 !: Number of data in the input field INTEGER, PUBLIC :: nsptint = 1 !: type of spatial interpolation CHARACTER(len=45) :: cfile_grid_T = 'dyna_grid_T.nc' ! name of the grid_T file CHARACTER(len=45) :: cfile_grid_U = 'dyna_grid_U.nc' ! name of the grid_U file CHARACTER(len=45) :: cfile_grid_V = 'dyna_grid_V.nc' ! name of the grid_V file CHARACTER(len=45) :: cfile_grid_W = 'dyna_grid_W.nc' ! name of the grid_W file REAL(wp) :: rnspdta ! number of time step per 2 consecutives data REAL(wp) :: rnspdta2 ! rnspdta * 0.5 INTEGER :: ndyn1, ndyn2 ! INTEGER :: nlecoff = 0 ! switch for the first read INTEGER :: numfl_t, numfl_u, numfl_v, numfl_w REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: tdta ! temperature at two consecutive times REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: sdta ! salinity at two consecutive times REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: udta ! zonal velocity at two consecutive times REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: vdta ! meridional velocity at two consecutive times REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: wdta ! vertical velocity at two consecutive times REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: avtdta ! vertical diffusivity coefficient REAL(wp), DIMENSION(jpi,jpj ,2) :: hmlddta ! mixed layer depth at two consecutive times REAL(wp), DIMENSION(jpi,jpj ,2) :: wspddta ! wind speed at two consecutive times REAL(wp), DIMENSION(jpi,jpj ,2) :: frlddta ! sea-ice fraction at two consecutive times REAL(wp), DIMENSION(jpi,jpj ,2) :: empdta ! E-P at two consecutive times REAL(wp), DIMENSION(jpi,jpj ,2) :: qsrdta ! short wave heat flux at two consecutive times REAL(wp), DIMENSION(jpi,jpj ,2) :: bblxdta ! frequency of bbl in the x direction at 2 consecutive times REAL(wp), DIMENSION(jpi,jpj ,2) :: bblydta ! frequency of bbl in the y direction at 2 consecutive times LOGICAL :: l_offbbl #if defined key_ldfslp REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: uslpdta ! zonal isopycnal slopes REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: vslpdta ! meridional isopycnal slopes REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: wslpidta ! zonal diapycnal slopes REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: wslpjdta ! meridional diapycnal slopes #endif #if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv REAL(wp), DIMENSION(jpi,jpj ,2) :: aeiwdta ! G&M coefficient #endif #if defined key_degrad REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: ahtudta, ahtvdta, ahtwdta ! Lateral diffusivity # if defined key_traldf_eiv REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: aeiudta, aeivdta, aeiwdta ! G&M coefficient # endif #endif !! * Substitutions # include "domzgr_substitute.h90" # include "vectopt_loop_substitute.h90" !!---------------------------------------------------------------------- !! NEMO/OFF 3.3 , NEMO Consortium (2010) !! $Id$ !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) !!---------------------------------------------------------------------- CONTAINS SUBROUTINE dta_dyn( kt ) !!---------------------------------------------------------------------- !! *** ROUTINE dta_dyn *** !! !! ** Purpose : Prepares dynamics and physics fields from an NEMO run !! for an off-line simulation of passive tracers !! !! ** Method : calculates the position of DATA to read READ DATA !! (example month changement) computes slopes IF needed !! interpolates DATA IF needed !!---------------------------------------------------------------------- INTEGER, INTENT(in) :: kt ! ocean time-step index !! INTEGER :: iper, iperm1, iswap, izt ! local integers REAL(wp) :: zt, zweigh ! local scalars !!---------------------------------------------------------------------- zt = ( REAL(kt,wp) + rnspdta2 ) / rnspdta izt = INT( zt ) zweigh = zt - REAL( INT(zt), wp ) IF( lperdyn ) THEN ; iperm1 = MOD( izt, ndtadyn ) ELSE ; iperm1 = MOD( izt, ndtatot - 1 ) + 1 ENDIF iper = iperm1 + 1 IF( iperm1 == 0 ) THEN IF( lperdyn ) THEN iperm1 = ndtadyn ELSE IF( lfirdyn ) THEN IF(lwp) WRITE (numout,*) 'dta_dyn: dynamic file is not periodic with or without interpolation & & we take the first value for the previous period iperm1 = 0 ' END IF END IF END IF iswap = 0 ! 1. First call lfirdyn = true ! ---------------------------- IF( lfirdyn ) THEN ndyn1 = iperm1 ! store the information of the period read ndyn2 = iper IF(lwp) THEN WRITE (numout,*) ' dynamics data read for the period ndyn1 =', ndyn1, & & ' and for the period ndyn2 = ', ndyn2 WRITE (numout,*) ' time step is : ', kt WRITE (numout,*) ' we have ndtadyn = ', ndtadyn, ' records in the dynamic file for one year' END IF ! CALL dynrea( kt, MAX( 1, iperm1) ) ! data read for the iperm1 period IF( lk_ldfslp .AND. .NOT. lk_c1d ) THEN ! Computes slopes (here tsn and avt are used as workspace) tsn (:,:,:,jp_tem) = tdta (:,:,:,2) tsn (:,:,:,jp_sal) = sdta (:,:,:,2) avt(:,:,:) = avtdta(:,:,:,2) CALL eos( tsn, rhd, rhop ) ! Time-filtered in situ density CALL bn2( tsn, rn2 ) ! before Brunt-Vaisala frequency IF( ln_zps ) & & CALL zps_hde( kt, jpts, tsn, gtsu, gtsv, & ! Partial steps: before Horizontal DErivative & rhd, gru , grv ) ! of t, s, rd at the bottom ocean level CALL zdf_mxl( kt ) ! mixed layer depth CALL ldf_slp( kt, rhd, rn2 ) uslpdta (:,:,:,2) = uslp (:,:,:) vslpdta (:,:,:,2) = vslp (:,:,:) wslpidta(:,:,:,2) = wslpi(:,:,:) wslpjdta(:,:,:,2) = wslpj(:,:,:) END IF ! CALL swap_dyn_data ! swap from record 2 to 1 ! iswap = 1 ! indicates swap ! CALL dynrea( kt, iper ) ! data read for the iper period ! IF( lk_ldfslp .AND. .NOT. lk_c1d ) THEN ! Computes slopes (here tsn and avt are used as workspace) tsn (:,:,:,jp_tem) = tdta (:,:,:,2) tsn (:,:,:,jp_sal) = sdta (:,:,:,2) avt(:,:,:) = avtdta(:,:,:,2) ! CALL eos( tsn, rhd, rhop ) ! now in situ density CALL bn2( tsn, rn2 ) ! now Brunt-Vaisala frequency IF( ln_zps ) CALL zps_hde( kt, jpts, tsn, gtsu, gtsv, & ! Partial steps: before Horizontal DErivative & rhd, gru , grv ) ! of t, s, rd at the bottom ocean level CALL zdf_mxl( kt ) ! mixed layer depth CALL ldf_slp( kt, rhd, rn2 ) ! slope of iso-neutral surfaces ! uslpdta (:,:,:,2) = uslp (:,:,:) vslpdta (:,:,:,2) = vslp (:,:,:) wslpidta(:,:,:,2) = wslpi(:,:,:) wslpjdta(:,:,:,2) = wslpj(:,:,:) END IF ! lfirdyn = .FALSE. ! trace the first call ENDIF ! ! And now what we have to do at every time step ! check the validity of the period in memory ! IF( iperm1 /= ndyn1 ) THEN ! IF( iperm1 == 0 ) THEN IF(lwp) THEN WRITE (numout,*) ' dynamic file is not periodic with periodic interpolation' WRITE (numout,*) ' we take the last value for the last period ' WRITE (numout,*) ' iperm1 = 12, iper = 13 ' ENDIF iperm1 = 12 iper = 13 ENDIF ! CALL swap_dyn_data ! We have to prepare a new read of data : swap from record 2 to 1 ! iswap = 1 ! indicates swap ! CALL dynrea( kt, iper ) ! data read for the iper period ! IF( lk_ldfslp .AND. .NOT. lk_c1d ) THEN ! Computes slopes. Caution : here tsn and avt are used as workspace tsn(:,:,:,jp_tem) = tdta (:,:,:,2) tsn(:,:,:,jp_sal) = sdta (:,:,:,2) avt(:,:,:) = avtdta(:,:,:,2) ! CALL eos( tsn, rhd, rhop ) ! now in situ density CALL bn2( tsn, rn2 ) ! now Brunt-Vaisala frequency IF( ln_zps ) CALL zps_hde( kt, jpts, tsn, gtsu, gtsv, & ! Partial steps: before Horizontal DErivative & rhd, gru , grv ) ! of t, s, rd at the bottom ocean level CALL zdf_mxl( kt ) ! mixed layer depth CALL ldf_slp( kt, rhd, rn2 ) ! slope of iso-neutral surfaces ! uslpdta (:,:,:,2) = uslp (:,:,:) vslpdta (:,:,:,2) = vslp (:,:,:) wslpidta(:,:,:,2) = wslpi(:,:,:) wslpjdta(:,:,:,2) = wslpj(:,:,:) END IF ! ndyn1 = ndyn2 ! store the information of the period read ndyn2 = iper ! IF(lwp) THEN WRITE (numout,*) ' dynamics data read for the period ndyn1 =', ndyn1, & & ' and for the period ndyn2 = ', ndyn2 WRITE (numout,*) ' time step is : ', kt END IF ! END IF ! ! Compute the data at the given time step !---------------------------------------- IF( nsptint == 0 ) THEN ! No space interpolation, data are probably correct ! ! We have to initialize data if we have changed the period CALL assign_dyn_data ELSEIF( nsptint == 1 ) THEN ! linear interpolation CALL linear_interp_dyn_data( zweigh ) ELSE ! other interpolation WRITE (numout,*) ' this kind of interpolation do not exist at the moment : we stop' STOP 'dtadyn' END IF ! CALL eos( tsn, rhd, rhop ) ! In any case, we need rhop ! #if ! defined key_degrad && defined key_traldf_c2d ! ! In case of 2D varying coefficients, we need aeiv and aeiu IF( lk_traldf_eiv ) CALL dta_eiv( kt ) ! eddy induced velocity coefficient #endif ! IF( .NOT. l_offbbl ) THEN ! Compute bbl coefficients if needed tsb(:,:,:,:) = tsn(:,:,:,:) CALL bbl( kt, 'TRC') END IF IF(ln_ctl) THEN CALL prt_ctl(tab3d_1=tsn(:,:,:,jp_tem), clinfo1=' tn - : ', mask1=tmask, ovlap=1, kdim=jpk ) CALL prt_ctl(tab3d_1=tsn(:,:,:,jp_sal), clinfo1=' sn - : ', mask1=tmask, ovlap=1, kdim=jpk ) CALL prt_ctl(tab3d_1=un , clinfo1=' un - : ', mask1=tmask, ovlap=1, kdim=jpk ) CALL prt_ctl(tab3d_1=vn , clinfo1=' vn - : ', mask1=tmask, ovlap=1, kdim=jpk ) CALL prt_ctl(tab3d_1=wn , clinfo1=' wn - : ', mask1=tmask, ovlap=1, kdim=jpk ) CALL prt_ctl(tab3d_1=avt , clinfo1=' kz - : ', mask1=tmask, ovlap=1, kdim=jpk ) CALL prt_ctl(tab2d_1=fr_i , clinfo1=' fr_i - : ', mask1=tmask, ovlap=1 ) CALL prt_ctl(tab2d_1=hmld , clinfo1=' hmld - : ', mask1=tmask, ovlap=1 ) CALL prt_ctl(tab2d_1=emps , clinfo1=' emps - : ', mask1=tmask, ovlap=1 ) CALL prt_ctl(tab2d_1=wndm , clinfo1=' wspd - : ', mask1=tmask, ovlap=1 ) CALL prt_ctl(tab2d_1=qsr , clinfo1=' qsr - : ', mask1=tmask, ovlap=1 ) ENDIF ! END SUBROUTINE dta_dyn SUBROUTINE dynrea( kt, kenr ) !!---------------------------------------------------------------------- !! *** ROUTINE dynrea *** !! !! ** Purpose : READ dynamics fiels from OPA9 netcdf output !! !! ** Method : READ the kenr records of DATA and store in udta(...,2), .... !!---------------------------------------------------------------------- INTEGER, INTENT(in) :: kt, kenr ! time index !! INTEGER :: jkenr REAL(wp), DIMENSION(jpi,jpj,jpk) :: zu, zv, zw, zt, zs, zavt , zhdiv ! 3D workspace REAL(wp), DIMENSION(jpi,jpj) :: zemp, zqsr, zmld, zice, zwspd, ztaux, ztauy ! 2D workspace REAL(wp), DIMENSION(jpi,jpj) :: zbblx, zbbly #if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv REAL(wp), DIMENSION(jpi,jpj) :: zaeiw #endif #if defined key_degrad REAL(wp), DIMENSION(jpi,jpj,jpk) :: zahtu, zahtv, zahtw ! Lateral diffusivity # if defined key_traldf_eiv REAL(wp), DIMENSION(jpi,jpj,jpk) :: zaeiu, zaeiv, zaeiw ! G&M coefficient # endif #endif !!---------------------------------------------------------------------- ! 0. Initialization ! cas d'un fichier non periodique : on utilise deux fois le premier et ! le dernier champ temporel jkenr = kenr IF(lwp) THEN WRITE(numout,*) WRITE(numout,*) 'Dynrea : read dynamical fields, kenr = ', jkenr WRITE(numout,*) '~~~~~~~' #if defined key_degrad WRITE(numout,*) ' Degraded fields' #endif WRITE(numout,*) ENDIF IF( kt == nit000 .AND. nlecoff == 0 ) THEN nlecoff = 1 CALL iom_open ( cfile_grid_T, numfl_t ) CALL iom_open ( cfile_grid_U, numfl_u ) CALL iom_open ( cfile_grid_V, numfl_v ) CALL iom_open ( cfile_grid_W, numfl_w ) ENDIF ! file grid-T !--------------- CALL iom_get( numfl_t, jpdom_data, 'votemper', zt (:,:,:), jkenr ) CALL iom_get( numfl_t, jpdom_data, 'vosaline', zs (:,:,:), jkenr ) CALL iom_get( numfl_t, jpdom_data, 'somixhgt', zmld (:,: ), jkenr ) CALL iom_get( numfl_t, jpdom_data, 'sowaflcd', zemp (:,: ), jkenr ) CALL iom_get( numfl_t, jpdom_data, 'soshfldo', zqsr (:,: ), jkenr ) CALL iom_get( numfl_t, jpdom_data, 'soicecov', zice (:,: ), jkenr ) IF( iom_varid( numfl_t, 'sowindsp', ldstop = .FALSE. ) > 0 ) THEN CALL iom_get( numfl_t, jpdom_data, 'sowindsp', zwspd(:,:), jkenr ) ELSE CALL iom_get( numfl_u, jpdom_data, 'sozotaux', ztaux(:,:), jkenr ) CALL iom_get( numfl_v, jpdom_data, 'sometauy', ztauy(:,:), jkenr ) CALL tau2wnd( ztaux, ztauy, zwspd ) ENDIF ! files grid-U / grid_V CALL iom_get( numfl_u, jpdom_data, 'vozocrtx', zu (:,:,:), jkenr ) CALL iom_get( numfl_v, jpdom_data, 'vomecrty', zv (:,:,:), jkenr ) #if defined key_trabbl IF( .NOT. lk_c1d .AND. nn_bbl_ldf == 1 ) THEN IF( iom_varid( numfl_u, 'sobblcox', ldstop = .FALSE. ) > 0 .AND. & & iom_varid( numfl_v, 'sobblcoy', ldstop = .FALSE. ) > 0 ) THEN CALL iom_get( numfl_u, jpdom_data, 'sobblcox', zbblx(:,:), jkenr ) CALL iom_get( numfl_v, jpdom_data, 'sobblcoy', zbbly(:,:), jkenr ) l_offbbl = .TRUE. ENDIF ENDIF #endif ! file grid-W !! CALL iom_get ( numfl_w, jpdom_data, 'vovecrtz', zw (:,:,:), jkenr ) ! Computation of vertical velocity using horizontal divergence CALL wzv( zu, zv, zw, zhdiv ) IF( iom_varid( numfl_w, 'voddmavs', ldstop = .FALSE. ) > 0 ) THEN ! avs exist: it is used CALL iom_get( numfl_w, jpdom_data, 'voddmavs', zavt (:,:,:), jkenr ) ELSE ! no avs: use avt CALL iom_get( numfl_w, jpdom_data, 'votkeavt', zavt (:,:,:), jkenr ) ENDIF #if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv CALL iom_get( numfl_w, jpdom_data, 'soleaeiw', zaeiw (:,: ), jkenr ) #endif #if defined key_degrad CALL iom_get( numfl_u, jpdom_data, 'vozoahtu', zahtu(:,:,:), jkenr ) CALL iom_get( numfl_v, jpdom_data, 'vomeahtv', zahtv(:,:,:), jkenr ) CALL iom_get( numfl_w, jpdom_data, 'voveahtw', zahtw(:,:,:), jkenr ) # if defined key_traldf_eiv CALL iom_get( numfl_u, jpdom_data, 'vozoaeiu', zaeiu(:,:,:), jkenr ) CALL iom_get( numfl_v, jpdom_data, 'vomeaeiv', zaeiv(:,:,:), jkenr ) CALL iom_get( numfl_w, jpdom_data, 'voveaeiw', zaeiw(:,:,:), jkenr ) # endif #endif udta(:,:,:,2) = zu(:,:,:) * umask(:,:,:) vdta(:,:,:,2) = zv(:,:,:) * vmask(:,:,:) wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:) tdta(:,:,:,2) = zt (:,:,:) * tmask(:,:,:) sdta(:,:,:,2) = zs (:,:,:) * tmask(:,:,:) avtdta(:,:,:,2) = zavt(:,:,:) * tmask(:,:,:) #if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv aeiwdta(:,:,2) = zaeiw(:,:) * tmask(:,:,1) #endif #if defined key_degrad ahtudta(:,:,:,2) = zahtu(:,:,:) * umask(:,:,:) ahtvdta(:,:,:,2) = zahtv(:,:,:) * vmask(:,:,:) ahtwdta(:,:,:,2) = zahtw(:,:,:) * tmask(:,:,:) # if defined key_traldf_eiv aeiudta(:,:,:,2) = zaeiu(:,:,:) * umask(:,:,:) aeivdta(:,:,:,2) = zaeiv(:,:,:) * vmask(:,:,:) aeiwdta(:,:,:,2) = zaeiw(:,:,:) * tmask(:,:,:) # endif #endif ! fluxes ! wspddta(:,:,2) = zwspd(:,:) * tmask(:,:,1) frlddta(:,:,2) = MIN( 1., zice(:,:) ) * tmask(:,:,1) empdta (:,:,2) = zemp(:,:) * tmask(:,:,1) qsrdta (:,:,2) = zqsr(:,:) * tmask(:,:,1) hmlddta(:,:,2) = zmld(:,:) * tmask(:,:,1) #if defined key_trabbl IF( l_offbbl ) THEN bblxdta(:,:,2) = MAX( 0., zbblx(:,:) ) bblydta(:,:,2) = MAX( 0., zbbly(:,:) ) WHERE( bblxdta(:,:,2) > 2. ) bblxdta(:,:,2) = 0. WHERE( bblydta(:,:,2) > 2. ) bblydta(:,:,2) = 0. ENDIF #endif IF( kt == nitend ) THEN CALL iom_close ( numfl_t ) CALL iom_close ( numfl_u ) CALL iom_close ( numfl_v ) CALL iom_close ( numfl_w ) ENDIF ! END SUBROUTINE dynrea SUBROUTINE dta_dyn_init !!---------------------------------------------------------------------- !! *** ROUTINE dta_dyn_init *** !! !! ** Purpose : initializations of parameters for the interpolation !! !! ** Method : !!---------------------------------------------------------------------- REAL(wp) :: znspyr !: number of time step per year !! NAMELIST/namdyn/ ndtadyn, ndtatot, nsptint, lperdyn, & & cfile_grid_T, cfile_grid_U, cfile_grid_V, cfile_grid_W !!---------------------------------------------------------------------- ! Define the dynamical input parameters ! ====================================== REWIND( numnam ) ! Read Namelist namdyn : Lateral physics on tracers READ ( numnam, namdyn ) IF(lwp) THEN ! control print WRITE(numout,*) WRITE(numout,*) 'namdyn : offline dynamical selection' WRITE(numout,*) '~~~~~~~' WRITE(numout,*) ' Namelist namdyn : set parameters for the lecture of the dynamical fields' WRITE(numout,*) WRITE(numout,*) ' number of elements in the FILE for a year ndtadyn = ' , ndtadyn WRITE(numout,*) ' total number of elements in the FILE ndtatot = ' , ndtatot WRITE(numout,*) ' type of interpolation nsptint = ' , nsptint WRITE(numout,*) ' loop on the same FILE lperdyn = ' , lperdyn WRITE(numout,*) ' ' WRITE(numout,*) ' name of grid_T file cfile_grid_T = ', TRIM(cfile_grid_T) WRITE(numout,*) ' name of grid_U file cfile_grid_U = ', TRIM(cfile_grid_U) WRITE(numout,*) ' name of grid_V file cfile_grid_V = ', TRIM(cfile_grid_V) WRITE(numout,*) ' name of grid_W file cfile_grid_W = ', TRIM(cfile_grid_W) WRITE(numout,*) ' ' ENDIF ! znspyr = nyear_len(1) * rday / rdt rnspdta = znspyr / FLOAT( ndtadyn ) rnspdta2 = rnspdta * 0.5 ! CALL dta_dyn( nit000 ) ! END SUBROUTINE dta_dyn_init SUBROUTINE wzv( pu, pv, pw, phdiv ) !!---------------------------------------------------------------------- !! *** ROUTINE wzv *** !! !! ** Purpose : Compute the now vertical velocity after the array swap !! !! ** Method : - compute the now divergence given by : !! * z-coordinate ONLY !!!! !! hdiv = 1/(e1t*e2t) [ di(e2u u) + dj(e1v v) ] !! - Using the incompressibility hypothesis, the vertical !! velocity is computed by integrating the horizontal divergence !! from the bottom to the surface. !! The boundary conditions are w=0 at the bottom (no flux). !!---------------------------------------------------------------------- REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in ) :: pu, pv !: horizontal velocities REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( out) :: pw !: verticla velocity REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) :: phdiv !: horizontal divergence !! INTEGER :: ji, jj, jk REAL(wp) :: zu, zu1, zv, zv1, zet !!---------------------------------------------------------------------- ! ! Computation of vertical velocity using horizontal divergence phdiv(:,:,:) = 0. DO jk = 1, jpkm1 DO jj = 2, jpjm1 DO ji = fs_2, fs_jpim1 ! vector opt. zu = pu(ji ,jj ,jk) * umask(ji ,jj ,jk) * e2u(ji ,jj ) * fse3u(ji ,jj ,jk) zu1 = pu(ji-1,jj ,jk) * umask(ji-1,jj ,jk) * e2u(ji-1,jj ) * fse3u(ji-1,jj ,jk) zv = pv(ji ,jj ,jk) * vmask(ji ,jj ,jk) * e1v(ji ,jj ) * fse3v(ji ,jj ,jk) zv1 = pv(ji ,jj-1,jk) * vmask(ji ,jj-1,jk) * e1v(ji ,jj-1) * fse3v(ji ,jj-1,jk) zet = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) phdiv(ji,jj,jk) = ( zu - zu1 + zv - zv1 ) * zet END DO END DO END DO CALL lbc_lnk( phdiv, 'T', 1. ) ! Lateral boundary conditions on phdiv ! ! computation of vertical velocity from the bottom pw(:,:,jpk) = 0._wp DO jk = jpkm1, 1, -1 pw(:,:,jk) = pw(:,:,jk+1) - fse3t(:,:,jk) * phdiv(:,:,jk) END DO ! END SUBROUTINE wzv SUBROUTINE dta_eiv( kt ) !!---------------------------------------------------------------------- !! *** ROUTINE dta_eiv *** !! !! ** Purpose : Compute the eddy induced velocity coefficient from the !! growth rate of baroclinic instability. !! !! ** Method : Specific to the offline model. Computes the horizontal !! values from the vertical value !!---------------------------------------------------------------------- INTEGER, INTENT( in ) :: kt ! ocean time-step inedx !! INTEGER :: ji, jj ! dummy loop indices !!---------------------------------------------------------------------- ! IF( kt == nit000 ) THEN IF(lwp) WRITE(numout,*) IF(lwp) WRITE(numout,*) 'dta_eiv : eddy induced velocity coefficients' IF(lwp) WRITE(numout,*) '~~~~~~~' ENDIF ! #if defined key_ldfeiv ! Average the diffusive coefficient at u- v- points DO jj = 2, jpjm1 DO ji = fs_2, fs_jpim1 ! vector opt. aeiu(ji,jj) = .5 * ( aeiw(ji,jj) + aeiw(ji+1,jj ) ) aeiv(ji,jj) = .5 * ( aeiw(ji,jj) + aeiw(ji ,jj+1) ) END DO END DO CALL lbc_lnk( aeiu, 'U', 1. ) ; CALL lbc_lnk( aeiv, 'V', 1. ) ! lateral boundary condition #endif ! END SUBROUTINE dta_eiv SUBROUTINE tau2wnd( ptaux, ptauy, pwspd ) !!--------------------------------------------------------------------- !! *** ROUTINE sbc_tau2wnd *** !! !! ** Purpose : Estimation of wind speed as a function of wind stress !! !! ** Method : |tau|=rhoa*Cd*|U|^2 !!--------------------------------------------------------------------- REAL(wp), DIMENSION(jpi,jpj), INTENT(in ) :: ptaux, ptauy ! wind stress in i-j direction resp. REAL(wp), DIMENSION(jpi,jpj), INTENT( out) :: pwspd ! wind speed !! REAL(wp) :: zrhoa = 1.22_wp ! Air density kg/m3 REAL(wp) :: zcdrag = 1.5e-3_wp ! drag coefficient REAL(wp) :: ztx, zty, ztau, zcoef ! temporary variables INTEGER :: ji, jj ! dummy indices !!--------------------------------------------------------------------- zcoef = 1. / ( zrhoa * zcdrag ) !CDIR NOVERRCHK DO jj = 2, jpjm1 !CDIR NOVERRCHK DO ji = fs_2, fs_jpim1 ! vector opt. ztx = ptaux(ji,jj) * umask(ji,jj,1) + ptaux(ji-1,jj ) * umask(ji-1,jj ,1) zty = ptauy(ji,jj) * vmask(ji,jj,1) + ptauy(ji ,jj-1) * vmask(ji ,jj-1,1) ztau = 0.5 * SQRT( ztx * ztx + zty * zty ) pwspd(ji,jj) = SQRT ( ztau * zcoef ) * tmask(ji,jj,1) END DO END DO CALL lbc_lnk( pwspd(:,:), 'T', 1. ) ! END SUBROUTINE tau2wnd SUBROUTINE swap_dyn_data !!---------------------------------------------------------------------- !! *** ROUTINE swap_dyn_data *** !! !! ** Purpose : swap array data !!---------------------------------------------------------------------- ! ! swap from record 2 to 1 tdta (:,:,:,1) = tdta (:,:,:,2) sdta (:,:,:,1) = sdta (:,:,:,2) avtdta (:,:,:,1) = avtdta (:,:,:,2) udta (:,:,:,1) = udta (:,:,:,2) vdta (:,:,:,1) = vdta (:,:,:,2) wdta (:,:,:,1) = wdta (:,:,:,2) #if defined key_ldfslp && ! defined key_c1d uslpdta (:,:,:,1) = uslpdta (:,:,:,2) vslpdta (:,:,:,1) = vslpdta (:,:,:,2) wslpidta(:,:,:,1) = wslpidta(:,:,:,2) wslpjdta(:,:,:,1) = wslpjdta(:,:,:,2) #endif hmlddta(:,:,1) = hmlddta(:,:,2) wspddta(:,:,1) = wspddta(:,:,2) frlddta(:,:,1) = frlddta(:,:,2) empdta (:,:,1) = empdta (:,:,2) qsrdta (:,:,1) = qsrdta (:,:,2) IF( l_offbbl ) THEN bblxdta(:,:,1) = bblxdta(:,:,2) bblydta(:,:,1) = bblydta(:,:,2) ENDIF #if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv aeiwdta(:,:,1) = aeiwdta(:,:,2) #endif #if defined key_degrad ahtudta(:,:,:,1) = ahtudta(:,:,:,2) ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) # if defined key_traldf_eiv aeiudta(:,:,:,1) = aeiudta(:,:,:,2) aeivdta(:,:,:,1) = aeivdta(:,:,:,2) aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) # endif #endif ! END SUBROUTINE swap_dyn_data SUBROUTINE assign_dyn_data !!---------------------------------------------------------------------- !! *** ROUTINE assign_dyn_data *** !! !! ** Purpose : Assign dynamical data to the data that have been read !! without time interpolation !! !!---------------------------------------------------------------------- tsn(:,:,:,jp_tem) = tdta (:,:,:,2) tsn(:,:,:,jp_sal) = sdta (:,:,:,2) avt(:,:,:) = avtdta(:,:,:,2) un (:,:,:) = udta (:,:,:,2) vn (:,:,:) = vdta (:,:,:,2) wn (:,:,:) = wdta (:,:,:,2) #if defined key_ldfslp && ! defined key_c1d uslp (:,:,:) = uslpdta (:,:,:,2) vslp (:,:,:) = vslpdta (:,:,:,2) wslpi(:,:,:) = wslpidta(:,:,:,2) wslpj(:,:,:) = wslpjdta(:,:,:,2) #endif hmld(:,:) = hmlddta(:,:,2) wndm(:,:) = wspddta(:,:,2) fr_i(:,:) = frlddta(:,:,2) emp (:,:) = empdta (:,:,2) emps(:,:) = emp(:,:) qsr (:,:) = qsrdta (:,:,2) #if defined key_trabbl IF( l_offbbl ) THEN ahu_bbl(:,:) = bblxdta(:,:,2) ahv_bbl(:,:) = bblydta(:,:,2) ENDIF #endif #if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv aeiw(:,:) = aeiwdta(:,:,2) #endif #if defined key_degrad ahtu(:,:,:) = ahtudta(:,:,:,2) ahtv(:,:,:) = ahtvdta(:,:,:,2) ahtw(:,:,:) = ahtwdta(:,:,:,2) # if defined key_traldf_eiv aeiu(:,:,:) = aeiudta(:,:,:,2) aeiv(:,:,:) = aeivdta(:,:,:,2) aeiw(:,:,:) = aeiwdta(:,:,:,2) # endif #endif ! END SUBROUTINE assign_dyn_data SUBROUTINE linear_interp_dyn_data( pweigh ) !!---------------------------------------------------------------------- !! *** ROUTINE linear_interp_dyn_data *** !! !! ** Purpose : linear interpolation of data !!---------------------------------------------------------------------- REAL(wp), INTENT(in) :: pweigh ! weigh !! REAL(wp) :: zweighm1 !!---------------------------------------------------------------------- zweighm1 = 1. - pweigh tsn(:,:,:,jp_tem) = zweighm1 * tdta (:,:,:,1) + pweigh * tdta (:,:,:,2) tsn(:,:,:,jp_sal) = zweighm1 * sdta (:,:,:,1) + pweigh * sdta (:,:,:,2) avt(:,:,:) = zweighm1 * avtdta(:,:,:,1) + pweigh * avtdta(:,:,:,2) un (:,:,:) = zweighm1 * udta (:,:,:,1) + pweigh * udta (:,:,:,2) vn (:,:,:) = zweighm1 * vdta (:,:,:,1) + pweigh * vdta (:,:,:,2) wn (:,:,:) = zweighm1 * wdta (:,:,:,1) + pweigh * wdta (:,:,:,2) #if defined key_ldfslp && ! defined key_c1d uslp (:,:,:) = zweighm1 * uslpdta (:,:,:,1) + pweigh * uslpdta (:,:,:,2) vslp (:,:,:) = zweighm1 * vslpdta (:,:,:,1) + pweigh * vslpdta (:,:,:,2) wslpi(:,:,:) = zweighm1 * wslpidta(:,:,:,1) + pweigh * wslpidta(:,:,:,2) wslpj(:,:,:) = zweighm1 * wslpjdta(:,:,:,1) + pweigh * wslpjdta(:,:,:,2) #endif hmld(:,:) = zweighm1 * hmlddta(:,:,1) + pweigh * hmlddta(:,:,2) wndm(:,:) = zweighm1 * wspddta(:,:,1) + pweigh * wspddta(:,:,2) fr_i(:,:) = zweighm1 * frlddta(:,:,1) + pweigh * frlddta(:,:,2) emp (:,:) = zweighm1 * empdta (:,:,1) + pweigh * empdta (:,:,2) emps(:,:) = emp(:,:) qsr (:,:) = zweighm1 * qsrdta (:,:,1) + pweigh * qsrdta (:,:,2) #if defined key_trabbl IF( l_offbbl ) THEN ahu_bbl(:,:) = zweighm1 * bblxdta(:,:,1) + pweigh * bblxdta(:,:,2) ahv_bbl(:,:) = zweighm1 * bblydta(:,:,1) + pweigh * bblydta(:,:,2) ENDIF #endif #if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv aeiw(:,:) = zweighm1 * aeiwdta(:,:,1) + pweigh * aeiwdta(:,:,2) #endif #if defined key_degrad ahtu(:,:,:) = zweighm1 * ahtudta(:,:,:,1) + pweigh * ahtudta(:,:,:,2) ahtv(:,:,:) = zweighm1 * ahtvdta(:,:,:,1) + pweigh * ahtvdta(:,:,:,2) ahtw(:,:,:) = zweighm1 * ahtwdta(:,:,:,1) + pweigh * ahtwdta(:,:,:,2) # if defined key_traldf_eiv aeiu(:,:,:) = zweighm1 * aeiudta(:,:,:,1) + pweigh * aeiudta(:,:,:,2) aeiv(:,:,:) = zweighm1 * aeivdta(:,:,:,1) + pweigh * aeivdta(:,:,:,2) aeiw(:,:,:) = zweighm1 * aeiwdta(:,:,:,1) + pweigh * aeiwdta(:,:,:,2) # endif #endif ! END SUBROUTINE linear_interp_dyn_data !!====================================================================== END MODULE dtadyn