MODULE nonoro_gwd_ran_mod IMPLICIT NONE REAL,allocatable,save :: du_nonoro_gwd(:,:) ! Zonal wind tendency due to GWD REAL,allocatable,save :: dv_nonoro_gwd(:,:) ! Meridional wind tendency due to GWD REAL,ALLOCATABLE,SAVE :: east_gwstress(:,:) ! Profile of eastward stress REAL,ALLOCATABLE,SAVE :: west_gwstress(:,:) ! Profile of westward stress !$OMP THREADPRIVATE(du_nonoro_gwd,dv_nonoro_gwd,east_gwstress,west_gwstress) CONTAINS SUBROUTINE NONORO_GWD_RAN(ngrid,nlayer,DTIME, cpnew, rnew, pp, & zmax_therm, pt, pu, pv, pdt, pdu, pdv, & zustr,zvstr,d_t, d_u, d_v) !-------------------------------------------------------------------------------- ! Parametrization of the momentum flux deposition due to a discrete ! number of gravity waves. ! F. Lott ! Version 14, Gaussian distribution of the source ! LMDz model online version ! ADAPTED FOR VENUS / F. LOTT + S. LEBONNOIS ! Version adapted on 03/04/2013: ! - input flux compensated in the deepest layers ! ! ADAPTED FOR MARS G.GILLI 02/2016 ! Revision with F.Forget 06/2016 Variable EP-flux according to ! PBL variation (max velocity thermals) ! UPDATED D.BARDET 01/2020 - reproductibility of the ! launching altitude calculation ! - wave characteristic ! calculation using MOD ! - adding east_gwstress and ! west_gwstress variables ! UPDATED J.LIU 12/2021 The rho (density) at the specific ! locations is introduced. The equation ! of EP-flux is corrected by adding the ! term of density at launch (source) ! altitude. ! !--------------------------------------------------------------------------------- use comcstfi_h, only: g, pi, r USE ioipsl_getin_p_mod, ONLY : getin_p use vertical_layers_mod, only : presnivs use geometry_mod, only: cell_area #ifdef CPP_XIOS use xios_output_mod, only: send_xios_field #endif implicit none include "yoegwd.h" include "callkeys.h" CHARACTER (LEN=20) :: modname='nonoro_gwd_ran' CHARACTER (LEN=80) :: abort_message ! 0. DECLARATIONS: ! 0.1 INPUTS INTEGER, intent(in):: ngrid ! number of atmospheric columns INTEGER, intent(in):: nlayer ! number of atmospheric layers REAL, intent(in):: DTIME ! Time step of the Physics(s) REAL, intent(in):: zmax_therm(ngrid) ! Altitude of max velocity thermals (m) REAL,INTENT(IN) :: cpnew(ngrid,nlayer)! Cp of the atmosphere REAL,INTENT(IN) :: rnew(ngrid,nlayer) ! R of the atmosphere REAL, intent(in):: pp(ngrid,nlayer) ! Pressure at full levels(Pa) REAL, intent(in):: pt(ngrid,nlayer) ! Temp at full levels(K) REAL, intent(in):: pu(ngrid,nlayer) ! Zonal wind at full levels(m/s) REAL, intent(in):: pv(ngrid,nlayer) ! Meridional winds at full levels(m/s) REAL,INTENT(in) :: pdt(ngrid,nlayer) ! Tendency on temperature (K/s) REAL,INTENT(in) :: pdu(ngrid,nlayer) ! Tendency on zonal wind (m/s/s) REAL,INTENT(in) :: pdv(ngrid,nlayer) ! Tendency on meridional wind (m/s/s) ! 0.2 OUTPUTS REAL, intent(out):: zustr(ngrid) ! Zonal surface stress REAL, intent(out):: zvstr(ngrid) ! Meridional surface stress REAL, intent(out):: d_t(ngrid, nlayer) ! Tendency on temperature (K/s) due to gravity waves (not used set to zero) REAL, intent(out):: d_u(ngrid, nlayer) ! Tendency on zonal wind (m/s/s) due to gravity waves REAL, intent(out):: d_v(ngrid, nlayer) ! Tendency on meridional wind (m/s/s) due to gravity waves ! 0.3 INTERNAL ARRAYS REAL :: TT(ngrid, nlayer) ! Temperature at full levels REAL :: RHO(ngrid, nlayer) ! Mass density at full levels REAL :: UU(ngrid, nlayer) ! Zonal wind at full levels REAL :: VV(ngrid, nlayer) ! Meridional winds at full levels REAL :: BVLOW(ngrid) ! initial N at each grid (not used) INTEGER II, JJ, LL ! 0.3.0 TIME SCALE OF THE LIFE CYCLE OF THE WAVES PARAMETERIZED REAL, parameter:: DELTAT = 24. * 3600. ! 0.3.1 GRAVITY-WAVES SPECIFICATIONS INTEGER, PARAMETER:: NK = 2 ! number of horizontal wavenumbers INTEGER, PARAMETER:: NP = 2 ! directions (eastward and westward) phase speed INTEGER, PARAMETER:: NO = 2 ! absolute values of phase speed INTEGER, PARAMETER:: NA = 5 ! number of realizations to get the phase speed INTEGER, PARAMETER:: NW = NK * NP * NO ! Total numbers of gravity waves INTEGER JK, JP, JO, JW ! Loop index for NK,NP,NO, and NW REAL, save :: kmax ! Max horizontal wavenumber (lambda min,lambda=2*pi/kmax),kmax=N/u, u(=30~50) zonal wind velocity at launch altitude !$OMP THREADPRIVATE(kmax) REAL, parameter:: kmin = 2.e-5 ! Min horizontal wavenumber (lambda max = 314 km,lambda=2*pi/kmin) REAL kstar ! Min horizontal wavenumber constrained by horizontal resolution REAL, save :: cmax ! Max horizontal absolute phase velocity (at the maginitide of zonal wind u at the launch altitude) !$OMP THREADPRIVATE(cmax) REAL, parameter:: cmin = 1. ! Min horizontal absolute phase velocity (not used) REAL CPHA ! absolute PHASE VELOCITY frequency REAL ZK(NW, ngrid) ! Horizontal wavenumber amplitude REAL ZP(NW, ngrid) ! Horizontal wavenumber angle REAL ZO(NW, ngrid) ! Absolute frequency REAL intr_freq_m(nw, ngrid) ! Waves Intr. freq. at the 1/2 lev below the full level (previous name: ZOM) REAL intr_freq_p(nw, ngrid) ! Waves Intr. freq. at the 1/2 lev above the full level (previous name: ZOP) REAL wwm(nw, ngrid) ! Wave EP-fluxes at the 1/2 level below the full level REAL wwp(nw, ngrid) ! Wave EP-fluxes at the 1/2 level above the full level REAL u_epflux_p(nw, ngrid) ! Partial zonal flux (=for each wave) at the 1/2 level above the full level (previous name: RUWP) REAL v_epflux_p(nw, ngrid) ! Partial meridional flux (=for each wave) at the 1/2 level above the full level (previous name: RVWP) REAL u_epflux_tot(ngrid, nlayer + 1) ! Total zonal flux (=for all waves (nw)) at the 1/2 level above the full level (3D) (previous name: RUW) REAL v_epflux_tot(ngrid, nlayer + 1) ! Total meridional flux (=for all waves (nw)) at the 1/2 level above the full level (3D) (previous name: RVW) REAL epflux_0(nw, ngrid) ! Fluxes at launching level (previous name: RUW0) REAL, save :: epflux_max ! Max EP flux value at launching altitude (previous name: RUWMAX, tunable) !$OMP THREADPRIVATE(epflux_max) INTEGER LAUNCH ! Launching altitude REAL, parameter:: xlaunch = 0.4 ! Control the launching altitude by pressure REAL, parameter:: zmaxth_top = 8000. ! Top of convective layer (approx. not used) ! 0.3.2 PARAMETERS OF WAVES DISSIPATIONS REAL, save :: sat ! saturation parameter(tunable) !$OMP THREADPRIVATE(sat) REAL, save :: rdiss ! dissipation coefficient (tunable) !$OMP THREADPRIVATE(rdiss) REAL, parameter:: zoisec = 1.e-10 ! security for intrinsic frequency ! 0.3.3 Background flow at 1/2 levels and vertical coordinate REAL H0bis(ngrid, nlayer) ! Characteristic Height of the atmosphere (specific locations) REAL, save:: H0 ! Characteristic Height of the atmosphere (constant) !$OMP THREADPRIVATE(H0) REAL, parameter:: pr = 250 ! Reference pressure [Pa] REAL, parameter:: tr = 220. ! Reference temperature [K] REAL ZH(ngrid, nlayer + 1) ! Log-pressure altitude (constant H0) REAL ZHbis(ngrid, nlayer + 1) ! Log-pressure altitude (varying H0bis) REAL UH(ngrid, nlayer + 1) ! zonal wind at 1/2 levels REAL VH(ngrid, nlayer + 1) ! meridional wind at 1/2 levels REAL PH(ngrid, nlayer + 1) ! Pressure at 1/2 levels REAL, parameter:: psec = 1.e-20 ! Security to avoid division by 0 pressure(!!IMPORTANT: should be lower than the topmost layer's pressure) REAL BV(ngrid, nlayer + 1) ! Brunt Vaisala freq. (N^2) at 1/2 levels REAL, parameter:: bvsec = 1.e-5 ! Security to avoid negative BV REAL HREF(nlayer + 1) ! Reference pressure for launching altitude REAL RAN_NUM_1,RAN_NUM_2,RAN_NUM_3 LOGICAL,SAVE :: firstcall = .true. !$OMP THREADPRIVATE(firstcall) !----------------------------------------------------------------------------------------------------------------------- ! 1. INITIALISATIONS !----------------------------------------------------------------------------------------------------------------------- IF (firstcall) THEN write(*,*) "nonoro_gwd_ran: FLott non-oro GW scheme is active!" epflux_max = 7.E-7 ! Mars' value !! call getin_p("nonoro_gwd_epflux_max", epflux_max) write(*,*) "nonoro_gwd_ran: epflux_max=", epflux_max sat = 1. ! default gravity waves saturation value !! call getin_p("nonoro_gwd_sat", sat) write(*,*) "nonoro_gwd_ran: sat=", sat cmax = 50. ! default gravity waves phase velocity value !! call getin_p("nonoro_gwd_cmax", cmax) write(*,*) "nonoro_gwd_ran: cmax=", cmax rdiss=0.01 ! default coefficient of dissipation !! call getin_p("nonoro_gwd_rdiss", rdiss) write(*,*) "nonoro_gwd_ran: rdiss=", rdiss kmax=7.E-4 ! default Max horizontal wavenumber !! call getin_p("nonoro_gwd_kmax", kmax) write(*,*) "nonoro_gwd_ran: kmax=", kmax ! Characteristic vertical scale height H0 = r * tr / g ! Control if (deltat .LT. dtime) THEN call abort_physic("nonoro_gwd_ran","gwd random: deltat lower than dtime!",1) endif if (nlayer .LT. nw) THEN call abort_physic("nonoro_gwd_ran","gwd random: nlayer lower than nw!",1) endif firstcall = .false. ENDIF ! Compute current values of temperature and winds tt(:,:)=pt(:,:)+dtime*pdt(:,:) uu(:,:)=pu(:,:)+dtime*pdu(:,:) vv(:,:)=pv(:,:)+dtime*pdv(:,:) ! Compute the real mass density by rho=p/R(T)T DO ll=1,nlayer DO ii=1,ngrid rho(ii,ll) = pp(ii,ll)/(rnew(ii,ll)*tt(ii,ll)) ENDDO ENDDO ! print*,'epflux_max just after firstcall:', epflux_max !----------------------------------------------------------------------------------------------------------------------- ! 2. EVALUATION OF THE BACKGROUND FLOW AT SEMI-LEVELS !----------------------------------------------------------------------------------------------------------------------- ! Pressure and inverse of pressure at 1/2 level DO LL = 2, nlayer PH(:, LL) = EXP((LOG(PP(:, LL)) + LOG(PP(:, LL - 1))) / 2.) end DO PH(:, nlayer + 1) = 0. PH(:, 1) = 2. * PP(:, 1) - PH(:, 2) call writediagfi(ngrid,'nonoro_pp','nonoro_pp', 'm',3,PP(:,1:nlayer)) call writediagfi(ngrid,'nonoro_ph','nonoro_ph', 'm',3,PH(:,1:nlayer)) ! Launching level for reproductible case !Pour revenir a la version non reproductible en changeant le nombre de !process ! Reprend la formule qui calcule PH en fonction de PP=play DO LL = 2, nlayer HREF(LL) = EXP((LOG(presnivs(LL))+ LOG(presnivs(LL - 1))) / 2.) end DO HREF(nlayer + 1) = 0. HREF(1) = 2. * presnivs(1) - HREF(2) LAUNCH=0 DO LL = 1, nlayer IF (HREF(LL) / HREF(1) > XLAUNCH) LAUNCH = LL ENDDO ! Log pressure vert. coordinate ZH(:,1) = 0. DO LL = 2, nlayer + 1 !ZH(:, LL) = H0 * LOG(PR / (PH(:, LL) + PSEC)) H0bis(:, LL-1) = r * TT(:, LL-1) / g ZH(:, LL) = ZH(:, LL-1) & + H0bis(:, LL-1)*(PH(:, LL-1)-PH(:,LL))/PP(:, LL-1) end DO ZH(:, 1) = H0 * LOG(PR / (PH(:, 1) + PSEC)) call writediagfi(ngrid,'nonoro_zh','nonoro_zh', 'm',3,ZH(:,2:nlayer+1)) ! Winds and Brunt Vaisala frequency DO LL = 2, nlayer UH(:, LL) = 0.5 * (UU(:, LL) + UU(:, LL - 1)) ! Zonal wind VH(:, LL) = 0.5 * (VV(:, LL) + VV(:, LL - 1)) ! Meridional wind ! Brunt Vaisala frequency (=g/T*[dT/dz + g/cp] ) BV(:, LL)= G/(0.5 * (TT(:, LL) + TT(:, LL - 1))) & *((TT(:, LL) - TT(:, LL - 1)) / (ZH(:, LL) - ZH(:, LL - 1))+ & G / cpnew(:,LL)) BV(:,LL) =MAX(1.E-12,BV(:,LL)) ! to ensure that BV is positive BV(:,LL) = MAX(BVSEC,SQRT(BV(:,LL))) ! make sure it is not too small end DO BV(:, 1) = BV(:, 2) UH(:, 1) = 0. VH(:, 1) = 0. BV(:, nlayer + 1) = BV(:, nlayer) UH(:, nlayer + 1) = UU(:, nlayer) VH(:, nlayer + 1) = VV(:, nlayer) call writediagfi(ngrid,'nonoro_bv','nonoro_bv', 'm',3,BV(:,2:nlayer+1)) !----------------------------------------------------------------------------------------------------------------------- ! 3. WAVES CHARACTERISTICS CHOSEN RANDOMLY !----------------------------------------------------------------------------------------------------------------------- ! The mod function of here a weird arguments are used to produce the waves characteristics in a stochastic way DO JW = 1, NW ! Angle DO II = 1, ngrid ! Angle (0 or PI so far) RAN_NUM_1=MOD(TT(II, JW) * 10., 1.) RAN_NUM_2= MOD(TT(II, JW) * 100., 1.) ZP(JW, II) = (SIGN(1., 0.5 - RAN_NUM_1) + 1.)* PI / 2. ! Horizontal wavenumber amplitude ! From Venus model: TN+GG April/June 2020 (rev2381) ! "Individual waves are not supposed to occupy the entire domain, but only a fraction of it" (Lott+2012) ! ZK(JW, II) = KMIN + (KMAX - KMIN) *RAN_NUM_2 KSTAR = PI/SQRT(cell_area(II)) ZK(JW, II) = MAX(KMIN,KSTAR) + (KMAX - MAX(KMIN,KSTAR)) *RAN_NUM_2 ! Horizontal phase speed ! this computation allows to get a gaussian distribution centered on 0 (right ?) ! then cmin is not useful, and it favors small values. CPHA = 0. DO JJ = 1, NA RAN_NUM_3=MOD(TT(II, JW+3*JJ)**2, 1.) CPHA = CPHA + 2.*CMAX* & (RAN_NUM_3 -0.5)* & SQRT(3.)/SQRT(NA*1.) END DO IF (CPHA.LT.0.) THEN CPHA = -1.*CPHA ZP(JW,II) = ZP(JW,II) + PI ENDIF ! otherwise, with the computation below, we get a uniform distribution between cmin and cmax. ! ran_num_3 = mod(tt(ii, jw)**2, 1.) ! cpha = cmin + (cmax - cmin) * ran_num_3 ! Intrinsic frequency ZO(JW, II) = CPHA * ZK(JW, II) ! Intrinsic frequency is imposed ZO(JW, II) = ZO(JW, II) & + ZK(JW, II) * COS(ZP(JW, II)) * UH(II, LAUNCH) & + ZK(JW, II) * SIN(ZP(JW, II)) * VH(II, LAUNCH) ! Momentum flux at launch level epflux_0(JW, II) = epflux_max & * MOD(100. * (UU(II, JW)**2 + VV(II, JW)**2), 1.) ENDDO end DO ! print*,'epflux_0 just after waves charac. ramdon:', epflux_0 !----------------------------------------------------------------------------------------------------------------------- ! 4. COMPUTE THE FLUXES !----------------------------------------------------------------------------------------------------------------------- ! 4.1 Vertical velocity at launching altitude to ensure the correct value to the imposed fluxes. !------------------------------------------------------ DO JW = 1, NW ! Evaluate intrinsic frequency at launching altitude: intr_freq_p(JW, :) = ZO(JW, :) & - ZK(JW, :) * COS(ZP(JW, :)) * UH(:, LAUNCH) & - ZK(JW, :) * SIN(ZP(JW, :)) * VH(:, LAUNCH) end DO ! VERSION WITHOUT CONVECTIVE SOURCE ! Vertical velocity at launch level, value to ensure the imposed ! mom flux: DO JW = 1, NW ! WW is directly a flux, here, not vertical velocity anymore WWP(JW, :) = epflux_0(JW,:) u_epflux_p(JW, :) = COS(ZP(JW, :)) * SIGN(1., intr_freq_p(JW, :)) * epflux_0(JW, :) v_epflux_p(JW, :) = SIN(ZP(JW, :)) * SIGN(1., intr_freq_p(JW, :)) * epflux_0(JW, :) end DO ! 4.2 Initial flux at launching altitude !------------------------------------------------------ u_epflux_tot(:, LAUNCH) = 0 v_epflux_tot(:, LAUNCH) = 0 DO JW = 1, NW u_epflux_tot(:, LAUNCH) = u_epflux_tot(:, LAUNCH) + u_epflux_p(JW, :) v_epflux_tot(:, LAUNCH) = v_epflux_tot(:, LAUNCH) + v_epflux_p(JW, :) end DO ! 4.3 Loop over altitudes, with passage from one level to the next done by: !---------------------------------------------------------------------------- ! i) conserving the EP flux, ! ii) dissipating a little, ! iii) testing critical levels, ! iv) testing the breaking. !---------------------------------------------------------------------------- DO LL = LAUNCH, nlayer - 1 ! W(KB)ARNING: ALL THE PHYSICS IS HERE (PASSAGE FROM ONE LEVEL TO THE NEXT) DO JW = 1, NW intr_freq_m(JW, :) = intr_freq_p(JW, :) WWM(JW, :) = WWP(JW, :) ! Intrinsic Frequency intr_freq_p(JW, :) = ZO(JW, :) - ZK(JW, :) * COS(ZP(JW,:)) * UH(:, LL + 1) & - ZK(JW, :) * SIN(ZP(JW,:)) * VH(:, LL + 1) ! Vertical velocity in flux formulation WWP(JW, :) = MIN( & ! No breaking (Eq.6) WWM(JW, :) & ! Dissipation (Eq. 8)(real rho used here rather than pressure ! parametration because the original code has a bug if the density of ! the planet at the launch altitude not approximate 1): ! * EXP(- RDISS*2./rho(:, LL + 1) & * ((BV(:, LL + 1) + BV(:, LL)) / 2.)**3 & / MAX(ABS(intr_freq_p(JW, :) + intr_freq_m(JW, :)) / 2., ZOISEC)**4 & * ZK(JW, :)**3 * (ZH(:, LL + 1) - ZH(:, LL))), & ! Critical levels (forced to zero if intrinsic frequency ! changes sign) MAX(0., SIGN(1., intr_freq_p(JW, :) * intr_freq_m(JW, :))) & ! Saturation (Eq. 12) (rho at launch altitude is imposed by J.Liu. ! Same reason with Eq. 8) * ABS(intr_freq_p(JW, :))**3 / (2.*BV(:, LL+1)) & * rho(:,launch)*exp(-zh(:, ll + 1) / H0) & * SAT**2 *KMIN**2 / ZK(JW, :)**4) end DO ! END OF W(KB)ARNING ! Evaluate EP-flux from Eq. 7 and give the right orientation to the stress DO JW = 1, NW u_epflux_p(JW, :) = SIGN(1.,intr_freq_p(JW, :)) * COS(ZP(JW, :)) * WWP(JW, :) v_epflux_p(JW, :) = SIGN(1.,intr_freq_p(JW, :)) * SIN(ZP(JW, :)) * WWP(JW, :) end DO u_epflux_tot(:, LL + 1) = 0. v_epflux_tot(:, LL + 1) = 0. DO JW = 1, NW u_epflux_tot(:, LL + 1) = u_epflux_tot(:, LL + 1) + u_epflux_p(JW, :) v_epflux_tot(:, LL + 1) = v_epflux_tot(:, LL + 1) + v_epflux_p(JW, :) EAST_GWSTRESS(:, LL)=EAST_GWSTRESS(:, LL)+MAX(0.,u_epflux_p(JW,:))/FLOAT(NW) WEST_GWSTRESS(:, LL)=WEST_GWSTRESS(:, LL)+MIN(0.,u_epflux_p(JW,:))/FLOAT(NW) end DO end DO ! DO LL = LAUNCH, nlayer - 1 !----------------------------------------------------------------------------------------------------------------------- ! 5. TENDENCY CALCULATIONS !----------------------------------------------------------------------------------------------------------------------- ! 5.1 Flow rectification at the top and in the low layers ! -------------------------------------------------------- ! Warning, this is the total on all GW u_epflux_tot(:, nlayer + 1) = 0. v_epflux_tot(:, nlayer + 1) = 0. ! Here, big change compared to FLott version: ! We compensate (u_epflux_tot(:, LAUNCH), ie total emitted upward flux ! over the layers max(1,LAUNCH-3) to LAUNCH-1 DO LL = 1, max(1,LAUNCH-3) u_epflux_tot(:, LL) = 0. v_epflux_tot(:, LL) = 0. end DO DO LL = max(2,LAUNCH-2), LAUNCH-1 u_epflux_tot(:, LL) = u_epflux_tot(:, LL - 1) + u_epflux_tot(:, LAUNCH) * & (PH(:,LL)-PH(:,LL-1)) / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3))) v_epflux_tot(:, LL) = v_epflux_tot(:, LL - 1) + v_epflux_tot(:, LAUNCH) * & (PH(:,LL)-PH(:,LL-1)) / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3))) EAST_GWSTRESS(:,LL) = EAST_GWSTRESS(:, LL - 1) + & EAST_GWSTRESS(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1))/ & (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3))) WEST_GWSTRESS(:,LL) = WEST_GWSTRESS(:, LL - 1) + & WEST_GWSTRESS(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1))/ & (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3))) end DO ! This way, the total flux from GW is zero, but there is a net transport ! (upward) that should be compensated by circulation ! and induce additional friction at the surface call writediagfi(ngrid,'nonoro_u_epflux_tot','nonoro_u_epflux_tot', '',3,u_epflux_tot(:,2:nlayer+1)) call writediagfi(ngrid,'nonoro_v_epflux_tot','nonoro_v_epflux_tot', '',3,v_epflux_tot(:,2:nlayer+1)) ! 5.2 AR-1 RECURSIVE FORMULA (13) IN VERSION 4 !--------------------------------------------- DO LL = 1, nlayer !Notice here the d_u and d_v are tendency (For Mars) but not increment(Venus). d_u(:, LL) = G * (u_epflux_tot(:, LL + 1) - u_epflux_tot(:, LL)) & / (PH(:, LL + 1) - PH(:, LL)) d_v(:, LL) = G * (v_epflux_tot(:, LL + 1) - v_epflux_tot(:, LL)) & / (PH(:, LL + 1) - PH(:, LL)) ENDDO d_t(:,:) = 0. call writediagfi(ngrid,'nonoro_d_u','nonoro_d_u', '',3,d_u) call writediagfi(ngrid,'nonoro_d_v','nonoro_d_v', '',3,d_v) ! 5.3 Update tendency of wind with the previous (and saved) values !----------------------------------------------------------------- d_u(:,:) = DTIME/DELTAT/REAL(NW) * d_u(:,:) & + (1.-DTIME/DELTAT) * du_nonoro_gwd(:,:) d_v(:,:) = DTIME/DELTAT/REAL(NW) * d_v(:,:) & + (1.-DTIME/DELTAT) * dv_nonoro_gwd(:,:) du_nonoro_gwd(:,:) = d_u(:,:) dv_nonoro_gwd(:,:) = d_v(:,:) call writediagfi(ngrid,'du_nonoro_gwd','du_nonoro_gwd', '',3,du_nonoro_gwd) call writediagfi(ngrid,'dv_nonoro_gwd','dv_nonoro_gwd', '',3,dv_nonoro_gwd) ! Cosmetic: evaluation of the cumulated stress ZUSTR(:) = 0. ZVSTR(:) = 0. DO LL = 1, nlayer ZUSTR(:) = ZUSTR(:) + D_U(:, LL) *rho(:, LL) * (ZH(:, LL + 1) - ZH(:, LL))*DTIME ZVSTR(:) = ZVSTR(:) + D_V(:, LL) *rho(:, LL) * (ZH(:, LL + 1) - ZH(:, LL))*DTIME ENDDO !#ifdef CPP_XIOS ! call send_xios_field("du_nonoro", d_u) ! call send_xios_field("dv_nonoro", d_v) ! call send_xios_field("east_gwstress", east_gwstress) ! call send_xios_field("west_gwstress", west_gwstress) !#endif END SUBROUTINE NONORO_GWD_RAN ! ======================================================== ! Subroutines used to allocate/deallocate module variables ! ======================================================== SUBROUTINE ini_nonoro_gwd_ran(ngrid,nlayer) IMPLICIT NONE INTEGER, INTENT (in) :: ngrid ! number of atmospheric columns INTEGER, INTENT (in) :: nlayer ! number of atmospheric layers allocate(du_nonoro_gwd(ngrid,nlayer)) allocate(dv_nonoro_gwd(ngrid,nlayer)) allocate(east_gwstress(ngrid,nlayer)) east_gwstress(:,:)=0 allocate(west_gwstress(ngrid,nlayer)) west_gwstress(:,:)=0 END SUBROUTINE ini_nonoro_gwd_ran ! ---------------------------------- SUBROUTINE end_nonoro_gwd_ran IMPLICIT NONE if (allocated(du_nonoro_gwd)) deallocate(du_nonoro_gwd) if (allocated(dv_nonoro_gwd)) deallocate(dv_nonoro_gwd) if (allocated(east_gwstress)) deallocate(east_gwstress) if (allocated(west_gwstress)) deallocate(west_gwstress) END SUBROUTINE end_nonoro_gwd_ran END MODULE nonoro_gwd_ran_mod