source: LMDZ6/trunk/libf/phylmd/lmdz_gwd_ini.f90

Last change on this file was 6061, checked in by evignon, 5 weeks ago

petite correction pour le 1D

File size: 18.8 KB
Line 
1MODULE lmdz_gwd_ini
2! -----------------------------------------------------
3! This module initializes all the parameters needed
4! for gravity wave drag (gwd) parameterizations
5!------------------------------------------------------
6
7   IMPLICIT NONE
8
9   REAL :: RCPD, RD, RG, RPI, RA, ROMEGA, RKAPPA, RLVTT
10   !$OMP THREADPRIVATE(RCPD, RD, RG, RPI, RA, ROMEGA, RKAPPA, RLVTT)
11
12   INTEGER, SAVE :: lunout, prt_level            ! Logical unit number and level for standard output
13   !$OMP THREADPRIVATE(lunout,prt_level)
14
15   ! for orographic gravity wave effects
16   LOGICAL :: ok_strato
17   !$OMP THREADPRIVATE(ok_strato)
18
19   LOGICAL, SAVE, PROTECTED :: ok_orodr = .true.  ! activate orographic gravity wave drag
20   !$OMP THREADPRIVATE(ok_orodr)
21
22   LOGICAL, SAVE, PROTECTED :: ok_orolf = .true.  ! activate orographic gravity wave lift
23   !$OMP THREADPRIVATE(ok_orolf)
24
25   REAL, SAVE, PROTECTED :: zpmm_orodr_t = 100.  ! threshold on zpeak-zmean [m] to activate sso gwd
26   !$OMP THREADPRIVATE(zpmm_orodr_t)
27
28   REAL, SAVE, PROTECTED :: zstd_orodr_t = 10.    ! threshold on zstd [m] to activate sso gwd
29   !$OMP THREADPRIVATE(zstd_orodr_t)
30
31   REAL, SAVE, PROTECTED :: zpmm_orolf_t = 100.  ! threshold on zpeak-zmean [m] to activate sso lift
32   !$OMP THREADPRIVATE(zpmm_orolf_t)
33
34   REAL, SAVE, PROTECTED :: nm_oro_t = -1      ! number of subgrid-scale mountains above which sso drag and lift are active
35   !$OMP THREADPRIVATE(nm_oro_t)
36
37   INTEGER, SAVE, PROTECTED :: NKTOPG ! Security value for blocked flow level
38   !$OMP THREADPRIVATE(NKTOPG)
39
40   INTEGER, SAVE, PROTECTED :: NSTRA ! An estimate to qualify the upper levels of the model where one wants to impose stress
41   !$OMP THREADPRIVATE(NSTRA)
42
43   REAL, SAVE, PROTECTED :: GFRCRIT = 1. ! Critical Non-dimensional mountain Height (HNC in (1), LOTT 1999)
44   !$OMP THREADPRIVATE(GFRCRIT)
45
46   REAL, SAVE, PROTECTED :: GKWAKE = 0.5 ! Bluff-body drag coefficient for low level wake (Cd in (2), LOTT 1999)
47   !$OMP THREADPRIVATE(GKWAKE)
48
49   REAL, SAVE, PROTECTED :: GRCRIT ! Critical Richardson Number (Ric, End of first column p791 of LOTT 1999)
50   !$OMP THREADPRIVATE(GRCRIT)
51
52   REAL, SAVE, PROTECTED :: GVCRIT
53   !$OMP THREADPRIVATE(GVCRIT)
54
55   REAL, SAVE, PROTECTED :: GKDRAG ! Gravity wave drag coefficient (G in (3), LOTT 1999)
56   !$OMP THREADPRIVATE(GKDRAG)
57
58   REAL, SAVE, PROTECTED :: GKLIFT ! Mountain Lift coefficient (Cl in (4), LOTT 1999)
59   !$OMP THREADPRIVATE(GKLIFT)
60
61   REAL, SAVE, PROTECTED :: GHMAX ! Not used
62   !$OMP THREADPRIVATE(GHMAX)
63
64   REAL, SAVE, PROTECTED :: GRAHILO ! Set-up the trapped waves fraction (Beta, End of first column,  LOTT 1999)
65   !$OMP THREADPRIVATE(GRAHILO)
66
67   REAL, SAVE, PROTECTED :: GSIGCR ! Security value for blocked flow depth
68   !$OMP THREADPRIVATE(GSIGCR)
69
70   REAL, SAVE, PROTECTED :: GSSEC ! Security min value for low-level B-V frequency
71   !$OMP THREADPRIVATE(GSSEC)
72
73   REAL, SAVE, PROTECTED :: GTSEC ! Security min value for anisotropy and GW stress.
74   !$OMP THREADPRIVATE(GTSEC)
75
76   REAL, SAVE, PROTECTED :: GVSEC ! Security min value for ulow
77   !$OMP THREADPRIVATE(GVSEC)
78
79   ! for non-orographic gravity wave drag
80   LOGICAL, SAVE, PROTECTED :: ok_hines = .false. ! activate non-orographic gravity wave drag from Hines
81   !$OMP THREADPRIVATE(ok_hines)
82
83   LOGICAL, SAVE, PROTECTED :: ok_gwd_rando = .false.   ! activate non-orographic stochastic gravity wave drag params
84   !$OMP THREADPRIVATE(ok_gwd_rando)
85
86   REAL, SAVE, PROTECTED :: GWD_RANDO_RUWMAX = 2.    ! Maximum Eliassen-Palm flux at launch level, in "FLOTT_GWD_rando"
87   !$OMP THREADPRIVATE(GWD_RANDO_RUWMAX)
88
89   REAL, SAVE, PROTECTED :: GWD_RANDO_SAT = 0.25     ! saturation parameter in "FLOTT_GWD_rando"  S_c in equation (12) of Lott (JGR, vol 118, page 8897, 2013)
90   !$OMP THREADPRIVATE(GWD_RANDO_SAT)
91
92   REAL, SAVE, PROTECTED ::  GWD_FRONT_RUWMAX = 2.5  ! Same as GWD_RANDO params but for fronal GWs
93   !$OMP THREADPRIVATE(GWD_FRONT_RUWMAX)
94
95   REAL, SAVE, PROTECTED ::  GWD_FRONT_SAT = 0.60  ! Same as GWD_RANDO params but for fronal GWs
96   !$OMP THREADPRIVATE(GWD_FRONT_SAT)
97
98   INTEGER, SAVE, PROTECTED :: NK = 2
99   !$OMP THREADPRIVATE(NK)
100
101   INTEGER, SAVE, PROTECTED :: NP = 2
102   !$OMP THREADPRIVATE(NP)
103
104   INTEGER, SAVE, PROTECTED :: NO = 2
105   !$OMP THREADPRIVATE(NO)
106
107   INTEGER, SAVE, PROTECTED :: NA = 5 ! number of realizations to get the phase speed
108   !$OMP THREADPRIVATE(NA)
109
110   INTEGER, SAVE, PROTECTED :: NW
111   !$OMP THREADPRIVATE(NW)
112
113   ! for TKE production by SSO drag
114   INTEGER, SAVE, PROTECTED :: addtkeoro = 0     ! activate TKE production by sso drag
115   !$OMP THREADPRIVATE(addtkeoro)
116
117   LOGICAL, SAVE, PROTECTED :: smallscales_tkeoro = .false.  ! considers all sso scales for TKE production by sso drag if addtkeoro >=2
118   !$OMP THREADPRIVATE(smallscales_tkeoro)
119
120   REAL, SAVE, PROTECTED :: alphatkeoro = 1.    ! tuning parameter for TKE production by sso drag
121   !$OMP THREADPRIVATE(alphatkeoro)
122
123   ! Others
124
125   LOGICAL, SAVE :: firstcall = .true.
126   !$OMP THREADPRIVATE(firstcall)
127
128   LOGICAL, SAVE :: gwd_reproductibilite_mpiomp = .true.
129   !$OMP THREADPRIVATE(gwd_reproductibilite_mpiomp)
130
131CONTAINS
132
133!===================================================================================================================================================
134   SUBROUTINE gwd_ini(klon, klev, pplay, paprs, lunout_in, prt_level_in, &
135                      RCPD_in, RD_in, RG_in, RPI_in, ROMEGA_in, RA_in, &
136                      RKAPPA_in, RLVTT_in, ok_strato_in)
137
138      USE ioipsl_getin_p_mod, ONLY: getin_p
139
140      IMPLICIT NONE
141
142      INTEGER, INTENT(IN) :: klon, klev, lunout_in, prt_level_in
143      REAL, DIMENSION(klon, klev) :: pplay
144      REAL, DIMENSION(klon, klev + 1) :: paprs
145      REAL, INTENT(IN) :: RCPD_in, RD_in, RG_in, RPI_in, ROMEGA_in
146      REAL, INTENT(IN) :: RLVTT_IN, RA_in, RKAPPA_in
147      LOGICAL, INTENT(IN) :: ok_strato_in
148      CHARACTER(len=80) :: abort_message
149
150      lunout = lunout_in
151      prt_level = prt_level_in
152      RCPD = RCPD_in
153      RD = RD_in
154      RG = RG_in
155      RPI = RPI_in
156      ROMEGA = ROMEGA_in
157      RA = RA_in
158      RKAPPA = RKAPPA_in
159      RLVTT = RLVTT_in
160      ok_strato = ok_strato_in
161
162! default value of orographic gravity wave drag depend on ok_strato
163      IF (ok_strato) THEN
164         gkdrag = 0.1875
165         grahilo = 0.1
166         grcrit = 1.
167         gklift = 0.25
168      ELSE
169         gkdrag = 0.2
170         grahilo = 1.
171         grcrit = 0.01
172         gklift = 0.50
173      END IF
174
175! getin from .def files
176
177      CALL getin_p('ok_orodr', ok_orodr)
178      CALL getin_p('ok_orolf', ok_orolf)
179      CALL getin_p('ok_hines', ok_hines)
180      CALL getin_p('ok_gwd_rando', ok_gwd_rando)
181      CALL getin_p('zpmm_orodr_t', zpmm_orodr_t)
182      CALL getin_p('zpmm_orolf_t', zpmm_orolf_t)
183      CALL getin_p('zstd_orodr_t', zstd_orodr_t)
184      CALL getin_p('nm_oro_t', nm_oro_t)
185      CALL getin_p('addtkeoro', addtkeoro)
186      CALL getin_p('alphatkeoro', alphatkeoro)
187      CALL getin_p('smallscales_tkeoro', smallscales_tkeoro)
188      CALL getin_p('gwd_rando_ruwmax', gwd_rando_ruwmax)
189      CALL getin_p('gwd_rando_sat', gwd_rando_sat)
190      CALL getin_p('gwd_front_ruwmax', gwd_front_ruwmax)
191      CALL getin_p('gwd_front_sat', gwd_front_sat)
192      CALL getin_p('sso_gkdrag', gkdrag)
193      CALL getin_p('sso_grahil', grahilo)
194      CALL getin_p('sso_grcrit', grcrit)
195      CALL getin_p('sso_gfrcri', gfrcrit)
196      CALL getin_p('sso_gkwake', gkwake)
197      CALL getin_p('sso_gklift', gklift)
198
199! write in used.def files
200      WRITE (lunout, *) 'gwd_ini, ok_orodr:', ok_orodr
201      WRITE (lunout, *) 'gwd_ini, ok_orolf:', ok_orolf
202      WRITE (lunout, *) 'gwd_ini, ok_hines:', ok_hines
203      WRITE (lunout, *) 'gwd_ini, ok_gwd_rando:', ok_gwd_rando
204      WRITE (lunout, *) 'gwd_ini, zpmm_orodr_t:', zpmm_orodr_t
205      WRITE (lunout, *) 'gwd_ini, zstd_orodr_t:', zstd_orodr_t
206      WRITE (lunout, *) 'gwd_ini, zpmm_orolf_t:', zpmm_orolf_t
207      WRITE (lunout, *) 'gwd_ini, nm_oro_t:', nm_oro_t
208      WRITE (lunout, *) 'gwd_ini, addtkeoro:', addtkeoro
209      WRITE (lunout, *) 'gwd_ini, alphatkeoro:', alphatkeoro
210      WRITE (lunout, *) 'gwd_ini, smallscales_tkeoro:', smallscales_tkeoro
211      WRITE (lunout, *) 'gwd_ini, gwd_rando_ruwmax:', gwd_rando_ruwmax
212      WRITE (lunout, *) 'gwd_ini, gwd_rando_sat:', gwd_rando_sat
213      WRITE (lunout, *) 'gwd_ini, gwd_front_ruwmax:', gwd_front_ruwmax
214      WRITE (lunout, *) 'gwd_ini, gwd_front_sat:', gwd_front_sat
215      WRITE (lunout, *) 'gwd_ini, gklift:', gklift
216      WRITE (lunout, *) 'gwd_ini, gkwake:', gkwake
217      WRITE (lunout, *) 'gwd_ini, gfrcrit:', gfrcrit
218      WRITE (lunout, *) 'gwd_ini, grcrit:', grcrit
219      WRITE (lunout, *) 'gwd_ini, grahilo:', grahilo
220
221! few checks:
222
223      IF (klon == 1 .AND. ok_gwd_rando) THEN
224         print*, 'stochastic non-orographic gravity wave drag param does not work in 1D'
225         print*, 'I set ok_gwd_rando to false'
226         ok_gwd_rando=.false.
227      END IF
228
229! initialisation of ogwd parameters
230
231      IF (ok_strato) THEN
232         CALL sugwd_strato(klon, klev, paprs, pplay)
233      ELSE
234         CALL sugwd(klon, klev, paprs, pplay)
235      END IF
236
237! initialisation of non-orographic gwd parameterisation
238      NW = NK*NP*NO
239      IF (ok_gwd_rando) THEN
240         CALL gwd_rando_first(klev)
241      END IF
242
243   END SUBROUTINE gwd_ini
244
245!===================================================================================================================================================
246   SUBROUTINE sugwd_strato(nlon, nlev, paprs, pplay)
247
248      ! **** *SUGWD* INITIALIZE PARAMETERS CONTROLLING GRAVITY WAVE DRAG
249
250      ! PURPOSE.
251      ! --------
252      ! INITIALIZE PARAMETERS THAT CONTROLS THE
253      ! GRAVITY WAVE DRAG PARAMETRIZATION.
254      ! VERY IMPORTANT:
255      ! ______________
256      ! THIS ROUTINE SET_UP THE "TUNABLE PARAMETERS" OF THE
257      ! VARIOUS SSO SCHEMES
258
259      ! **   INTERFACE.
260      ! ----------
261      ! CALL *SUGWD* FROM *SUPHEC*
262      ! -----        ------
263
264      ! EXPLICIT ARGUMENTS :
265      ! --------------------
266      ! PAPRS,PPLAY : Pressure at semi and full model levels
267      ! NLEV        : number of model levels
268      ! NLON        : number of points treated in the physics
269
270      ! IMPLICIT ARGUMENTS :
271      ! --------------------
272      ! -GFRCRIT-R:  Critical Non-dimensional mountain Height
273      ! (HNC in (1),    LOTT 1999)
274      ! -GKWAKE--R:  Bluff-body drag coefficient for low level wake
275      ! (Cd in (2),     LOTT 1999)
276      ! -GRCRIT--R:  Critical Richardson Number
277      ! (Ric, End of first column p791 of LOTT 1999)
278      ! -GKDRAG--R:  Gravity wave drag coefficient
279      ! (G in (3),      LOTT 1999)
280      ! -GKLIFT--R:  Mountain Lift coefficient
281      ! (Cl in (4),     LOTT 1999)
282      ! -GHMAX---R:  Not used
283      ! -GRAHILO-R:  Set-up the trapped waves fraction
284      ! (Beta , End of first column,  LOTT 1999)
285
286      ! -GSIGCR--R:  Security value for blocked flow depth
287      ! -NKTOPG--I:  Security value for blocked flow level
288      ! -nstra----I:  An estimate to qualify the upper levels of
289      ! the model where one wants to impose strees
290      ! profiles
291      ! -GSSECC--R:  Security min value for low-level B-V frequency
292      ! -GTSEC---R:  Security min value for anisotropy and GW stress.
293      ! -GVSEC---R:  Security min value for ulow
294
295      ! METHOD.
296      ! -------
297      ! SEE DOCUMENTATION
298
299      ! EXTERNALS.
300      ! ----------
301      ! NONE
302
303      ! REFERENCE.
304      ! ----------
305      ! Lott, 1999: Alleviation of stationary biases in a GCM through...
306      ! Monthly Weather Review, 127, pp 788-801.
307
308      ! AUTHOR.
309      ! -------
310      ! FRANCOIS LOTT        *LMD*
311
312      ! MODIFICATIONS.
313      ! --------------
314      ! ORIGINAL : 90-01-01 (MARTIN MILLER, ECMWF)
315      ! LAST:  99-07-09     (FRANCOIS LOTT,LMD)
316      ! ------------------------------------------------------------------
317
318      USE dimphy
319      USE mod_phys_lmdz_para
320      USE mod_grid_phy_lmdz
321      USE geometry_mod
322
323      IMPLICIT NONE
324
325      ! -----------------------------------------------------------------
326      ! ----------------------------------------------------------------
327
328      ! ARGUMENTS
329      INTEGER nlon, nlev
330      REAL paprs(nlon, nlev + 1)
331      REAL pplay(nlon, nlev)
332
333      INTEGER jk
334      REAL zpr, ztop, zsigt, zpm1r
335      INTEGER :: cell, ij, nstra_tmp, nktopg_tmp
336      REAL :: current_dist, dist_min, dist_min_glo
337
338      ! *       1.    SET THE VALUES OF THE PARAMETERS
339      ! --------------------------------
340
341      ghmax = 10000.
342      zpr = 100000.
343      ZTOP = 0.00005
344      zsigt = 0.94
345
346!ym Take the point at equator close to (0,0) coordinates.
347      dist_min = 360
348      dist_min_glo = 360.
349      cell = -1
350      DO ij = 1, klon
351         current_dist = sqrt(longitude_deg(ij)**2 + latitude_deg(ij)**2)
352         current_dist = current_dist*(1 + (1e-10*ind_cell_glo(ij))/klon_glo) ! For point unicity
353         IF (dist_min > current_dist) THEN
354            dist_min = current_dist
355            cell = ij
356         END IF
357      END DO
358
359      CALL reduce_min(dist_min, dist_min_glo)
360      CALL bcast(dist_min_glo)
361      IF (dist_min /= dist_min_glo) cell = -1
362!ym in future find the point at equator close to (0,0) coordinates.
363      PRINT *, 'SUGWD distmin dist_min_glo cell=', dist_min, dist_min_glo, cell
364
365      nktopg_tmp = nktopg
366      nstra_tmp = nstra
367
368      IF (cell /= -1) THEN
369
370         DO jk = 1, nlev
371            !zpm1r = pplay(cell+1, jk)/paprs(cell+1, 1)
372            zpm1r = pplay(cell, jk)/paprs(cell, 1)
373            IF (zpm1r >= zsigt) THEN
374               nktopg_tmp = jk
375            END IF
376            IF (zpm1r >= ztop) THEN
377               nstra_tmp = jk
378            END IF
379         END DO
380      ELSE
381         nktopg_tmp = 0
382         nstra_tmp = 0
383      END IF
384
385      CALL reduce_sum(nktopg_tmp, nktopg)
386      CALL bcast(nktopg)
387      CALL reduce_sum(nstra_tmp, nstra)
388      CALL bcast(nstra)
389
390      ! inversion car dans orodrag on compte les niveaux a l'envers
391      nktopg = nlev - nktopg + 1
392      nstra = nlev - nstra
393      PRINT *, ' DANS SUGWD nktopg=', nktopg
394      PRINT *, ' DANS SUGWD nstra=', nstra
395      if (nstra == 0) call abort_physic("sugwd_strato", "no level in stratosphere", 1)
396
397      gsigcr = 0.80 ! Top of low level flow
398      gvcrit = 0.1
399
400      WRITE (UNIT=6, FMT='('' *** SSO essential constants ***'')')
401      WRITE (UNIT=6, FMT='('' *** SPECIFIED IN SUGWD ***'')')
402      WRITE (UNIT=6, FMT='('' Gravity wave ct '',E13.7,'' '')') gkdrag
403      WRITE (UNIT=6, FMT='('' Trapped/total wave dag '',E13.7,'' '')') grahilo
404      WRITE (UNIT=6, FMT='('' Critical Richardson   = '',E13.7,'' '')') grcrit
405      WRITE (UNIT=6, FMT='('' Critical Froude'',e13.7)') gfrcrit
406      WRITE (UNIT=6, FMT='('' Low level Wake bluff cte'',e13.7)') gkwake
407      WRITE (UNIT=6, FMT='('' Low level lift  cte'',e13.7)') gklift
408
409      ! ----------------------------------------------------------------
410
411      ! *       2.    SET VALUES OF SECURITY PARAMETERS
412      ! ---------------------------------
413
414      gvsec = 0.10
415      gssec = 0.0001
416      gtsec = 0.00001
417
418      RETURN
419   END SUBROUTINE sugwd_strato
420
421!===================================================================================================================================================
422   SUBROUTINE sugwd(nlon, nlev, paprs, pplay)
423
424      USE dimphy
425      USE mod_phys_lmdz_para
426      USE mod_grid_phy_lmdz
427      ! USE parallel
428
429      ! **** *SUGWD* INITIALIZE PARAMETERS CONTROLLING GRAVITY WAVE DRAG
430
431      ! PURPOSE.
432      ! --------
433      ! INITIALIZE PARAMETERS THAT CONTROLS THE
434      ! GRAVITY WAVE DRAG PARAMETRIZATION.
435
436      ! **   INTERFACE.
437      ! ----------
438      ! CALL *SUGWD* FROM *SUPHEC*
439      ! -----        ------
440
441      ! EXPLICIT ARGUMENTS :
442      ! --------------------
443      ! PSIG        : VERTICAL COORDINATE TABLE
444      ! NLEV        : NUMBER OF MODEL LEVELS
445
446      ! IMPLICIT ARGUMENTS :
447      ! --------------------
448      ! METHOD.
449      ! -------
450      ! SEE DOCUMENTATION
451
452      ! EXTERNALS.
453      ! ----------
454      ! NONE
455
456      ! REFERENCE.
457      ! ----------
458      ! ECMWF Research Department documentation of the IFS
459
460      ! AUTHOR.
461      ! -------
462      ! MARTIN MILLER             *ECMWF*
463
464      ! MODIFICATIONS.
465      ! --------------
466      ! ORIGINAL : 90-01-01
467      ! ------------------------------------------------------------------
468      IMPLICIT NONE
469
470      ! -----------------------------------------------------------------
471      ! ----------------------------------------------------------------
472
473      INTEGER nlon, nlev, jk
474      REAL paprs(nlon, nlev + 1)
475      REAL pplay(nlon, nlev)
476      REAL zpr, zstra, zsigt, zpm1r
477      REAL :: pplay_glo(klon_glo, nlev)
478      REAL :: paprs_glo(klon_glo, nlev + 1)
479
480      ! *       1.    SET THE VALUES OF THE PARAMETERS
481      ! --------------------------------
482
483      ghmax = 10000.
484      zpr = 100000.
485      zstra = 0.1
486      zsigt = 0.94
487
488      CALL gather(pplay, pplay_glo)
489      CALL bcast(pplay_glo)
490      CALL gather(paprs, paprs_glo)
491      CALL bcast(paprs_glo)
492
493      DO jk = 1, nlev
494         zpm1r = pplay_glo((klon_glo/2) + 1, jk)/paprs_glo((klon_glo/2) + 1, 1)
495         IF (zpm1r >= zsigt) THEN
496            nktopg = jk
497         END IF
498         zpm1r = pplay_glo((klon_glo/2) + 1, jk)/paprs_glo((klon_glo/2) + 1, 1)
499         IF (zpm1r >= zstra) THEN
500            nstra = jk
501         END IF
502      END DO
503
504      ! inversion car dans orodrag on compte les niveaux a l'envers
505      nktopg = nlev - nktopg + 1
506      nstra = nlev - nstra
507      PRINT *, ' DANS SUGWD nktopg=', nktopg
508      PRINT *, ' DANS SUGWD nstra=', nstra
509
510      gsigcr = 0.80
511      gvcrit = 0.0
512
513      ! ----------------------------------------------------------------
514
515      ! *       2.    SET VALUES OF SECURITY PARAMETERS
516      ! ---------------------------------
517
518      gvsec = 0.10
519      gssec = 1.E-12
520
521      gtsec = 1.E-07
522
523      ! ----------------------------------------------------------------
524
525      RETURN
526   END SUBROUTINE sugwd
527!===================================================================================================================================================
528
529   SUBROUTINE gwd_rando_first(klev)
530
531      USE ioipsl_getin_p_mod, ONLY: getin_p
532
533      IMPLICIT NONE
534      INTEGER, INTENT(IN) :: klev
535      CHARACTER(LEN=20), PARAMETER :: modname = 'gwd_rando_first'
536      CHARACTER(LEN=80) :: abort_message
537
538      IF (firstcall) THEN
539         ! Key to solve a non-reproductability issue. The aim is to test to get back to previous version
540         ! to remove asap
541         CALL getin_p('gwd_reproductibilite_mpiomp', gwd_reproductibilite_mpiomp)
542
543         IF (NW + 3*NA >= KLEV) THEN
544            abort_message = 'NW+3*NA>=KLEV Problem to generate waves associated with precip.'
545            CALL abort_physic(modname, abort_message, 1)
546         END IF
547
548         IF (.NOT. ok_hines) THEN
549            IF (NW + 4*(NA - 1) + NA >= KLEV) THEN
550               abort_message = 'NW+3*NA>=KLEV Problem to generate waves associated with fronts'
551               CALL abort_physic(modname, abort_message, 1)
552            END IF
553         END IF
554
555         firstcall = .false.
556      END IF
557   END SUBROUTINE gwd_rando_first
558
559!===================================================================================================================================================
560END MODULE lmdz_gwd_ini
Note: See TracBrowser for help on using the repository browser.