source: trunk/LMDZ.GENERIC/libf/phystd/nonoro_gwd_ran_mod.F90 @ 2590

Last change on this file since 2590 was 2550, checked in by emillour, 4 years ago

Generic GCM:
Some OpenMP fixes in routines initracer.F, nonoro_gwd_ran_mod.F90,
phys_state_var_mod.F90 and sugas_corrk.F90
EM

File size: 23.0 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                   zmax_therm, pt, pu, pv, pdt, pdu, pdv, &
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!---------------------------------------------------------------------------------
41
42
43
44     use comcstfi_mod, only: g, pi, cpp, r
45     USE ioipsl_getin_p_mod, ONLY : getin_p
46     use assert_m, only : assert
47     use vertical_layers_mod, only : presnivs
48     use callkeys_mod, only : epflux_max, & ! Max EP flux value at launching altitude (previous name: RUWMAX)
49                              gwd_convective_source
50     use inifis_mod
51     use geometry_mod, only: cell_area
52#ifdef CPP_XIOS
53     use xios_output_mod, only: send_xios_field
54#endif
55 implicit none
56
57! 0. DECLARATIONS:
58
59     ! 0.1 INPUTS
60
61     INTEGER, intent(in):: ngrid, nlayer
62     REAL, intent(in):: dtime              ! Time step of the Physics
63     REAL, intent(in):: zmax_therm(ngrid)  ! Altitude of max velocity thermals (m)
64     REAL, intent(in):: pp(ngrid, nlayer)  ! Pressure at full levels
65     REAL, intent(in):: pt(ngrid, nlayer)  ! Temperature at full levels
66     REAL, intent(in):: pu(ngrid, nlayer)  ! Zonal wind at full levels
67     REAL, intent(in):: pv(ngrid, nlayer)  ! Meridional wind at full levels
68     REAL, intent(in):: pdt(ngrid, nlayer) ! Tendency on temperature
69     REAL, intent(in):: pdu(ngrid, nlayer) ! Tendency on zonal wind
70     REAL, intent(in):: pdv(ngrid, nlayer) ! Tendency on meridional wind
71     
72     ! 0.2 OUTPUTS
73
74     REAL, intent(out):: zustr(ngrid)         ! Zonal surface stress
75     REAL, intent(out):: zvstr(ngrid)         ! Meridional surface stress
76     REAL, intent(out):: d_t(ngrid, nlayer)   ! Tendency on temperature
77     REAL, intent(out):: d_u(ngrid, nlayer)   ! Tendency on zonal wind
78     REAL, intent(out):: d_v(ngrid, nlayer)   ! Tendency on meridional wind
79
80     ! 0.3 INTERNAL ARRAYS
81
82     REAL :: tt(ngrid, nlayer) ! Temperature at full levels
83     REAL :: uu(ngrid, nlayer) ! Zonal wind at full levels
84     REAL :: vv(ngrid, nlayer) ! Meridional wind at full levels
85     REAL :: bvlow(ngrid)      ! Qu est-ce ?
86     REAL :: dz
87
88     INTEGER ii, jj, ll
89
90     ! 0.3.0 Time scale of the like cycle of the waves parametrized
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
100     REAL kstar                  ! Control value to constrain the min horizontal wavenumber by the horizontal resolution
101     REAL, parameter:: kmax = 1.e-4 ! Max horizontal wavenumber=N/u, u(=30) zonal wind velocity at launch
102     !REAL, parameter:: kmax = 4.e-5 ! Max horizontal wavenumber=N/u, u(=70) zonal wind velocity at launch
103     REAL, parameter:: kmin = 2.e-6 ! Min horizontal wavenumber=1/sqrt(DxDy) Dx and Dy horizontal grid
104     REAL, parameter:: cmax = 30.   ! Max horizontal absolute phase velocity=zonal wind at launch
105     !REAL, parameter:: cmax = 100.  ! Test for Saturn: Max horizontal absolute phase velocity
106     !REAL, parameter:: cmax = 70.   ! Test for Saturn: Max horizontal absolute phase velocity
107     REAL, parameter:: cmin = 1.    ! Min horizontal absolute phase velocity
108     REAL cpha                      ! absolute phase velocity frequency
109     REAL zk(nw, ngrid)             ! horizontal wavenumber amplitude
110     REAL zp(nw, ngrid)             ! horizontal wavenumber angle
111     REAL zo(nw, ngrid)             ! absolute frequency
112
113     REAL intr_freq_m(nw, ngrid)          ! Waves Intr. freq. at the 1/2 lev below the full level (previous name: ZOM)
114     REAL intr_freq_p(nw, ngrid)          ! Waves Intr. freq. at the 1/2 lev above the full level (previous name: ZOP)
115     REAL wwm(nw, ngrid)                  ! Wave EP-fluxes at the 1/2 level below the full level
116     REAL wwp(nw, ngrid)                  ! Wave EP-fluxes at the 1/2 level above the full level
117     REAL u_epflux_p(nw, ngrid)           ! Partial zonal flux (=for each wave) at the 1/2 level above the full level (previous name: RUWP)
118     REAL v_epflux_p(nw, ngrid)           ! Partial meridional flux (=for each wave) at the 1/2 level above the full level (previous name: RVWP)
119     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) 
120     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)
121     REAL epflux_0(nw, ngrid)             ! Fluxes at launching level (previous name: RUW0)
122     INTEGER launch                       ! Launching altitude
123     REAL, parameter:: xlaunch = 0.2      ! Control the launching altitude
124     REAL, parameter:: zmaxth_top = 8000. ! Top of convective layer in m (approx.)
125
126     REAL prec(ngrid)     ! precipitations
127     REAL prmax           ! Max value of precipitation
128
129     ! 0.3.2 Parameters of waves dissipations
130     REAL, parameter:: sat   = 1.     ! saturation parameter
131     REAL, parameter:: rdiss = 1.     ! coefficient of dissipation
132     REAL, parameter:: zoisec = 1.e-6 ! security for intrinsic freguency
133
134     ! 0.3.3 Background flow at 1/2 levels and vertical coordinate
135     REAL H0bis(ngrid, nlayer)       ! characteristic height of the atmosphere
136     REAL, save::  H0                ! characteristic height of the atmosphere
137     !REAL, parameter:: pr = 250      ! Reference pressure [Pa]
138     !REAL, parameter:: tr = 220.     ! Reference temperature [K]
139     ! TEST FOR SATURN STRATOSPHERE
140     REAL, parameter:: pr = 110      ! Reference pressure [Pa]
141     REAL, parameter:: tr = 130.     ! Reference temperature [K]
142     REAL zh(ngrid, nlayer + 1)      ! Log-pressure altitude (constant H0)
143     REAL zhbis(ngrid, nlayer + 1)   ! Log-pressure altitude (varying H)
144     REAL uh(ngrid, nlayer + 1)      ! Zonal wind at 1/2 levels
145     REAL vh(ngrid, nlayer + 1)      ! Meridional wind at 1/2 levels
146     REAL ph(ngrid, nlayer + 1)      ! Pressure at 1/2 levels
147     REAL, parameter:: psec = 1.e-6  ! Security to avoid division by 0 pressure
148     REAL bv(ngrid, nlayer + 1)      ! Brunt Vaisala freq. at 1/2 levels
149     REAL, parameter:: bvsec = 1.e-5 ! Security to avoid negative BV
150     REAL href(nlayer + 1)           ! Reference altitude for launching altitude
151
152
153
154     REAL ran_num_1, ran_num_2, ran_num_3
155   
156
157     LOGICAL, save :: firstcall = .true.
158
159!-------------------------------------------------------------------------------------------------- 
160! 1. INITIALISATION
161!------------------
162     IF (firstcall) THEN
163        write(*,*) "nonoro_gwd_ran: FLott non-oro GW scheme is active!"
164        ! Characteristic vertical scale height
165        H0 = r * tr / g
166        ! Control
167        if (deltat .LT. dtime) THEN
168             call abort_physic("nonoro_gwd_ran","gwd random: deltat lower than dtime!",1)
169        endif
170        if (nlayer .LT. nw) THEN
171             call abort_physic("nonoro_gwd_ran","gwd random: nlayer lower than nw!",1)
172        endif
173        firstcall = .false.
174     ENDIF
175     gwd_convective_source = .false.
176
177     ! Compute subroutine's current values of temperature and winds
178     tt(:,:) = pt(:,:) + dtime * pdt(:,:)
179     uu(:,:) = pu(:,:) + dtime * pdu(:,:)
180     vv(:,:) = pv(:,:) + dtime * pdv(:,:)
181!    print*,'epflux_max just after firstcall:', epflux_max
182
183
184! 2. EVALUATION OF THE BACKGROUND FLOW AT SEMI-LEVELS
185!----------------------------------------------------
186     ! Pressure and Inv of pressure, temperature at 1/2 level
187     DO ll = 2, nlayer
188        ph(:, ll) = exp((log(pp(:, ll)) + log(pp(:, ll - 1))) / 2.)
189     end DO
190
191     ph(:, nlayer + 1) = 0.
192     ph(:, 1) = 2. * pp(:, 1) - ph(:, 2)
193   
194     ! Launching altitude for reproductible case
195     DO ll = 2, nlayer
196        href(ll) = exp((log(presnivs(ll)) + log(presnivs(ll - 1))) / 2.)
197     end DO
198     href(nlayer + 1) = 0.
199     href(1) = 2. * presnivs(1) - href(2)
200
201     launch = 0.
202     DO ll =1, nlayer
203        IF (href (ll) / href(1) > xlaunch) launch = ll
204     end DO
205   
206     ! Log-pressure vertical coordinate
207     DO ll = 1, nlayer + 1
208        zh(:,ll) = H0 * log(pr / (ph(:,ll) + psec))
209     end DO
210     
211     ! Winds and Brunt Vaisala frequency
212     DO ll = 2, nlayer
213        uh(:, ll) = 0.5 * (uu(:, ll) + uu(:, ll - 1)) ! Zonal wind
214        vh(:, ll) = 0.5 * (vv(:, ll) + vv(:, ll - 1)) ! Meridional wind
215
216        bv(:, ll) = 0.5 * (tt(:, ll) + tt(:, ll - 1)) &
217                    * r**2 / cpp / H0**2              &
218                    + (tt(:, ll) - tt(:, ll - 1))     &
219                    / (zh(:, ll) - zh(:, ll - 1))     &
220                    * r / H0
221        bv(:, ll) = sqrt(max(bvsec, bv(:,ll)))
222     end DO
223
224     bv(:, 1) = bv(:, 2)
225     uh(:, 1) = 0.
226     vh(:, 1) = 0.
227     bv(:, nlayer + 1) = bv(:, nlayer)
228     uh(:, nlayer + 1) = uu(:, nlayer)
229     vh(:, nlayer + 1) = vv(:, nlayer)
230
231! 3. WAVES CHARACTERISTICS CHOSEN RANDOMLY
232!-----------------------------------------
233! The mod function of here a weird arguments
234! are used to produce the waves characteristics
235! in a stochastic way
236     DO jw = 1, nw
237        DO ii = 1, ngrid
238           
239           ran_num_1 = mod(tt(ii, jw) * 10., 1.)
240           ran_num_2 = mod(tt(ii, jw) * 100., 1.)
241
242           ! angle (0 or pi so far)
243           zp(jw, ii) = (sign(1., 0.5 - ran_num_1) &
244                        + 1.) * pi / 2.
245           ! horizontal wavenumber amplitude
246           ! TN+GG April/June 2020
247           ! "Individual waves are not supposed
248           ! to occupy the entire domain, but
249           ! only a faction of it" Lott+2012
250
251           !zk(jw, ii) = kmin + (kmax - kmin) * ran_num_2
252           kstar = pi / sqrt(cell_area(ii))
253           zk(jw, ii) = max(kmin,kstar) + (kmax - max(kmin,kstar)) * ran_num_2
254           ! horizontal phase speed
255           cpha = 0.
256           DO jj = 1, na
257              ran_num_3 = mod(tt(ii, jw + 3 * jj)**2, 1.)
258              cpha = cpha + 2. * cmax *            &
259                     (ran_num_3 - 0.5) *           &
260                     sqrt(3.) / sqrt(na * 1.)
261           end DO
262           IF (cpha < 0.) THEN
263              cpha = - 1. * cpha
264              zp (jw, ii) = zp(jw, ii) + pi
265           ENDIF
266           ! Intrinsic frequency
267           zo(jw, ii) = cpha * zk(jw, ii)
268           ! Intrinsic frequency is imposed
269           zo(jw, ii) = zo(jw, ii)                                    &
270                      + zk(jw, ii) * cos(zp(jw, ii)) * uh(ii, launch) &
271                      + zk(jw, ii) * sin(zp(jw, ii)) * vh(ii, launch)
272           ! Momentum flux at launch level
273           epflux_0(jw, ii) = epflux_max                              &
274                        * mod(100. * (uu(ii, jw)**2 + vv(ii, jw)**2), 1.)
275
276        end DO
277     end DO
278
279!    print*,'epflux_max just after waves charac. ramdon:', epflux_max
280!    print*,'epflux_0 just after waves charac. ramdon:', epflux_0
281! 4. COMPUTE THE FLUXES
282!----------------------
283     ! 4.1 Vertical velocity at launching altitude to ensure
284     !     the correct value to the imposed fluxes.
285     !------------------------------------------------------
286     DO jw = 1, nw
287        ! Evaluate intrinsic frequency at launching altutide:
288        intr_freq_p(jw, :) = zo(jw, :)                                &
289                   - zk(jw, :) * cos(zp(jw, :)) * uh(:, launch)       &
290                   - zk(jw, :) * sin(zp(jw, :)) * vh(:, launch)
291     end DO
292 
293     IF (gwd_convective_source) THEN
294        DO jw = 1, nw
295       ! VERSION WITH CONVECTIVE SOURCE ! designed for Earth
296
297       ! Vertical velocity at launch level, value to ensure the
298       ! imposed mmt flux factor related to the convective forcing:
299       ! precipitations.
300
301       ! tanh limitation to values above prmax:
302!       WWP(JW, :) = epflux_0(JW, :) &
303!            * (r / cpp / H0 * RLVTT * PRMAX * TANH(PREC(:) / PRMAX))**2
304!       Here, we neglected the kinetic energy providing of the thermodynamic
305!       phase change
306
307!
308
309       ! Factor related to the characteristics of the waves:
310            wwp(jw, :) = wwp(jw, :) * zk(jw, :)**3 / kmin / bvlow(:)  &
311                       / MAX(ABS(intr_freq_p(jw, :)), zoisec)**3
312
313      ! Moderation by the depth of the source (dz here):
314            wwp(jw, :) = wwp(jw, :)                                          &
315                       * exp(- bvlow(:)**2 / max(abs(intr_freq_p(jw, :)), zoisec)**2 &
316                       * zk(jw, :)**2 * dz**2)
317
318      ! Put the stress in the right direction:
319            u_epflux_p(jw, :) = intr_freq_p(jw, :) / max(abs(intr_freq_p(jw, :)), zoisec)**2       &
320                        * bv(:, launch) * cos(zp(jw, :)) * wwp(jw, :)**2
321            v_epflux_p(JW, :) = intr_freq_p(jw, :) / max(abs(intr_freq_p(jw, :)), zoisec)**2       &
322                        * bv(:, launch) * sin(zp(jw, :)) * wwp(jw, :)**2
323        end DO
324     ELSE ! VERSION WITHOUT CONVECTIVE SOURCE
325       ! Vertical velocity at launch level, value to ensure the imposed
326       ! mom flux:
327        DO jw = 1, nw
328       ! WW is directly a flux, here, not vertical velocity anymore
329            wwp(jw, :) = epflux_0(JW,:)
330            u_epflux_p(jw, :) = cos(zp(jw, :)) * sign(1., intr_freq_p(jw, :)) * epflux_0(jw, :)
331            v_epflux_p(jw, :) = sin(zp(jw, :)) * sign(1., intr_freq_p(jw, :)) * epflux_0(jw, :)
332
333        end DO
334     ENDIF
335     
336     ! 4.2 Initial flux at launching altitude
337     !---------------------------------------
338     u_epflux_tot(:, launch) = 0.
339     v_epflux_tot(:, launch) = 0.
340     DO jw = 1, nw
341        u_epflux_tot(:, launch) = u_epflux_tot(:, launch) + u_epflux_p(jw, :)
342        v_epflux_tot(:, launch) = v_epflux_tot(:, launch) + v_epflux_p(jw, :)
343     end DO
344
345     ! 4.3 Loop over altitudes, with passage from one level to the next done by:
346     !--------------------------------------------------------------------------
347     !     i) conserving the EP flux,
348     !     ii) dissipating a little,
349     !     iii) testing critical levels,
350     !     iv) testing the breaking.
351     !--------------------------------------------------------------------------
352     DO ll = launch, nlayer - 1
353        ! Warning! all the physics is here (passage from one level to the next
354        DO jw = 1, nw
355           intr_freq_m(jw, :) = intr_freq_p(jw, :)
356           wwm(jw, :) = wwp(jw, :)
357           ! Intrinsic frequency
358           intr_freq_p(jw, :) = zo(jw, :) - zk(jw, :) * cos(zp(jw, :)) * uh(:, ll + 1)    &
359                      - zk(jw, :) * sin(zp(jw, :)) * vh(:, ll + 1)
360           ! Vertical velocity in flux formulation
361           wwp(jw, :) = min(                                                              &
362                      ! No breaking (Eq. 6):
363                      wwm(jw, :)                                                          &
364                      ! Dissipation (Eq. 8):
365                      * exp(-rdiss * pr / (ph(:, ll + 1) + ph (:, ll))                    &
366                      * ((bv(:, ll + 1) + bv (:, ll)) / 2.)**3                            &
367                      / max(abs(intr_freq_p(jw, :) + intr_freq_m(jw, :)) / 2., zoisec)**4 &
368                      * zk(jw, :)**3 * (zh(:, ll + 1) - zh(:, ll))) ,                     &
369                      ! Critical levels (forced to zero if intrinsic frequency
370                      ! changes sign)
371                      max(0., sign(1., intr_freq_p(jw, :) * intr_freq_m(jw, :)))          &
372                      ! Saturation (Eq. 12)
373                      * abs(intr_freq_p(jw, :))**3 / bv(:, ll + 1)                        &
374                      * exp(-zh(:, ll + 1) / H0)                                          &
375                      * sat**2 * kmin**2 / zk(jw, :)**4)
376        end DO
377
378     ! Evaluate EP-flux from Eq. 7 and give the right orientation to the stress
379        DO jw = 1, nw
380           u_epflux_p(jw, :) = sign(1., intr_freq_p(jw, :)) * cos(zp(jw, :)) * wwp(jw, :)
381           v_epflux_p(jw, :) = sign(1., intr_freq_p(jw, :)) * sin(zp(jw, :)) * wwp(jw, :)
382        end DO
383
384        u_epflux_tot(:, ll + 1) = 0.
385        v_epflux_tot(:, ll + 1) = 0.
386        DO jw = 1, nw
387           u_epflux_tot(:, ll + 1) = u_epflux_tot(:, ll + 1) + u_epflux_p(jw, :)
388           v_epflux_tot(:, ll + 1) = v_epflux_tot(:, ll + 1) + v_epflux_p(jw, :)
389           east_gwstress(:, ll) = east_gwstress(:, ll)                         & 
390                                + max(0., u_epflux_p(jw, :)) / float(nw)
391           west_gwstress(:, ll) = west_gwstress(:, ll)                         &
392                                + min(0., u_epflux_p(jw, ll)) / float(nw)
393        end DO
394     end DO ! DO LL = LAUNCH, nlayer - 1
395
396!    print*,'u_epflux_tot just after launching:', u_epflux_tot
397!    print*,'v_epflux_tot just after launching:', v_epflux_tot
398!    print*,'u_epflux_p just after launching:', maxval(u_epflux_p), minval(u_epflux_p)
399!    print*,'v_epflux_p just after launching:', maxval(v_epflux_p), minval(v_epflux_p)
400
401! 5. TENDENCY CALCULATIONS
402!-------------------------
403     ! 5.1 Flow rectification at the top and in the low layers
404     ! --------------------------------------------------------
405     ! Warning, this is the total on all GW...
406
407     u_epflux_tot(:, nlayer + 1) = 0.
408     v_epflux_tot(:, nlayer + 1) = 0.
409
410     ! Here, big change compared to FLott version:
411     ! We compensate (u_epflux_(:, LAUNCH), ie total emitted upward flux
412     ! over the layers max(1,LAUNCH-3) to LAUNCH-1
413     DO LL = 1, max(1,LAUNCH-3)
414        u_epflux_tot(:, LL) = 0.
415        v_epflux_tot(:, LL) = 0.
416     end DO
417     DO LL = max(2,LAUNCH-2), LAUNCH-1
418        u_epflux_tot(:, LL) = u_epflux_tot(:, LL - 1)                          &
419                            + u_epflux_tot(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1))  &
420                            / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
421        v_epflux_tot(:, LL) = v_epflux_tot(:, LL - 1)                          &
422                            + v_epflux_tot(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1))  &
423                            / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
424        east_gwstress(:,LL) = east_gwstress(:, LL - 1)                         &
425                            + east_gwstress(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1)) &
426                            / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
427        west_gwstress(:,LL) = west_gwstress(:, LL - 1)                         &
428                            + west_gwstress(:, LAUNCH) * (PH(:,LL)-PH(:,LL-1)) &
429                            / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
430     end DO
431     ! This way, the total flux from GW is zero, but there is a net transport
432     ! (upward) that should be compensated by circulation
433     ! and induce additional friction at the surface
434
435     ! 5.2 AR-1 RECURSIVE FORMULA (13) IN VERSION 4
436     !---------------------------------------------
437     DO LL = 1, nlayer
438       d_u(:, LL) = (u_epflux_tot(:, LL + 1) - u_epflux_tot(:, LL))! &
439       !d_u(:, LL) = G * (u_epflux_tot(:, LL + 1) - u_epflux_tot(:, LL))! &
440        !          / (PH(:, LL + 1) - PH(:, LL)) * DTIME
441       d_v(:, LL) = (v_epflux_tot(:, LL + 1) - v_epflux_tot(:, LL))! &
442       !d_v(:, LL) = G * (v_epflux_tot(:, LL + 1) - v_epflux_tot(:, LL))! &
443         !         / (PH(:, LL + 1) - PH(:, LL)) * DTIME
444     ENDDO
445     d_t(:,:) = 0.
446
447     ! 5.3 Update tendency of wind with the previous (and saved) values
448     !-----------------------------------------------------------------
449     d_u(:,:) = DTIME/DELTAT/REAL(NW) * d_u(:,:)                       &
450              + (1.-DTIME/DELTAT) * du_nonoro_gwd(:,:)
451     d_v(:,:) = DTIME/DELTAT/REAL(NW) * d_v(:,:)                       &
452              + (1.-DTIME/DELTAT) * dv_nonoro_gwd(:,:)
453     du_nonoro_gwd(:,:) = d_u(:,:)
454     dv_nonoro_gwd(:,:) = d_v(:,:)
455
456!    print*,'u_epflux_tot just after tendency:', u_epflux_tot
457!    print*,'v_epflux_tot just after tendency:', v_epflux_tot
458!    print*,'d_u just after tendency:', maxval(d_u(:,:)), minval(d_u)
459!    print*,'d_v just after tendency:', maxval(d_v(:,:)), minval(d_v)
460     ! Cosmetic: evaluation of the cumulated stress
461
462     ZUSTR(:) = 0.
463     ZVSTR(:) = 0.
464     DO LL = 1, nlayer
465        ZUSTR(:) = ZUSTR(:) + D_U(:, LL) !/ g * (PH(:, LL + 1) - PH(:, LL))
466        ZVSTR(:) = ZVSTR(:) + D_V(:, LL) !/ g * (PH(:, LL + 1) - PH(:, LL))
467     ENDDO
468
469#ifdef CPP_XIOS
470     call send_xios_field("du_nonoro", d_u)
471     call send_xios_field("dv_nonoro", d_v)
472     call send_xios_field("east_gwstress", east_gwstress)
473     call send_xios_field("west_gwstress", west_gwstress)
474#endif
475
476      END SUBROUTINE nonoro_gwd_ran
477
478! ===================================================================
479! Subroutines used to write variables of memory in start files       
480! ===================================================================
481
482      SUBROUTINE ini_nonoro_gwd_ran(ngrid,nlayer)
483
484      IMPLICIT NONE
485
486      INTEGER, INTENT (in) :: ngrid  ! number of atmospheric columns
487      INTEGER, INTENT (in) :: nlayer ! number of atmospheric layers
488
489         allocate(du_nonoro_gwd(ngrid, nlayer))
490         allocate(dv_nonoro_gwd(ngrid, nlayer))
491         allocate(east_gwstress(ngrid, nlayer))
492         allocate(west_gwstress(ngrid, nlayer))
493
494      END SUBROUTINE ini_nonoro_gwd_ran
495! ----------------------------------
496      SUBROUTINE end_nonoro_gwd_ran
497
498      IMPLICIT NONE
499
500         if (allocated(du_nonoro_gwd)) deallocate(du_nonoro_gwd)
501         if (allocated(dv_nonoro_gwd)) deallocate(dv_nonoro_gwd)
502         if (allocated(east_gwstress)) deallocate(east_gwstress)
503         if (allocated(west_gwstress)) deallocate(west_gwstress)
504
505      END SUBROUTINE end_nonoro_gwd_ran
506
507END MODULE nonoro_gwd_ran_mod
508
Note: See TracBrowser for help on using the repository browser.