source: trunk/LMDZ.VENUS/libf/phyvenus/nonoro_gwd_ran_mod.F90 @ 2745

Last change on this file since 2745 was 2686, checked in by slebonnois, 3 years ago

SL+AM : various corrections prior to VCD 2.0

File size: 22.1 KB
Line 
1MODULE nonoro_gwd_ran_mod
2
3IMPLICIT NONE
4
5      REAL, allocatable, save :: du_nonoro_gwd(:, :) ! Zonal wind tendency due to GWD
6      REAL, allocatable, save :: dv_nonoro_gwd(:, :) ! Meridional wind tendency due to GWD
7      REAL, allocatable, save :: east_gwstress(:, :) ! Eastward stress profile
8      REAL, allocatable, save :: west_gwstress(:, :) ! Westward stress profile
9!$OMP THREADPRIVATE(du_nonoro_gwd,dv_nonoro_gwd,east_gwstress,west_gwstress)
10CONTAINS
11
12      SUBROUTINE nonoro_gwd_ran(ngrid, nlayer, dtime, pp, &
13                   pn2, presnivs, pt, pu, pv, &
14                   zustr, zvstr, d_t, d_u, d_v)
15
16!--------------------------------------------------------------------------------
17! Parametrization of the momentum flux deposition due to a discrete
18! number of gravity waves.
19! F. Lott
20! Version 14, Gaussian distribution of the source
21! LMDz model online version     
22! ADAPTED FOR VENUS /  F. LOTT + S. LEBONNOIS
23! Version adapted on 03/04/2013:
24!      - input flux compensated in the deepest layers
25!                           
26! ADAPTED FOR MARS     G.GILLI     02/2016
27!        Revision with F.Forget    06/2016  Variable EP-flux
28!        according to
29!                                           PBL variation (max
30!                                           velocity thermals)
31!                      D.BARDET    01/2020  - reproductibility of
32!                                           the launching altitude
33!                                           calculation
34!                                           - wave characteristic
35!                                           calculation using MOD
36!                                           - adding east_gwstress
37!                                           and west_gwstress variables   
38!
39! ADAPTED FOR GENERIC  D.BARDET    01/2020
40! READAPTED FOR VENUS  S.LEBONNOIS 08/2021
41!---------------------------------------------------------------------------------
42
43     USE ioipsl_getin_p_mod, ONLY : getin_p
44     use geometry_mod, only: cell_area
45!#ifdef CPP_XIOS
46!     use xios_output_mod, only: send_xios_field
47!#endif
48 implicit none
49#include "YOMCST.h"
50
51! 0. DECLARATIONS:
52
53     ! 0.1 INPUTS
54
55     INTEGER, intent(in):: ngrid, nlayer
56     REAL, intent(in):: dtime              ! Time step of the Physics
57     REAL, intent(in):: pp(ngrid, nlayer)  ! Pressure at full levels
58     REAL, intent(in):: pn2(ngrid, nlayer) ! N2 (BV^2) at 1/2 levels
59     REAL, intent(in):: presnivs(nlayer)   ! Approximate pressure for reproductible launching altitude
60     REAL, intent(in):: pt(ngrid, nlayer)  ! Temperature at full levels
61     REAL, intent(in):: pu(ngrid, nlayer)  ! Zonal wind at full levels
62     REAL, intent(in):: pv(ngrid, nlayer)  ! Meridional wind at full levels
63     
64     ! 0.2 OUTPUTS
65
66     REAL, intent(out):: zustr(ngrid)         ! Zonal surface stress
67     REAL, intent(out):: zvstr(ngrid)         ! Meridional surface stress
68     REAL, intent(out):: d_t(ngrid, nlayer)   ! Tendency on temperature
69     REAL, intent(out):: d_u(ngrid, nlayer)   ! Tendency on zonal wind
70     REAL, intent(out):: d_v(ngrid, nlayer)   ! Tendency on meridional wind
71
72     ! 0.3 INTERNAL ARRAYS
73
74     REAL :: tt(ngrid, nlayer) ! Temperature at full levels
75     REAL :: uu(ngrid, nlayer) ! Zonal wind at full levels
76     REAL :: vv(ngrid, nlayer) ! Meridional wind at full levels
77
78     INTEGER ii, jj, ll
79
80     ! 0.3.0 Time scale of the like cycle of the waves parametrized
81     REAL, parameter:: deltat = 24. * 3600.
82
83     ! 0.3.1 Gravity waves specifications
84     INTEGER, parameter:: nk = 2 ! number of horizontal wavenumbers
85     INTEGER, parameter:: np = 2 ! directions (eastward and westward) phase speed
86     INTEGER, parameter:: no = 2 ! absolute values of phase speed
87     INTEGER, parameter:: na = 5 ! Number of realizations to get the phase speed
88     INTEGER, parameter:: nw = nk * np *no ! Total numbers of gravity waves
89     INTEGER jk, jp, jo, jw
90     REAL kstar                  ! Control value to constrain the min horizontal wavenumber by the horizontal resolution
91
92!!! TEST Kobs  and BestFit Diogo
93     REAL, parameter:: kmax = 1.e-4 ! Max horizontal wavenumber (lambda min 50 km)
94! generic : kmax=N/u, u(=30) zonal wind velocity at launch
95     REAL, parameter:: kmin = 1.e-5 ! Min horizontal wavenumber (lambda max 628 km)
96! generic : kmin=1/sqrt(DxDy) Dx and Dy horizontal grid
97
98!----------------------------------------
99! VCD 1.1 tuning
100!     REAL, parameter:: cmax = 61.   ! Max horizontal absolute phase velocity
101!----------------------------------------
102! VCD 2.0 tuning
103      REAL, parameter:: cmax = 61.   ! Max horizontal absolute phase velocity
104!----------------------------------------
105   
106! generic : cmax=zonal wind at launch
107     REAL, parameter:: cmin = 1.    ! Min horizontal absolute phase velocity
108
109!!! Nominal as Gilli2017
110     REAL, parameter:: epflux_max = 0.005  ! Max EP flux value at launching altitude (previous name: RUWMAX)
111
112     REAL cpha                      ! absolute phase velocity frequency
113     REAL zk(nw, ngrid)             ! horizontal wavenumber amplitude
114     REAL zp(nw, ngrid)             ! horizontal wavenumber angle
115     REAL zo(nw, ngrid)             ! absolute frequency
116
117     REAL intr_freq_m(nw, ngrid)          ! Waves Intr. freq. at the 1/2 lev below the full level (previous name: ZOM)
118     REAL intr_freq_p(nw, ngrid)          ! Waves Intr. freq. at the 1/2 lev above the full level (previous name: ZOP)
119     REAL wwm(nw, ngrid)                  ! Wave EP-fluxes at the 1/2 level below the full level
120     REAL wwp(nw, ngrid)                  ! Wave EP-fluxes at the 1/2 level above the full level
121     REAL u_epflux_p(nw, ngrid)           ! Partial zonal flux (=for each wave) at the 1/2 level above the full level (previous name: RUWP)
122     REAL v_epflux_p(nw, ngrid)           ! Partial meridional flux (=for each wave) at the 1/2 level above the full level (previous name: RVWP)
123     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) 
124     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)
125     REAL epflux_0(nw, ngrid)             ! Fluxes at launching level (previous name: RUW0)
126     INTEGER launch                       ! Launching altitude
127     REAL, parameter:: xlaunch = 5e-3     ! Value for top of cloud convective region
128
129     ! 0.3.2 Parameters of waves dissipations
130!----------------------------------------
131! VCD 1.1 tuning
132!    REAL, parameter:: sat   = 0.85   ! saturation parameter
133!    REAL, parameter:: rdiss = 0.1    ! coefficient of dissipation
134!----------------------------------------
135! VCD 2.0 tuning
136     REAL, parameter:: sat   = 0.6    ! saturation parameter
137     REAL, parameter:: rdiss = 8.e-4  ! coefficient of dissipation
138!----------------------------------------
139     REAL, parameter:: zoisec = 1.e-8 ! security for intrinsic freguency
140
141     ! 0.3.3 Background flow at 1/2 levels and vertical coordinate
142     REAL H0bis(ngrid, nlayer)       ! characteristic height of the atmosphere
143     REAL min_k(ngrid)               ! min(kstar,kmin)
144     REAL, save::  H0                ! characteristic height of the atmosphere
145     REAL, save::  MDR               ! characteristic mass density 
146!----------------------------------------
147! VCD 1.1 tuning
148!    REAL, parameter:: pr = 5e5      ! Reference pressure [Pa]   ! VENUS: cloud layer
149!----------------------------------------
150! VCD 2.0 tuning
151     REAL, parameter:: pr = 5e4      ! Reference pressure [Pa]   ! VENUS: cloud layer
152!----------------------------------------
153     REAL, parameter:: tr = 300.     ! Reference temperature [K] ! VENUS: cloud layer
154     REAL zh(ngrid, nlayer + 1)      ! Log-pressure altitude (constant H0)
155     REAL zhbis(ngrid, nlayer + 1)   ! Log-pressure altitude (varying H)
156     REAL uh(ngrid, nlayer + 1)      ! Zonal wind at 1/2 levels
157     REAL vh(ngrid, nlayer + 1)      ! Meridional wind at 1/2 levels
158     REAL ph(ngrid, nlayer + 1)      ! Pressure at 1/2 levels
159     REAL, parameter:: psec = 1.e-9  ! Security to avoid division by 0 pressure
160     REAL bv(ngrid, nlayer + 1)      ! Brunt Vaisala freq. at 1/2 levels
161     REAL, parameter:: bvsec = 1.e-5 ! Security to avoid negative BV
162     REAL href(nlayer + 1)           ! Reference altitude for launching altitude
163
164     REAL ran_num_1, ran_num_2, ran_num_3
165   
166     LOGICAL, save :: firstcall = .true.
167
168!-------------------------------------------------------------------------------------------------- 
169! 1. INITIALISATION
170!------------------
171     IF (firstcall) THEN
172        write(*,*) "nonoro_gwd_ran: FLott non-oro GW scheme is active!"
173        ! Characteristic vertical scale height
174        H0 = RD * tr / RG
175        ! Characteristic mass density at launch altitude
176        MDR = pr / ( RD * tr )
177        ! Control
178        if (deltat .LT. dtime) THEN
179             call abort_physic("nonoro_gwd_ran","gwd random: deltat lower than dtime!",1)
180        endif
181        if (nlayer .LT. nw) THEN
182             call abort_physic("nonoro_gwd_ran","gwd random: nlayer lower than nw!",1)
183        endif
184        firstcall = .false.
185     ENDIF
186
187! For VENUS, we use *_seri variables, that already integrate the different previous tendencies
188     tt(:,:) = pt(:,:)
189     uu(:,:) = pu(:,:)
190     vv(:,:) = pv(:,:)
191
192!    print*,'epflux_max just after firstcall:', epflux_max
193
194! 2. EVALUATION OF THE BACKGROUND FLOW AT SEMI-LEVELS
195!----------------------------------------------------
196     ! Pressure and Inv of pressure, temperature at 1/2 level
197     DO ll = 2, nlayer
198        ph(:, ll) = exp((log(pp(:, ll)) + log(pp(:, ll - 1))) / 2.)
199     end DO
200
201     ph(:, nlayer + 1) = 0.
202     ph(:, 1) = 2. * pp(:, 1) - ph(:, 2)
203   
204     ! Launching altitude for reproductible case
205     DO ll = 2, nlayer
206        href(ll) = exp((log(presnivs(ll)) + log(presnivs(ll - 1))) / 2.)
207     end DO
208     href(nlayer + 1) = 0.
209     href(1) = 2. * presnivs(1) - href(2)
210
211     launch = 0.
212     DO ll =1, nlayer
213        IF (href (ll) / href(1) > xlaunch) launch = ll
214     end DO
215   
216     ! Log-pressure vertical coordinate
217     DO ll = 1, nlayer + 1
218        zh(:,ll) = H0 * log(pr / (ph(:,ll) + psec))
219     end DO
220     
221     ! Winds and Brunt Vaisala frequency
222     DO ll = 2, nlayer
223        uh(:, ll) = 0.5 * (uu(:, ll) + uu(:, ll - 1)) ! Zonal wind
224        vh(:, ll) = 0.5 * (vv(:, ll) + vv(:, ll - 1)) ! Meridional wind
225! VENUS ATTENTION: CP VARIABLE PSTAB CALCULE EN AMONT DES PARAMETRISATIONS
226        bv(:, ll) = max(bvsec,sqrt(pn2(:,ll)))
227     end DO
228
229     bv(:, 1) = bv(:, 2)
230     uh(:, 1) = 0.
231     vh(:, 1) = 0.
232     bv(:, nlayer + 1) = bv(:, nlayer)
233     uh(:, nlayer + 1) = uu(:, nlayer)
234     vh(:, nlayer + 1) = vv(:, nlayer)
235
236     ! TN+GG April/June 2020
237     ! "Individual waves are not supposed
238     ! to occupy the entire domain, but
239     ! only a faction of it" Lott+2012
240     ! minimum value of k between kmin and kstar
241     DO ii = 1, ngrid
242        kstar = RPI / sqrt(cell_area(ii))
243        min_k(ii) = max(kmin,kstar)
244     end DO
245
246! 3. WAVES CHARACTERISTICS CHOSEN RANDOMLY
247!-----------------------------------------
248! The mod function of here a weird arguments
249! are used to produce the waves characteristics
250! in a stochastic way
251     DO jw = 1, nw
252        DO ii = 1, ngrid
253           
254           ran_num_1 = mod(tt(ii, jw) * 10., 1.)
255           ran_num_2 = mod(tt(ii, jw) * 100., 1.)
256
257!! OPTIONS GENERIC DIFF VENUS !!
258           ! angle (random)
259           zp(jw, ii) = (sign(1., 0.5 - ran_num_1) &
260                        + 1.) * RPI / 2.
261           ! angle (0 or RPI)
262           !zp(jw, ii) = RPI*mod(jw,2)
263
264           ! horizontal wavenumber amplitude
265!           zk(jw, ii) = kmin + (kmax - kmin) * ran_num_2
266           zk(jw, ii) = min_k(ii) + (kmax - min_k(ii)) * ran_num_2
267
268           ! horizontal phase speed
269! this computation allows to get a gaussian distribution centered on 0 (right ?)
270! then cmin is not useful, and it favors small values.
271! VCD 2.0 tuning
272           cpha = 0.
273           DO jj = 1, na
274              ran_num_3 = mod(tt(ii, jw + 3 * jj)**2, 1.)
275              cpha = cpha + 2. * cmax *            &
276                     (ran_num_3 - 0.5) *           &
277                     sqrt(3.) / sqrt(na * 1.)
278           end DO
279           IF (cpha < 0.) THEN
280              cpha = - 1. * cpha
281              zp (jw, ii) = zp(jw, ii) + RPI
282           ENDIF
283! otherwise, with the computation below, we get a uniform distribution between cmin and cmax.
284!           ran_num_3 = mod(tt(ii, jw)**2, 1.)
285!           cpha = cmin + (cmax - cmin) * ran_num_3
286
287           ! Intrinsic frequency
288           zo(jw, ii) = cpha * zk(jw, ii)
289           ! Intrinsic frequency is imposed
290           zo(jw, ii) = zo(jw, ii)                                    &
291                      + zk(jw, ii) * cos(zp(jw, ii)) * uh(ii, launch) &
292                      + zk(jw, ii) * sin(zp(jw, ii)) * vh(ii, launch)
293
294           ! Momentum flux at launch level
295           epflux_0(jw, ii) = epflux_max                              &
296                        * mod(100. * (uu(ii, jw)**2 + vv(ii, jw)**2), 1.)
297
298        end DO
299     end DO
300!     print*,'epflux_0 just after waves charac. ramdon:', epflux_0
301
302! 4. COMPUTE THE FLUXES
303!----------------------
304     ! 4.1 Vertical velocity at launching altitude to ensure
305     !     the correct value to the imposed fluxes.
306     !------------------------------------------------------
307     DO jw = 1, nw
308        ! Evaluate intrinsic frequency at launching altutide:
309        intr_freq_p(jw, :) = zo(jw, :)                                &
310                   - zk(jw, :) * cos(zp(jw, :)) * uh(:, launch)       &
311                   - zk(jw, :) * sin(zp(jw, :)) * vh(:, launch)
312     end DO
313 
314       ! Vertical velocity at launch level, value to ensure the imposed
315       ! mom flux:
316        DO jw = 1, nw
317       ! WW is directly a flux, here, not vertical velocity anymore
318            wwp(jw, :) = epflux_0(jw,:)
319            u_epflux_p(jw, :) = cos(zp(jw, :)) * sign(1., intr_freq_p(jw, :)) * epflux_0(jw, :)
320            v_epflux_p(jw, :) = sin(zp(jw, :)) * sign(1., intr_freq_p(jw, :)) * epflux_0(jw, :)
321        end DO
322     
323     ! 4.2 Initial flux at launching altitude
324     !---------------------------------------
325     u_epflux_tot(:, launch) = 0.
326     v_epflux_tot(:, launch) = 0.
327     DO jw = 1, nw
328        u_epflux_tot(:, launch) = u_epflux_tot(:, launch) + u_epflux_p(jw, :)
329        v_epflux_tot(:, launch) = v_epflux_tot(:, launch) + v_epflux_p(jw, :)
330     end DO
331
332     ! 4.3 Loop over altitudes, with passage from one level to the next done by:
333     !--------------------------------------------------------------------------
334     !     i) conserving the EP flux,
335     !     ii) dissipating a little,
336     !     iii) testing critical levels,
337     !     iv) testing the breaking.
338     !--------------------------------------------------------------------------
339     DO ll = launch, nlayer - 1
340        ! Warning! all the physics is here (passage from one level to the next
341        DO jw = 1, nw
342           intr_freq_m(jw, :) = intr_freq_p(jw, :)
343           wwm(jw, :) = wwp(jw, :)
344           ! Intrinsic frequency
345           intr_freq_p(jw, :) = zo(jw, :) - zk(jw, :) * cos(zp(jw, :)) * uh(:, ll + 1)    &
346                                          - zk(jw, :) * sin(zp(jw, :)) * vh(:, ll + 1)
347           ! Vertical velocity in flux formulation
348           wwp(jw, :) = min(                                                              &
349                      ! No breaking (Eq. 6):
350                      wwm(jw, :)                                                          &
351                      ! Dissipation (Eq. 8):
352                      * exp(-rdiss * pr / (ph(:, ll + 1) + ph (:, ll))                    &
353                      * ((bv(:, ll + 1) + bv (:, ll)) / 2.)**3                            &
354                      / max(abs(intr_freq_p(jw, :) + intr_freq_m(jw, :)) / 2., zoisec)**4 &
355                      * zk(jw, :)**3 * (zh(:, ll + 1) - zh(:, ll))) ,                     &
356                      ! Critical levels (forced to zero if intrinsic frequency
357                      ! changes sign)
358                      max(0., sign(1., intr_freq_p(jw, :) * intr_freq_m(jw, :)))          &
359                      ! Saturation (Eq. 12)
360                      * abs(intr_freq_p(jw, :))**3 / bv(:, ll + 1)                        &
361                      * exp(-zh(:, ll + 1) / H0)                                          &
362!! Correction for VCD 2.0
363!                      * sat**2 * kmin**2 / zk(jw, :)**4)
364                      * sat**2 * min_k(:)**2 * MDR / zk(jw, :)**4)                       
365        end DO
366
367     ! Evaluate EP-flux from Eq. 7 and give the right orientation to the stress
368        DO jw = 1, nw
369           u_epflux_p(jw, :) = sign(1., intr_freq_p(jw, :)) * cos(zp(jw, :)) * wwp(jw, :)
370           v_epflux_p(jw, :) = sign(1., intr_freq_p(jw, :)) * sin(zp(jw, :)) * wwp(jw, :)
371        end DO
372
373        u_epflux_tot(:, ll + 1) = 0.
374        v_epflux_tot(:, ll + 1) = 0.
375        DO jw = 1, nw
376           u_epflux_tot(:, ll + 1) = u_epflux_tot(:, ll + 1) + u_epflux_p(jw, :)
377           v_epflux_tot(:, ll + 1) = v_epflux_tot(:, ll + 1) + v_epflux_p(jw, :)
378           east_gwstress(:, ll) = east_gwstress(:, ll)                         & 
379                                + max(0., u_epflux_p(jw, :)) / float(nw)
380           west_gwstress(:, ll) = west_gwstress(:, ll)                         &
381                                + min(0., u_epflux_p(jw, ll)) / float(nw)
382        end DO
383     end DO ! DO LL = LAUNCH, nlayer - 1
384
385!    print*,'u_epflux_tot just after launching:', u_epflux_tot
386!    print*,'v_epflux_tot just after launching:', v_epflux_tot
387!    print*,'u_epflux_p just after launching:', maxval(u_epflux_p), minval(u_epflux_p)
388!    print*,'v_epflux_p just after launching:', maxval(v_epflux_p), minval(v_epflux_p)
389
390! 5. TENDENCY CALCULATIONS
391!-------------------------
392     ! 5.1 Flow rectification at the top and in the low layers
393     ! --------------------------------------------------------
394     ! Warning, this is the total on all GW...
395
396     u_epflux_tot(:, nlayer + 1) = 0.
397     v_epflux_tot(:, nlayer + 1) = 0.
398
399     ! Here, big change compared to FLott version:
400     ! We compensate (u_epflux_(:, LAUNCH), ie total emitted upward flux
401     ! over the layers max(1,LAUNCH-3) to LAUNCH-1
402     DO LL = 1, max(1,LAUNCH-3)
403        u_epflux_tot(:, LL) = 0.
404        v_epflux_tot(:, LL) = 0.
405     end DO
406     DO LL = max(2,LAUNCH-2), LAUNCH-1
407        u_epflux_tot(:, LL) = u_epflux_tot(:, LL - 1)                          &
408                            + u_epflux_tot(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1))  &
409                            / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
410        v_epflux_tot(:, LL) = v_epflux_tot(:, LL - 1)                          &
411                            + v_epflux_tot(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1))  &
412                            / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
413        east_gwstress(:,LL) = east_gwstress(:, LL - 1)                         &
414                            + east_gwstress(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1)) &
415                            / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
416        west_gwstress(:,LL) = west_gwstress(:, LL - 1)                         &
417                            + west_gwstress(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1)) &
418                            / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
419     end DO
420     ! This way, the total flux from GW is zero, but there is a net transport
421     ! (upward) that should be compensated by circulation
422     ! and induce additional friction at the surface
423
424     ! 5.2 AR-1 RECURSIVE FORMULA (13) IN VERSION 4
425     !---------------------------------------------
426     DO LL = 1, nlayer
427!       d_u(:, LL) = (u_epflux_tot(:, LL + 1) - u_epflux_tot(:, LL))
428       d_u(:, LL) = RG * (u_epflux_tot(:, LL + 1) - u_epflux_tot(:, LL)) &
429                  / (PH(:, LL + 1) - PH(:, LL)) * DTIME
430!       d_v(:, LL) = (v_epflux_tot(:, LL + 1) - v_epflux_tot(:, LL))
431       d_v(:, LL) = RG * (v_epflux_tot(:, LL + 1) - v_epflux_tot(:, LL)) &
432                  / (PH(:, LL + 1) - PH(:, LL)) * DTIME
433     ENDDO
434     d_t(:,:) = 0.
435
436     ! 5.3 Update tendency of wind with the previous (and saved) values
437     !-----------------------------------------------------------------
438     d_u(:,:) = DTIME/DELTAT/REAL(NW) * d_u(:,:)                       &
439              + (1.-DTIME/DELTAT) * du_nonoro_gwd(:,:)
440     d_v(:,:) = DTIME/DELTAT/REAL(NW) * d_v(:,:)                       &
441              + (1.-DTIME/DELTAT) * dv_nonoro_gwd(:,:)
442     du_nonoro_gwd(:,:) = d_u(:,:)
443     dv_nonoro_gwd(:,:) = d_v(:,:)
444
445!    print*,'u_epflux_tot just after tendency:', u_epflux_tot
446!    print*,'v_epflux_tot just after tendency:', v_epflux_tot
447!    print*,'d_u just after tendency:', maxval(d_u(:,:)), minval(d_u)
448!    print*,'d_v just after tendency:', maxval(d_v(:,:)), minval(d_v)
449
450     ! Cosmetic: evaluation of the cumulated stress
451     ZUSTR(:) = 0.
452     ZVSTR(:) = 0.
453     DO LL = 1, nlayer
454        ZUSTR(:) = ZUSTR(:) + D_U(:, LL) / RG * (PH(:, LL + 1) - PH(:, LL))
455        ZVSTR(:) = ZVSTR(:) + D_V(:, LL) / RG * (PH(:, LL + 1) - PH(:, LL))
456     ENDDO
457
458!#ifdef CPP_XIOS
459!     call send_xios_field("du_nonoro", d_u)
460!     call send_xios_field("dv_nonoro", d_v)
461!     call send_xios_field("east_gwstress", east_gwstress)
462!     call send_xios_field("west_gwstress", west_gwstress)
463!#endif
464
465      END SUBROUTINE nonoro_gwd_ran
466
467! ===================================================================
468! Subroutines used to write variables of memory in start files       
469! ===================================================================
470
471      SUBROUTINE ini_nonoro_gwd_ran(ngrid,nlayer)
472
473      IMPLICIT NONE
474
475      INTEGER, INTENT (in) :: ngrid  ! number of atmospheric columns
476      INTEGER, INTENT (in) :: nlayer ! number of atmospheric layers
477
478         allocate(du_nonoro_gwd(ngrid, nlayer))
479         allocate(dv_nonoro_gwd(ngrid, nlayer))
480         allocate(east_gwstress(ngrid, nlayer))
481         allocate(west_gwstress(ngrid, nlayer))
482
483      END SUBROUTINE ini_nonoro_gwd_ran
484! ----------------------------------
485      SUBROUTINE end_nonoro_gwd_ran
486
487      IMPLICIT NONE
488
489         if (allocated(du_nonoro_gwd)) deallocate(du_nonoro_gwd)
490         if (allocated(dv_nonoro_gwd)) deallocate(dv_nonoro_gwd)
491         if (allocated(east_gwstress)) deallocate(east_gwstress)
492         if (allocated(west_gwstress)) deallocate(west_gwstress)
493
494      END SUBROUTINE end_nonoro_gwd_ran
495
496END MODULE nonoro_gwd_ran_mod
497
498
Note: See TracBrowser for help on using the repository browser.