source: trunk/LMDZ.MARS/libf/phymars/nonoro_gwd_ran_mod.F90 @ 3799

Last change on this file since 3799 was 3726, checked in by emillour, 2 months ago

Mars PCM:
Turn "callkeys.h" into module "callkeys_mod.F90"
EM

File size: 25.4 KB
Line 
1MODULE nonoro_gwd_ran_mod
2
3IMPLICIT NONE
4
5REAL,allocatable,save :: du_nonoro_gwd(:,:) ! Zonal wind tendency due to GWD
6REAL,allocatable,save :: dv_nonoro_gwd(:,:) ! Meridional wind tendency due to GWD
7REAL,ALLOCATABLE,SAVE :: east_gwstress(:,:) ! Profile of eastward stress
8REAL,ALLOCATABLE,SAVE :: west_gwstress(:,:) ! Profile of westward stress
9
10!$OMP THREADPRIVATE(du_nonoro_gwd,dv_nonoro_gwd,east_gwstress,west_gwstress)
11
12CONTAINS
13
14      SUBROUTINE NONORO_GWD_RAN(ngrid,nlayer,DTIME, cpnew, rnew, pp,  &
15                  zmax_therm, pt, pu, pv, pdt, pdu, pdv, &
16                  zustr,zvstr,d_t, d_u, d_v)
17
18    !--------------------------------------------------------------------------------
19    ! Parametrization of the momentum flux deposition due to a discrete
20    ! number of gravity waves.
21    ! F. Lott
22    ! Version 14, Gaussian distribution of the source
23    ! LMDz model online version     
24    ! ADAPTED FOR VENUS /  F. LOTT + S. LEBONNOIS
25    ! Version adapted on 03/04/2013:
26    !      - input flux compensated in the deepest layers
27    !                           
28    ! ADAPTED FOR MARS     G.GILLI     02/2016
29    !        Revision with F.Forget    06/2016  Variable EP-flux according to
30    !                                           PBL variation (max velocity thermals)
31    ! UPDATED              D.BARDET    01/2020  - reproductibility of the
32    !                                           launching altitude calculation
33    !                                           - wave characteristic
34    !                                           calculation using MOD
35    !                                           - adding east_gwstress and
36    !                                           west_gwstress variables
37    ! UPDATED              J.LIU       12/2021  The rho (density) at the specific
38    !                                           locations is introduced. The equation
39    !                                           of EP-flux is corrected by adding the
40    !                                           term of density at launch (source)
41    !                                           altitude.
42    ! 
43    !---------------------------------------------------------------------------------
44
45      use comcstfi_h, only: g, pi, r
46      USE ioipsl_getin_p_mod, ONLY : getin_p
47      use vertical_layers_mod, only : presnivs
48      use geometry_mod, only: cell_area
49      use write_output_mod, only: write_output
50     
51      implicit none
52
53      CHARACTER (LEN=20) :: modname='nonoro_gwd_ran'
54      CHARACTER (LEN=80) :: abort_message
55
56
57    ! 0. DECLARATIONS:
58
59    ! 0.1 INPUTS
60    INTEGER, intent(in):: ngrid          ! number of atmospheric columns
61    INTEGER, intent(in):: nlayer         ! number of atmospheric layers
62    REAL, intent(in):: DTIME             ! Time step of the Physics(s)
63    REAL, intent(in):: zmax_therm(ngrid) ! Altitude of max velocity thermals (m)
64    REAL,INTENT(IN) :: cpnew(ngrid,nlayer)! Cp of the atmosphere
65    REAL,INTENT(IN) :: rnew(ngrid,nlayer) ! R of the atmosphere
66    REAL, intent(in):: pp(ngrid,nlayer)  ! Pressure at full levels(Pa)
67    REAL, intent(in):: pt(ngrid,nlayer)  ! Temp at full levels(K)
68    REAL, intent(in):: pu(ngrid,nlayer)  ! Zonal wind at full levels(m/s)
69    REAL, intent(in):: pv(ngrid,nlayer)  ! Meridional winds at full levels(m/s)
70    REAL,INTENT(in) :: pdt(ngrid,nlayer) ! Tendency on temperature (K/s)
71    REAL,INTENT(in) :: pdu(ngrid,nlayer) ! Tendency on zonal wind (m/s/s)
72    REAL,INTENT(in) :: pdv(ngrid,nlayer) ! Tendency on meridional wind (m/s/s)
73
74    ! 0.2 OUTPUTS
75    REAL, intent(out):: zustr(ngrid)       ! Zonal surface stress
76    REAL, intent(out):: zvstr(ngrid)       ! Meridional surface stress
77    REAL, intent(out):: d_t(ngrid, nlayer) ! Tendency on temperature (K/s) due to gravity waves (not used set to zero)
78    REAL, intent(out):: d_u(ngrid, nlayer) ! Tendency on zonal wind (m/s/s) due to gravity waves
79    REAL, intent(out):: d_v(ngrid, nlayer) ! Tendency on meridional wind (m/s/s) due to gravity waves
80
81    ! 0.3 INTERNAL ARRAYS
82    REAL :: TT(ngrid, nlayer)   ! Temperature at full levels
83    REAL :: RHO(ngrid, nlayer)  ! Mass density at full levels
84    REAL :: UU(ngrid, nlayer)   ! Zonal wind at full levels
85    REAL :: VV(ngrid, nlayer)   ! Meridional winds at full levels
86    REAL :: BVLOW(ngrid)        ! initial N at each grid (not used)
87
88    INTEGER II, JJ, LL
89
90    ! 0.3.0 TIME SCALE OF THE LIFE CYCLE OF THE WAVES PARAMETERIZED
91    REAL, parameter:: DELTAT = 24. * 3600.
92   
93    ! 0.3.1 GRAVITY-WAVES SPECIFICATIONS
94    INTEGER, PARAMETER:: NK = 2 ! number of horizontal wavenumbers
95    INTEGER, PARAMETER:: NP = 2 ! directions (eastward and westward) phase speed
96    INTEGER, PARAMETER:: NO = 2 ! absolute values of phase speed
97    INTEGER, PARAMETER:: NA = 5 ! number of realizations to get the phase speed
98    INTEGER, PARAMETER:: NW = NK * NP * NO ! Total numbers of gravity waves
99    INTEGER JK, JP, JO, JW      ! Loop index for NK,NP,NO, and NW
100
101    REAL, save :: kmax             ! Max horizontal wavenumber (lambda min,lambda=2*pi/kmax),kmax=N/u, u(=30~50) zonal wind velocity at launch altitude
102!$OMP THREADPRIVATE(kmax)
103    REAL, save :: kmin             ! Min horizontal wavenumber (lambda max = 314 km,lambda=2*pi/kmin)
104!$OMP THREADPRIVATE(kmin)
105    REAL kstar                     ! Min horizontal wavenumber constrained by horizontal resolution
106
107    REAL :: max_k(ngrid)           ! max_k=max(kstar,kmin)
108    REAL, parameter:: cmin = 1.    ! Min horizontal absolute phase velocity (not used)   
109    REAL CPHA(ngrid)                      ! absolute PHASE VELOCITY frequency
110    REAL ZK(NW, ngrid)             ! Horizontal wavenumber amplitude
111    REAL ZP(NW, ngrid)             ! Horizontal wavenumber angle       
112    REAL ZO(NW, ngrid)             ! Absolute frequency
113
114    REAL intr_freq_m(nw, ngrid)          ! Waves Intr. freq. at the 1/2 lev below the full level (previous name: ZOM)
115    REAL intr_freq_p(nw, ngrid)          ! Waves Intr. freq. at the 1/2 lev above the full level (previous name: ZOP)
116    REAL wwm(nw, ngrid)                  ! Wave EP-fluxes at the 1/2 level below the full level
117    REAL wwp(nw, ngrid)                  ! Wave EP-fluxes at the 1/2 level above the full level
118    REAL u_epflux_p(nw, ngrid)           ! Partial zonal flux (=for each wave) at the 1/2 level above the full level (previous name: RUWP)
119    REAL v_epflux_p(nw, ngrid)           ! Partial meridional flux (=for each wave) at the 1/2 level above the full level (previous name: RVWP)
120    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) 
121    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)
122    REAL epflux_0(nw, ngrid)             ! Fluxes at launching level (previous name: RUW0)
123    REAL, save :: epflux_max             ! Max EP flux value at launching altitude (previous name: RUWMAX, tunable)
124!$OMP THREADPRIVATE(epflux_max)
125    INTEGER LAUNCH                       ! Launching altitude
126    REAL, save :: xlaunch                ! Control the launching altitude by pressure
127!$OMP THREADPRIVATE(xlaunch)
128    REAL, parameter:: zmaxth_top = 8000. ! Top of convective layer (approx. not used)
129    REAL cmax(ngrid,nlayer)             ! Max horizontal absolute phase velocity (at the maginitide of zonal wind u at the launch altitude)
130
131    ! 0.3.2 PARAMETERS OF WAVES DISSIPATIONS
132    REAL, save :: sat                 ! saturation parameter(tunable)
133!$OMP THREADPRIVATE(sat)
134    REAL, save :: rdiss               ! dissipation coefficient (tunable)
135!$OMP THREADPRIVATE(rdiss)
136    REAL, parameter:: zoisec = 1.e-10 ! security for intrinsic frequency
137
138    ! 0.3.3 Background flow at 1/2 levels and vertical coordinate
139    REAL H0bis(ngrid, nlayer)          ! Characteristic Height of the atmosphere (specific locations)
140    REAL, save:: H0                    ! Characteristic Height of the atmosphere (constant)
141!$OMP THREADPRIVATE(H0)
142    REAL, parameter:: pr = 250         ! Reference pressure [Pa]
143    REAL, parameter:: tr = 220.        ! Reference temperature [K]
144    REAL ZH(ngrid, nlayer + 1)         ! Log-pressure altitude (constant H0)
145    REAL ZHbis(ngrid, nlayer + 1)      ! Log-pressure altitude (varying H0bis)
146    REAL UH(ngrid, nlayer + 1)         ! zonal wind at 1/2 levels
147    REAL VH(ngrid, nlayer + 1)         ! meridional wind at 1/2 levels
148    REAL PH(ngrid, nlayer + 1)         ! Pressure at 1/2 levels
149    REAL, parameter:: psec = 1.e-20    ! Security to avoid division by 0 pressure(!!IMPORTANT: should be lower than the topmost layer's pressure)
150    REAL BV(ngrid, nlayer + 1)         ! Brunt Vaisala freq. (N^2) at 1/2 levels
151    REAL, parameter:: bvsec = 1.e-5    ! Security to avoid negative BV 
152    REAL HREF(nlayer + 1)              ! Reference pressure for launching altitude
153
154
155    REAL RAN_NUM_1,RAN_NUM_2,RAN_NUM_3
156
157
158    LOGICAL,SAVE :: firstcall = .true.
159!$OMP THREADPRIVATE(firstcall)
160
161!-----------------------------------------------------------------------------------------------------------------------
162!  1. INITIALISATIONS
163!-----------------------------------------------------------------------------------------------------------------------
164     IF (firstcall) THEN
165        write(*,*) "nonoro_gwd_ran: FLott non-oro GW scheme is active!"
166        epflux_max = 5.E-4 ! Mars' value !!
167        call getin_p("nonoro_gwd_epflux_max", epflux_max)
168        write(*,*) "nonoro_gwd_ran: epflux_max=", epflux_max
169        sat = 1.5 ! default gravity waves saturation value !!
170        call getin_p("nonoro_gwd_sat", sat)
171        write(*,*) "nonoro_gwd_ran: sat=", sat     
172!        cmax = 50. ! default gravity waves phase velocity value !!
173!        call getin_p("nonoro_gwd_cmax", cmax)
174!        write(*,*) "nonoro_gwd_ran: cmax=", cmax
175        rdiss=0.07 ! default coefficient of dissipation !!
176        call getin_p("nonoro_gwd_rdiss", rdiss)
177        write(*,*) "nonoro_gwd_ran: rdiss=", rdiss
178        kmax=1.E-4 ! default Max horizontal wavenumber !!
179        call getin_p("nonoro_gwd_kmax", kmax)
180        write(*,*) "nonoro_gwd_ran: kmax=", kmax
181        kmin=7.E-6 ! default Max horizontal wavenumber !!
182        call getin_p("nonoro_gwd_kmin", kmin)
183        write(*,*) "nonoro_gwd_ran: kmin=", kmin
184        xlaunch=0.6 ! default GW luanch altitude controller !!
185        call getin_p("nonoro_gwd_xlaunch", xlaunch)
186        write(*,*) "nonoro_gwd_ran: xlaunch=", xlaunch
187        ! Characteristic vertical scale height
188        H0 = r * tr / g
189        ! Control
190        if (deltat .LT. dtime) THEN
191             call abort_physic("nonoro_gwd_ran","gwd random: deltat lower than dtime!",1)
192        endif
193        if (nlayer .LT. nw) THEN
194             call abort_physic("nonoro_gwd_ran","gwd random: nlayer lower than nw!",1)
195        endif
196        firstcall = .false.
197     ENDIF
198
199! Compute current values of temperature and winds
200    tt(:,:)=pt(:,:)+dtime*pdt(:,:)
201    uu(:,:)=pu(:,:)+dtime*pdu(:,:)
202    vv(:,:)=pv(:,:)+dtime*pdv(:,:)
203! Compute the real mass density by rho=p/R(T)T
204     DO ll=1,nlayer
205       DO ii=1,ngrid
206            rho(ii,ll) = pp(ii,ll)/(rnew(ii,ll)*tt(ii,ll))
207       ENDDO
208     ENDDO
209!    print*,'epflux_max just after firstcall:', epflux_max
210
211!-----------------------------------------------------------------------------------------------------------------------
212!  2. EVALUATION OF THE BACKGROUND FLOW AT SEMI-LEVELS
213!-----------------------------------------------------------------------------------------------------------------------
214    ! Pressure and inverse of pressure at 1/2 level
215    DO LL = 2, nlayer
216       PH(:, LL) = EXP((LOG(PP(:, LL)) + LOG(PP(:, LL - 1))) / 2.)
217    end DO
218    PH(:, nlayer + 1) = 0.
219    PH(:, 1) = 2. * PP(:, 1) - PH(:, 2)
220
221!    call write_output('nonoro_pp','nonoro_pp', 'm',PP(:,:))
222!    call write_output('nonoro_ph','nonoro_ph', 'm',PH(:,:))
223
224    ! Launching level for reproductible case
225    !Pour revenir a la version non reproductible en changeant le nombre de
226    !process
227    ! Reprend la formule qui calcule PH en fonction de PP=play
228    DO LL = 2, nlayer
229       HREF(LL) = EXP((LOG(presnivs(LL))+ LOG(presnivs(LL - 1))) / 2.)
230    end DO
231    HREF(nlayer + 1) = 0.
232    HREF(1) = 2. * presnivs(1) - HREF(2)
233
234    LAUNCH=0
235    DO LL = 1, nlayer
236       IF (HREF(LL) / HREF(1) > XLAUNCH) LAUNCH = LL
237    ENDDO
238       
239    ! Log pressure vert. coordinate
240       ZH(:,1) = 0.
241    DO LL = 2, nlayer + 1
242       !ZH(:, LL) = H0 * LOG(PR / (PH(:, LL) + PSEC))
243        H0bis(:, LL-1) = r * TT(:, LL-1) / g
244          ZH(:, LL) = ZH(:, LL-1) &
245           + H0bis(:, LL-1)*(PH(:, LL-1)-PH(:,LL))/PP(:, LL-1)
246    end DO
247        ZH(:, 1) = H0 * LOG(PR / (PH(:, 1) + PSEC))
248
249!    call write_output('nonoro_zh','nonoro_zh', 'm',ZH(:,2:nlayer+1))
250
251    ! Winds and Brunt Vaisala frequency
252    DO LL = 2, nlayer
253       UH(:, LL) = 0.5 * (UU(:, LL) + UU(:, LL - 1))                          ! Zonal wind
254       VH(:, LL) = 0.5 * (VV(:, LL) + VV(:, LL - 1))                          ! Meridional wind
255       ! Brunt Vaisala frequency (=g/T*[dT/dz + g/cp] )
256       BV(:, LL)= G/(0.5 * (TT(:, LL) + TT(:, LL - 1)))                     &
257        *((TT(:, LL) - TT(:, LL - 1)) / (ZH(:, LL) - ZH(:, LL - 1))+ &
258        G / cpnew(:,LL))
259
260       BV(:,LL) =MAX(1.E-12,BV(:,LL)) ! to ensure that BV is positive
261       BV(:,LL) = MAX(BVSEC,SQRT(BV(:,LL))) ! make sure it is not too small
262    end DO
263    BV(:, 1) = BV(:, 2)
264    UH(:, 1) = 0.
265    VH(:, 1) = 0.
266    BV(:, nlayer + 1) = BV(:, nlayer)
267    UH(:, nlayer + 1) = UU(:, nlayer)
268    VH(:, nlayer + 1) = VV(:, nlayer)
269    cmax(:,launch)=UU(:, launch)
270    DO ii=1,ngrid
271       KSTAR = PI/SQRT(cell_area(II))
272       MAX_K(II)=MAX(kmin,kstar)
273    ENDDO
274    call write_output('nonoro_bv','Brunt Vaisala frequency in nonoro', 'Hz',BV(:,2:nlayer+1))
275
276!-----------------------------------------------------------------------------------------------------------------------
277! 3. WAVES CHARACTERISTICS CHOSEN RANDOMLY
278!-----------------------------------------------------------------------------------------------------------------------
279! The mod function of here a weird arguments are used to produce the waves characteristics in a stochastic way
280    DO JW = 1, NW
281             !  Angle
282             DO II = 1, ngrid
283                ! Angle (0 or PI so far)
284                RAN_NUM_1=MOD(TT(II, JW) * 10., 1.)
285                RAN_NUM_2= MOD(TT(II, JW) * 100., 1.)
286                ZP(JW, II) = (SIGN(1., 0.5 - RAN_NUM_1) + 1.)* PI / 2.
287               
288                ! Horizontal wavenumber amplitude
289                ! From Venus model: TN+GG April/June 2020 (rev2381)
290                ! "Individual waves are not supposed to occupy the entire domain, but only a fraction of it" (Lott+2012)
291                ! ZK(JW, II) = KMIN + (KMAX - KMIN) *RAN_NUM_2
292                KSTAR = PI/SQRT(cell_area(II))         
293                ZK(JW, II) = MAX_K(II) + (KMAX - MAX_K(II)) *RAN_NUM_2
294               
295               ! Horizontal phase speed
296! this computation allows to get a gaussian distribution centered on 0 (right ?)
297! then cmin is not useful, and it favors small values.
298                CPHA(:) = 0.
299                DO JJ = 1, NA
300                    RAN_NUM_3=MOD(TT(II, JW+3*JJ)**2, 1.)
301                    CPHA(ii) = CPHA(ii) + 2.*CMAX(ii,launch)*      &
302                           (RAN_NUM_3 -0.5)*                       &
303                           SQRT(3.)/SQRT(NA*1.)
304                END DO
305                IF (CPHA(ii).LT.0.)  THEN
306                   CPHA(ii) = -1.*CPHA(ii)
307                   ZP(JW,II) = ZP(JW,II) + PI
308                ENDIF
309! otherwise, with the computation below, we get a uniform distribution between cmin and cmax.
310!           ran_num_3 = mod(tt(ii, jw)**2, 1.)
311!           cpha = cmin + (cmax - cmin) * ran_num_3
312
313                ! Intrinsic frequency
314                ZO(JW, II) = CPHA(II) * ZK(JW, II)
315                ! Intrinsic frequency  is imposed
316                ZO(JW, II) = ZO(JW, II)                                            &
317                            + ZK(JW, II) * COS(ZP(JW, II)) * UH(II, LAUNCH)        &
318                            + ZK(JW, II) * SIN(ZP(JW, II)) * VH(II, LAUNCH)
319               
320                ! Momentum flux at launch level
321                epflux_0(JW, II) = epflux_max                                      &
322                                  * MOD(100. * (UU(II, JW)**2 + VV(II, JW)**2), 1.)
323              ENDDO
324   end DO
325!     print*,'epflux_0 just after waves charac. ramdon:', epflux_0
326
327!-----------------------------------------------------------------------------------------------------------------------
328! 4. COMPUTE THE FLUXES
329!-----------------------------------------------------------------------------------------------------------------------
330    !  4.1  Vertical velocity at launching altitude to ensure the correct value to the imposed fluxes.
331    !------------------------------------------------------
332    DO JW = 1, NW
333       ! Evaluate intrinsic frequency at launching altitude:
334       intr_freq_p(JW, :) = ZO(JW, :)                                    &
335                            - ZK(JW, :) * COS(ZP(JW, :)) * UH(:, LAUNCH) &
336                            - ZK(JW, :) * SIN(ZP(JW, :)) * VH(:, LAUNCH)
337    end DO
338
339! VERSION WITHOUT CONVECTIVE SOURCE
340       ! Vertical velocity at launch level, value to ensure the imposed
341       ! mom flux:
342         DO JW = 1, NW
343       ! WW is directly a flux, here, not vertical velocity anymore
344            WWP(JW, :) = epflux_0(JW,:)
345            u_epflux_p(JW, :) = COS(ZP(JW, :)) * SIGN(1., intr_freq_p(JW, :)) * epflux_0(JW, :)
346            v_epflux_p(JW, :) = SIN(ZP(JW, :)) * SIGN(1., intr_freq_p(JW, :)) * epflux_0(JW, :)
347         end DO
348   
349    !  4.2 Initial flux at launching altitude
350    !------------------------------------------------------
351    u_epflux_tot(:, LAUNCH) = 0
352    v_epflux_tot(:, LAUNCH) = 0
353    DO JW = 1, NW
354       u_epflux_tot(:, LAUNCH) = u_epflux_tot(:, LAUNCH) + u_epflux_p(JW, :)
355       v_epflux_tot(:, LAUNCH) = v_epflux_tot(:, LAUNCH) + v_epflux_p(JW, :)
356    end DO
357
358    !  4.3 Loop over altitudes, with passage from one level to the next done by:
359    !----------------------------------------------------------------------------
360    !    i) conserving the EP flux,
361    !    ii) dissipating a little,
362    !    iii) testing critical levels,
363    !    iv) testing the breaking.
364    !----------------------------------------------------------------------------
365    DO LL = LAUNCH, nlayer - 1
366       !  W(KB)ARNING: ALL THE PHYSICS IS HERE (PASSAGE FROM ONE LEVEL TO THE NEXT)
367       DO JW = 1, NW
368          intr_freq_m(JW, :) = intr_freq_p(JW, :)
369          WWM(JW, :) = WWP(JW, :)
370          ! Intrinsic Frequency
371          intr_freq_p(JW, :) = ZO(JW, :) - ZK(JW, :) * COS(ZP(JW,:)) * UH(:, LL + 1)     &
372                               - ZK(JW, :) * SIN(ZP(JW,:)) * VH(:, LL + 1)
373          ! Vertical velocity in flux formulation
374          WWP(JW, :) = MIN(                                                              &
375                     ! No breaking (Eq.6)
376                     WWM(JW, :)                                                          &
377                     ! Dissipation (Eq. 8)(real rho used here rather than pressure
378                     ! parametration because the original code has a bug if the density of
379                     ! the planet at the launch altitude not approximate 1):                     !
380                     * EXP(- RDISS*2./rho(:, LL + 1)                                     &
381                     * ((BV(:, LL + 1) + BV(:, LL)) / 2.)**3                             &
382                     / MAX(ABS(intr_freq_p(JW, :) + intr_freq_m(JW, :)) / 2., ZOISEC)**4 &
383                     * ZK(JW, :)**3 * (ZH(:, LL + 1) - ZH(:, LL))),                      &
384                     ! Critical levels (forced to zero if intrinsic frequency
385                     ! changes sign)
386                     MAX(0., SIGN(1., intr_freq_p(JW, :) * intr_freq_m(JW, :)))          &
387                     ! Saturation (Eq. 12) (rho at launch altitude is imposed by J.Liu.
388                     ! Same reason with Eq. 8)
389                     * ABS(intr_freq_p(JW, :))**3 / (2.*BV(:, LL+1))                     &
390                     * rho(:,launch)*exp(-zh(:, ll + 1) / H0)                            &
391                     * SAT**2 *KMIN**2 / ZK(JW, :)**4) 
392       end DO
393       ! END OF W(KB)ARNING
394
395       ! Evaluate EP-flux from Eq. 7 and give the right orientation to the stress
396       DO JW = 1, NW
397          u_epflux_p(JW, :) = SIGN(1.,intr_freq_p(JW, :)) * COS(ZP(JW, :)) *  WWP(JW, :)
398          v_epflux_p(JW, :) = SIGN(1.,intr_freq_p(JW, :)) * SIN(ZP(JW, :)) *  WWP(JW, :)
399       end DO
400 
401       u_epflux_tot(:, LL + 1) = 0.
402       v_epflux_tot(:, LL + 1) = 0.
403       DO JW = 1, NW
404          u_epflux_tot(:, LL + 1) = u_epflux_tot(:, LL + 1) + u_epflux_p(JW, :)
405          v_epflux_tot(:, LL + 1) = v_epflux_tot(:, LL + 1) + v_epflux_p(JW, :)
406          EAST_GWSTRESS(:, LL)=EAST_GWSTRESS(:, LL)+MAX(0.,u_epflux_p(JW,:))/FLOAT(NW)
407          WEST_GWSTRESS(:, LL)=WEST_GWSTRESS(:, LL)+MIN(0.,u_epflux_p(JW,:))/FLOAT(NW)
408       end DO       
409    end DO ! DO LL = LAUNCH, nlayer - 1
410
411!-----------------------------------------------------------------------------------------------------------------------
412! 5. TENDENCY CALCULATIONS
413!-----------------------------------------------------------------------------------------------------------------------
414
415    ! 5.1 Flow rectification at the top and in the low layers
416    ! --------------------------------------------------------
417    ! Warning, this is the total on all GW
418    u_epflux_tot(:, nlayer + 1) = 0.
419    v_epflux_tot(:, nlayer + 1) = 0.
420
421    ! Here, big change compared to FLott version:
422    ! We compensate (u_epflux_tot(:, LAUNCH), ie total emitted upward flux
423    !  over the layers max(1,LAUNCH-3) to LAUNCH-1
424    DO LL = 1, max(1,LAUNCH-3)
425      u_epflux_tot(:, LL) = 0.
426      v_epflux_tot(:, LL) = 0.
427    end DO
428    DO LL = max(2,LAUNCH-2), LAUNCH-1
429       u_epflux_tot(:, LL) = u_epflux_tot(:, LL - 1) + u_epflux_tot(:, LAUNCH) *          &
430                             (PH(:,LL)-PH(:,LL-1)) / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
431       v_epflux_tot(:, LL) = v_epflux_tot(:, LL - 1) + v_epflux_tot(:, LAUNCH) *          &
432                             (PH(:,LL)-PH(:,LL-1)) / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
433       EAST_GWSTRESS(:,LL) = EAST_GWSTRESS(:, LL - 1) +                                   &
434                             EAST_GWSTRESS(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1))/            &
435                             (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
436       WEST_GWSTRESS(:,LL) = WEST_GWSTRESS(:, LL - 1) +                                   &
437                             WEST_GWSTRESS(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1))/            &
438                             (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
439    end DO
440    ! This way, the total flux from GW is zero, but there is a net transport
441    ! (upward) that should be compensated by circulation
442    ! and induce additional friction at the surface
443    call write_output('nonoro_u_epflux_tot','Total EP Flux along U in nonoro', '',u_epflux_tot(:,2:nlayer+1))
444    call write_output('nonoro_v_epflux_tot','Total EP Flux along V in nonoro', '',v_epflux_tot(:,2:nlayer+1))
445
446    ! 5.2 AR-1 RECURSIVE FORMULA (13) IN VERSION 4
447    !---------------------------------------------
448    DO LL = 1, nlayer
449       !Notice here the d_u and d_v are tendency (For Mars) but not increment(Venus).
450       d_u(:, LL) = G * (u_epflux_tot(:, LL + 1) - u_epflux_tot(:, LL)) &
451                    / (PH(:, LL + 1) - PH(:, LL))
452       d_v(:, LL) = G * (v_epflux_tot(:, LL + 1) - v_epflux_tot(:, LL)) &
453                    / (PH(:, LL + 1) - PH(:, LL))
454    ENDDO
455    d_t(:,:) = 0.
456!    call write_output('nonoro_d_u','nonoro_d_u', '',d_u(:,:))
457!    call write_output('nonoro_d_v','nonoro_d_v', '',d_v(:,:))
458
459    ! 5.3 Update tendency of wind with the previous (and saved) values
460    !-----------------------------------------------------------------
461    d_u(:,:) = DTIME/DELTAT/REAL(NW) * d_u(:,:)                       &
462               + (1.-DTIME/DELTAT) * du_nonoro_gwd(:,:)
463    d_v(:,:) = DTIME/DELTAT/REAL(NW) * d_v(:,:)                       &
464               + (1.-DTIME/DELTAT) * dv_nonoro_gwd(:,:)
465    du_nonoro_gwd(:,:) = d_u(:,:)
466    dv_nonoro_gwd(:,:) = d_v(:,:)
467
468    call write_output('du_nonoro_gwd','Tendency on U due to nonoro GW', 'm.s-2',du_nonoro_gwd(:,:))
469    call write_output('dv_nonoro_gwd','Tendency on V due to nonoro GW', 'm.s-2',dv_nonoro_gwd(:,:))
470   
471    ! Cosmetic: evaluation of the cumulated stress
472    ZUSTR(:) = 0.
473    ZVSTR(:) = 0.
474    DO LL = 1, nlayer
475       ZUSTR(:) = ZUSTR(:) + D_U(:, LL) *rho(:, LL) * (ZH(:, LL + 1) - ZH(:, LL))*DTIME
476       ZVSTR(:) = ZVSTR(:) + D_V(:, LL) *rho(:, LL) * (ZH(:, LL + 1) - ZH(:, LL))*DTIME
477    ENDDO
478
479  END SUBROUTINE NONORO_GWD_RAN
480
481
482
483! ========================================================
484! Subroutines used to allocate/deallocate module variables       
485! ========================================================
486  SUBROUTINE ini_nonoro_gwd_ran(ngrid,nlayer)
487
488  IMPLICIT NONE
489
490      INTEGER, INTENT (in) :: ngrid  ! number of atmospheric columns
491      INTEGER, INTENT (in) :: nlayer ! number of atmospheric layers
492
493         allocate(du_nonoro_gwd(ngrid,nlayer))
494         allocate(dv_nonoro_gwd(ngrid,nlayer))
495         allocate(east_gwstress(ngrid,nlayer))
496         east_gwstress(:,:)=0
497         allocate(west_gwstress(ngrid,nlayer))
498         west_gwstress(:,:)=0
499
500  END SUBROUTINE ini_nonoro_gwd_ran
501! ----------------------------------
502  SUBROUTINE end_nonoro_gwd_ran
503
504  IMPLICIT NONE
505
506         if (allocated(du_nonoro_gwd)) deallocate(du_nonoro_gwd)
507         if (allocated(dv_nonoro_gwd)) deallocate(dv_nonoro_gwd)
508         if (allocated(east_gwstress)) deallocate(east_gwstress)
509         if (allocated(west_gwstress)) deallocate(west_gwstress)
510
511  END SUBROUTINE end_nonoro_gwd_ran
512
513END MODULE nonoro_gwd_ran_mod
Note: See TracBrowser for help on using the repository browser.