source: trunk/LMDZ.VENUS/libf/phyvenus/flott_gwd_ran.F90 @ 3461

Last change on this file since 3461 was 2381, checked in by slebonnois, 4 years ago

GG+TN+SL : includes a limitation in the horizontal wavelength of non-orographic gravity waves, related to GCM resolution

File size: 21.1 KB
Line 
1      SUBROUTINE FLOTT_GWD_RAN(NLON,NLEV,DTIME, pp, pn2,  &
2                  tt,uu,vv, plevmoy, &
3                  zustr,zvstr,d_t, d_u, d_v)
4
5    !----------------------------------------------------------------------
6    ! Parametrization of the momentum flux deposition due to a discrete
7    ! number of gravity waves.
8    ! F. Lott
9    ! Version 14, Gaussian distribution of the source
10    ! LMDz model online version     
11    ! ADAPTED FOR VENUS /  F. LOTT + S. LEBONNOIS
12    ! Version adapted on 03/04/2013:
13    !      - input flux compensated in the deepest layers
14    !---------------------------------------------------------------------
15
16      use dimphy
17      USE geometry_mod, only: cell_area, latitude_deg
18      implicit none
19
20#include "YOEGWD.h"
21#include "YOMCST.h"
22
23    ! 0. DECLARATIONS:
24
25    ! 0.1 INPUTS
26    INTEGER, intent(in):: NLON, NLEV
27    REAL, intent(in):: DTIME ! Time step of the Physics
28    REAL, intent(in):: pp(NLON, NLEV)   ! Pressure at full levels
29! VENUS ATTENTION: CP VARIABLE PN2 CALCULE EN AMONT DES PARAMETRISATIONS
30    REAL, intent(in):: pn2(NLON,NLEV)   ! N2 (BV^2) at 1/2 levels
31    REAL, intent(in):: TT(NLON, NLEV)   ! Temp at full levels
32
33    REAL, intent(in):: UU(NLON, NLEV) , VV(NLON, NLEV)
34    ! Hor winds at full levels
35    REAL, intent(in) :: plevmoy(NLEV+1) ! press (Pa) at interlayers, at klon/2+1
36
37    ! 0.2 OUTPUTS
38    REAL, intent(out):: zustr(NLON), zvstr(NLON) ! Surface Stresses
39    REAL, intent(inout):: d_t(NLON, NLEV)        ! Tendency on Temp.
40
41    REAL, intent(inout):: d_u(NLON, NLEV), d_v(NLON, NLEV)
42    ! Tendencies on winds
43
44    ! O.3 INTERNAL ARRAYS
45
46    INTEGER II, LL
47
48    ! 0.3.0 TIME SCALE OF THE LIFE CYCLE OF THE WAVES PARAMETERIZED
49
50    REAL DELTAT
51
52    ! 0.3.1 GRAVITY-WAVES SPECIFICATIONS
53
54!VENUS
55    INTEGER, PARAMETER:: NK = 2, NP = 2, NO = 2, NW = NK * NP * NO
56    !Online output: change NO
57!    INTEGER, PARAMETER:: NK = 1, NP = 2, NO = 11, NW = NK * NP * NO
58    INTEGER JK, JP, JO, JW
59    REAL KMIN, KMAX ! Min and Max horizontal wavenumbers
60    REAL KSTAR ! Min horizontal wavenumber constrained by resolution
61    REAL CMIN, CMAX ! Min and Max absolute ph. vel.
62    REAL CPHA ! absolute PHASE VELOCITY frequency
63    REAL ZK(NW, KLON) ! Horizontal wavenumber amplitude
64    REAL ZP(NW)       ! Horizontal wavenumber angle       
65    REAL ZO(NW, KLON) ! Absolute frequency      !
66
67    ! Waves Intr. freq. at the 1/2 lev surrounding the full level
68    REAL ZOM(NW, KLON), ZOP(NW, KLON)
69
70    ! Wave EP-fluxes at the 2 semi levels surrounding the full level
71    REAL WWM(NW, KLON), WWP(NW, KLON)
72    ! Fluxes X and Y for each waves at 1/2 Levels
73    REAL RUWP(NW, KLON), RVWP(NW, KLON)
74    REAL RUW(KLON, KLEV + 1) ! Flux x at semi levels
75    REAL RVW(KLON, KLEV + 1) ! Flux y at semi levels
76
77    REAL RUW0(NW, KLON) ! Fluxes at launching level
78    REAL RUWMAX         ! Max value
79    INTEGER LAUNCH      ! Launching altitude
80    REAL XLAUNCH        ! Control the launching altitude
81
82    ! 0.3.2 PARAMETERS OF WAVES DISSIPATIONS
83
84    REAL SAT  ! saturation parameter
85    REAL RDISS, ZOISEC ! COEFF DE DISSIPATION, SECURITY FOR INTRINSIC FREQ
86
87    ! 0.3.3 BACKGROUND FLOW AT 1/2 LEVELS AND VERTICAL COORDINATE
88
89    REAL H0bis(KLON, KLEV) ! Characteristic Height of the atmosphere
90    REAL H0 ! Characteristic Height of the atmosphere
91    REAL PR, TR ! Reference Pressure and Temperature
92
93    REAL ZH(KLON, KLEV + 1) ! Log-pressure altitude (constant H0)
94    REAL ZHbis(KLON, KLEV + 1) ! Log-pressure altitude (varying H)
95
96    REAL UH(KLON, KLEV + 1), VH(KLON, KLEV + 1) ! Winds at 1/2 levels
97    REAL PH(KLON, KLEV + 1) ! Pressure at 1/2 levels
98    REAL PSEC ! Security to avoid division by 0 pressure
99    REAL BV(KLON, KLEV + 1) ! Brunt Vaisala freq. (BVF) at 1/2 levels
100    REAL BVSEC ! Security to avoid negative BVF
101
102! COSMETICS TO DIAGNOSE EACH WAVES CONTRIBUTION.
103    logical output
104    data output/.false./
105!    data output/.true./
106 ! CAUTION ! IF output is .true. THEN change NO to 10 at least !
107    character*14 outform
108    character*2  str2
109    integer      ieq
110
111! ON CONSERVE LA MEMOIRE un certain temps AVEC UN SAVE
112    real,save,allocatable :: d_u_sav(:,:),d_v_sav(:,:)
113    LOGICAL firstcall
114    SAVE firstcall
115    DATA firstcall/.true./
116
117      REAL ALEAS
118      EXTERNAL ALEAS
119
120    !-----------------------------------------------------------------
121    !  1. INITIALISATIONS
122
123      IF (firstcall) THEN
124        allocate(d_u_sav(NLON,NLEV),d_v_sav(NLON,NLEV))
125        d_u_sav = 0.
126        d_v_sav = 0.
127        firstcall=.false.
128      ENDIF
129   
130    !    1.1 Basic parameter
131
132    !  PARAMETERS CORRESPONDING TO V3:
133    RUWMAX = 0.005      ! Max EP-Flux at Launch altitude (Nominal as Gilli2017)
134!    RUWMAX = 0.0025    ! Half-nominal
135!    RUWMAX = 0.01       ! Double-nominal
136    SAT    = 0.85       ! Saturation parameter: Sc in (12)
137    RDISS  = 0.1        ! Diffusion parameter
138
139    DELTAT=24.*3600.    ! Time scale of the waves (first introduced in 9b)
140
141!!!! TESTS
142!    KMIN = 1.E-5         ! Min horizontal wavenumber    (lambda max 628)
143!    KMIN = 6.3E-6       ! Min horizontal wavenumber    (lambda max 1000 km)
144!    KMIN = 1.2E-6       ! Min horizontal wavenumber     (lambda max 5000 km)
145!    KMAX = 1.E-4        ! Max horizontal wavenumber     (lambda min 50 km)
146
147!!! TEST Kobs  and BestFit Diogo
148    KMIN = 1.E-5         ! Min horizontal wavenumber    (lambda max 628)
149    KMAX = 1.E-4        ! Max horizontal wavenumber     (lambda min 50 km)
150
151!!!  TEST GW8
152
153!    KMIN = 1.25E-5       ! Min horizontal wavenumber     (lambda max 500 km)
154!    KMAX = 1.E-4        ! Max horizontal wavenumber     (lambda min 50 km)
155
156!!!  TEST GW9
157
158!    KMIN = 1.E-6       ! Min horizontal wavenumber     (lambda max 6000 km)
159!    KMAX = 2.E-5       ! Max horizontal wavenumber     (lambda min 300 km)
160
161!!!  TEST BestFit GG
162!    KMIN = 6.3E-6       ! Min horizontal wavenumber    (lambda max 1000 km)
163!    KMAX = 1.E-4        ! Max horizontal wavenumber     (lambda min 50 km)
164
165    !Online output: one value only
166    if (output) then
167      KMIN = 3E-5
168      KMAX = 3E-5
169    endif
170
171!!
172    CMIN = 1.           ! Min phase velocity
173    CMAX = 61.          ! Max phase speed velocity  TestGW6
174!     CMAX= 111.         ! TestGW2
175!    XLAUNCH=0.6        ! Parameter that control launching altitude
176     XLAUNCH=5e-3        ! Value for top of cloud convective region
177!    XLAUNCH=2.e-4         ! Value for ~ 70 km? Review this
178
179!    PR = 9.2e6          ! Reference pressure    ! VENUS!!
180    PR = 5e5            ! Reference pressure    ! VENUS: cloud layer
181    TR = 300.           ! Reference Temperature ! VENUS: cloud layer
182    H0 = RD * TR / RG   ! Characteristic vertical scale height
183
184    BVSEC  = 1.E-5      ! Security to avoid negative BVF 
185    PSEC   = 1.E-8      ! Security to avoid division by 0 pressure
186    ZOISEC = 1.E-8      ! Security FOR 0 INTRINSIC FREQ
187
188    IF(DELTAT.LT.DTIME)THEN
189       PRINT *,'GWD RANDO: DELTAT LT DTIME!'
190       STOP
191    ENDIF
192
193
194    IF (NLEV < NW) THEN
195       PRINT *, 'YOU WILL HAVE PROBLEM WITH RANDOM NUMBERS'
196       PRINT *, 'FLOTT GWD STOP'
197       STOP 1
198    ENDIF
199
200    !  2. EVALUATION OF THE BACKGROUND FLOW AT SEMI-LEVELS
201    !-------------------------------------------------------------
202
203    !Online output
204    if (output) OPEN(11,file="impact-gwno.dat")
205
206    ! Pressure and Inv of pressure, Temperature / at 1/2 level
207    DO LL = 2, KLEV
208       PH(:, LL) = EXP((LOG(PP(:, LL)) + LOG(PP(:, LL - 1))) / 2.)
209    end DO
210
211    PH(:, KLEV + 1) = 0.
212    PH(:, 1) = 2. * PP(:, 1) - PH(:, 2)
213
214    ! Launching altitude
215
216    DO LL = 1, NLEV
217       IF (plevmoy(LL) / plevmoy(1) > XLAUNCH) LAUNCH = LL
218    ENDDO
219! test
220!    print*,"launch=",LAUNCH
221!    print*,"launch p,N2=",plevmoy(LAUNCH),pn2(nlon/2+1,LAUNCH)
222
223    ! Log pressure vert. coordinate
224    DO LL = 1, KLEV + 1
225       ZH(:, LL) = H0 * LOG(PR / (PH(:, LL) + PSEC))
226    end DO
227
228    if (output) then
229    ! altitude above surface
230       ZHbis(:,1) = 0.   
231       DO LL = 2, KLEV + 1
232          H0bis(:, LL-1) = RD * TT(:, LL-1) / RG
233          ZHbis(:, LL) = ZHbis(:, LL-1) &
234           + H0bis(:, LL-1)*(PH(:, LL-1)-PH(:,LL))/PP(:, LL-1)
235       end DO
236    endif
237
238    ! Winds and BV frequency
239    DO LL = 2, KLEV
240       UH(:, LL) = 0.5 * (UU(:, LL) + UU(:, LL - 1)) ! Zonal wind
241       VH(:, LL) = 0.5 * (VV(:, LL) + VV(:, LL - 1)) ! Meridional wind
242       ! BVSEC: BV Frequency
243! VENUS ATTENTION: CP VARIABLE PSTAB CALCULE EN AMONT DES PARAMETRISATIONS
244       BV(:, LL) = MAX(BVSEC,SQRT(pn2(:,LL)))
245    end DO
246    BV(:, 1) = BV(:, 2)
247    UH(:, 1) = 0.
248    VH(:, 1) = 0.
249    BV(:, KLEV + 1) = BV(:, KLEV)
250    UH(:, KLEV + 1) = UU(:, KLEV)
251    VH(:, KLEV + 1) = VV(:, KLEV)
252
253
254    ! 3. WAVES CHARACTERISTICS CHOSEN RANDOMLY
255    !-------------------------------------------
256
257    ! The mod function of here a weird arguments
258    ! are used to produce the waves characteristics
259    ! in a stochastic way
260
261!! A REVOIR:
262!! - utilisation de MOD ou bien de aleas ?
263!! - distribution gaussienne des CPHA ? (avec signe ZP qui est ajuste apres)
264
265    JW = 0
266    DO JP = 1, NP
267       DO JK = 1, NK
268          DO JO = 1, NO
269             JW = JW + 1
270             !  Angle
271             ZP(JW) = 2. * RPI * REAL(JP - 1) / REAL(NP)
272             DO II = 1, KLON
273                ! Horizontal wavenumber amplitude
274!                ZK(JW, II) = KMIN + (KMAX - KMIN) * MOD(TT(II, JW) * 100., 1.)
275                ! TN+GG April/June 2020 - "Individual waves are not supposed to occupy
276                ! the entire domain, but only a fraction of it" [Lott+2012]:
277                KSTAR = RPI/SQRT(cell_area(II))     ! KSTAR should be < KMAX...
278                ZK(JW, II) = MAX(KMIN,KSTAR) + (KMAX - MAX(KMIN,KSTAR)) * ALEAS(0.)
279!               ZK(JW, II) = KMIN + (KMAX - KMIN) * ALEAS(0.)
280                ! Horizontal phase speed
281!                CPHA = CMIN + (CMAX - CMIN) * MOD(TT(II, JW)**2, 1.)
282                CPHA = CMIN + (CMAX - CMIN) * ALEAS(0.)
283       !Online output: linear
284                if (output) CPHA = CMIN + (CMAX - CMIN) * (JO-1)/(NO-1)
285                ! Intrinsic frequency
286                ZO(JW, II) = CPHA * ZK(JW, II)
287                ! Intrinsic frequency  is imposed
288                    ZO(JW, II) = ZO(JW, II)      &
289                  + ZK(JW, II) * COS(ZP(JW)) * UH(II, LAUNCH) &
290                  + ZK(JW, II) * SIN(ZP(JW)) * VH(II, LAUNCH)
291                ! Momentum flux at launch lev
292                ! RUW0(JW, II) = RUWMAX / REAL(NW) &
293                RUW0(JW, II) = RUWMAX &
294!                     * MOD(100. * (UU(II, JW)**2 + VV(II, JW)**2), 1.)
295                     * ALEAS(0.)
296       !Online output: fixed to max
297                if (output) RUW0(JW, II) = RUWMAX
298             ENDDO
299          end DO
300       end DO
301    end DO
302
303    ! 4. COMPUTE THE FLUXES
304    !--------------------------
305
306    !  4.1  Vertical velocity at launching altitude to ensure
307    !       the correct value to the imposed fluxes.
308    !
309    DO JW = 1, NW
310
311       ! Evaluate intrinsic frequency at launching altitude:
312       ZOP(JW, :) = ZO(JW, :) &
313            - ZK(JW, :) * COS(ZP(JW)) * UH(:, LAUNCH) &
314            - ZK(JW, :) * SIN(ZP(JW)) * VH(:, LAUNCH)
315       ! WW is directly a flux, here, not vertical velocity anymore
316       WWP(JW, :) = RUW0(JW,:)
317       RUWP(JW, :) = COS(ZP(JW)) * SIGN(1., ZOP(JW, :)) * RUW0(JW, :)
318       RVWP(JW, :) = SIN(ZP(JW)) * SIGN(1., ZOP(JW, :)) * RUW0(JW, :)
319
320    end DO
321
322    !  4.2 Initial flux at launching altitude
323
324    RUW(:, LAUNCH) = 0
325    RVW(:, LAUNCH) = 0
326    DO JW = 1, NW
327       RUW(:, LAUNCH) = RUW(:, LAUNCH) + RUWP(JW, :)
328       RVW(:, LAUNCH) = RVW(:, LAUNCH) + RVWP(JW, :)
329    end DO
330
331    !  4.3 Loop over altitudes, with passage from one level to the
332    !      next done by i) conserving the EP flux, ii) dissipating
333    !      a little, iii) testing critical levels, and vi) testing
334    !      the breaking.
335
336    !Online output
337    if (output) then
338        ieq=nlon/2+1
339        write(str2,'(i2)') NW+2
340        outform="("//str2//"(E12.4,1X))"
341        WRITE(11,outform) ZH(IEQ, 1) / 1000., ZHbis(IEQ, 1) / 1000., &
342               (ZO(JW, IEQ)/ZK(JW, IEQ)*COS(ZP(JW)), JW = 1, NW)
343    endif
344
345    DO LL = LAUNCH, KLEV - 1
346
347
348       !  W(KB)ARNING: ALL THE PHYSICS IS HERE (PASSAGE FROM ONE LEVEL
349       ! TO THE NEXT)
350       DO JW = 1, NW
351          ZOM(JW, :) = ZOP(JW, :)
352          WWM(JW, :) = WWP(JW, :)
353          ! Intrinsic Frequency
354          ZOP(JW, :) = ZO(JW, :) - ZK(JW, :) * COS(ZP(JW)) * UH(:, LL + 1) &
355               - ZK(JW, :) * SIN(ZP(JW)) * VH(:, LL + 1)
356
357          WWP(JW, :) = MIN( &
358       ! No breaking (Eq.6)
359               WWM(JW, :) &
360       ! Dissipation (Eq. 8):
361               * EXP(- RDISS * PR / (PH(:, LL + 1) + PH(:, LL)) &
362               * ((BV(:, LL + 1) + BV(:, LL)) / 2.)**3 &
363               / MAX(ABS(ZOP(JW, :) + ZOM(JW, :)) / 2., ZOISEC)**4 &
364               * ZK(JW, :)**3 * (ZH(:, LL + 1) - ZH(:, LL))), &
365       ! Critical levels (forced to zero if intrinsic frequency changes sign)
366               MAX(0., SIGN(1., ZOP(JW, :) * ZOM(JW, :))) &
367       ! Saturation (Eq. 12)
368               * ABS(ZOP(JW, :))**3 /BV(:, LL+1) &
369               * EXP(-ZH(:, LL + 1)/H0) * SAT**2*KMIN**2/ZK(JW, :)**4) 
370       end DO
371
372       ! END OF W(KB)ARNING
373       ! Evaluate EP-flux from Eq. 7 and
374       ! Give the right orientation to the stress
375
376       DO JW = 1, NW
377          RUWP(JW, :) = SIGN(1.,ZOP(JW, :))*COS(ZP(JW))*WWP(JW, :)
378          RVWP(JW, :) = SIGN(1.,ZOP(JW, :))*SIN(ZP(JW))*WWP(JW, :)
379       end DO
380       !
381       RUW(:, LL + 1) = 0.
382       RVW(:, LL + 1) = 0.
383
384       DO JW = 1, NW
385          RUW(:, LL + 1) = RUW(:, LL + 1) + RUWP(JW, :)
386          RVW(:, LL + 1) = RVW(:, LL + 1) + RVWP(JW, :)
387       end DO
388       !Online output
389       if (output) then
390         do JW=1,NW
391            if(RUWP(JW, IEQ).gt.0.) then
392              RUWP(JW, IEQ) = max(RUWP(JW, IEQ), 1.e-99)
393            else
394              RUWP(JW, IEQ) = min(RUWP(JW, IEQ), -1.e-99)
395            endif
396         enddo
397                   WRITE(11,outform) ZH(IEQ, LL+1) / 1000., &
398                                  ZHbis(IEQ, LL+1) / 1000., &
399                                  (RUWP(JW, IEQ), JW = 1, NW)
400       endif
401
402    end DO
403
404    ! 5 CALCUL DES TENDANCES:
405    !------------------------
406
407    ! 5.1 Rectification des flux au sommet et dans les basses couches:
408! MODIF SL
409 
410! Attention, ici c'est le total sur toutes les ondes...
411
412    RUW(:, KLEV + 1) = 0.
413    RVW(:, KLEV + 1) = 0.
414
415    ! Here, big change compared to FLott version:
416    ! We compensate (RUW(:, LAUNCH), ie total emitted upward flux
417    !  over the layers max(1,LAUNCH-3) to LAUNCH-1
418    DO LL = 1, max(1,LAUNCH-3)
419      RUW(:, LL) = 0.
420      RVW(:, LL) = 0.
421    end DO
422    DO LL = max(2,LAUNCH-2), LAUNCH-1
423       RUW(:, LL) = RUW(:, LL - 1) + RUW(:, LAUNCH) * &
424            (PH(:,LL)-PH(:,LL-1)) / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
425       RVW(:, LL) = RVW(:, LL - 1) + RVW(:, LAUNCH) * &
426            (PH(:,LL)-PH(:,LL-1)) / (PH(:,LAUNCH)-PH(:,max(1,LAUNCH-3)))
427    end DO
428    ! This way, the total flux from GW is zero, but there is a net transport
429    ! (upward) that should be compensated by circulation
430    ! and induce additional friction at the surface
431
432    !Online output
433    if (output) then
434       DO LL = 1, KLEV - 1
435           WRITE(11,*) ZHbis(IEQ, LL)/1000.,RUW(IEQ,LL)
436       end DO
437       CLOSE(11)
438       stop
439    endif
440
441    ! AR-1 RECURSIVE FORMULA (13) IN VERSION 4
442    DO LL = 1, KLEV
443       d_u(:, LL) = RG * (RUW(:, LL + 1) - RUW(:, LL)) &
444            / (PH(:, LL + 1) - PH(:, LL)) * DTIME
445       d_v(:, LL) = RG * (RVW(:, LL + 1) - RVW(:, LL)) &
446            / (PH(:, LL + 1) - PH(:, LL)) * DTIME
447    ENDDO
448        d_t = 0.
449    ! ON CONSERVE LA MEMOIRE un certain temps AVEC UN SAVE
450        d_u = DTIME/DELTAT/REAL(NW) * d_u + (1.-DTIME/DELTAT) * d_u_sav
451        d_v = DTIME/DELTAT/REAL(NW) * d_v + (1.-DTIME/DELTAT) * d_v_sav
452        d_u_sav = d_u
453        d_v_sav = d_v
454
455    ! Cosmetic: evaluation of the cumulated stress
456
457    ZUSTR(:) = 0.
458    ZVSTR(:) = 0.
459    DO LL = 1, KLEV
460       ZUSTR(:) = ZUSTR(:) + D_U(:, LL) / RG * (PH(:, LL + 1) - PH(:, LL))
461       ZVSTR(:) = ZVSTR(:) + D_V(:, LL) / RG * (PH(:, LL + 1) - PH(:, LL))
462    ENDDO
463
464  END SUBROUTINE FLOTT_GWD_RAN
465
466!===================================================================
467!===================================================================
468!===================================================================
469!===================================================================
470
471      FUNCTION ALEAS (R)
472!***BEGIN PROLOGUE  ALEAS
473!***PURPOSE  Generate a uniformly distributed random number.
474!***LIBRARY   SLATEC (FNLIB)
475!***CATEGORY  L6A21
476!***TYPE      SINGLE PRECISION (ALEAS-S)
477!***KEYWORDS  FNLIB, ALEAS NUMBER, SPECIAL FUNCTIONS, UNIFORM
478!***AUTHOR  Fullerton, W., (LANL)
479!***DESCRIPTION
480!
481!      This pseudo-random number generator is portable among a wide
482! variety of computers.  RAND(R) undoubtedly is not as good as many
483! readily available installation dependent versions, and so this
484! routine is not recommended for widespread usage.  Its redeeming
485! feature is that the exact same random numbers (to within final round-
486! off error) can be generated from machine to machine.  Thus, programs
487! that make use of random numbers can be easily transported to and
488! checked in a new environment.
489!
490!      The random numbers are generated by the linear congruential
491! method described, e.g., by Knuth in Seminumerical Methods (p.9),
492! Addison-Wesley, 1969.  Given the I-th number of a pseudo-random
493! sequence, the I+1 -st number is generated from
494!             X(I+1) = (A*X(I) + C) MOD M,
495! where here M = 2**22 = 4194304, C = 1731 and several suitable values
496! of the multiplier A are discussed below.  Both the multiplier A and
497! random number X are represented in double precision as two 11-bit
498! words.  The constants are chosen so that the period is the maximum
499! possible, 4194304.
500!
501!      In order that the same numbers be generated from machine to
502! machine, it is necessary that 23-bit integers be reducible modulo
503! 2**11 exactly, that 23-bit integers be added exactly, and that 11-bit
504! integers be multiplied exactly.  Furthermore, if the restart option
505! is used (where R is between 0 and 1), then the product R*2**22 =
506! R*4194304 must be correct to the nearest integer.
507!
508!      The first four random numbers should be .0004127026,
509! .6750836372, .1614754200, and .9086198807.  The tenth random number
510! is .5527787209, and the hundredth is .3600893021 .  The thousandth
511! number should be .2176990509 .
512!
513!      In order to generate several effectively independent sequences
514! with the same generator, it is necessary to know the random number
515! for several widely spaced calls.  The I-th random number times 2**22,
516! where I=K*P/8 and P is the period of the sequence (P = 2**22), is
517! still of the form L*P/8.  In particular we find the I-th random
518! number multiplied by 2**22 is given by
519! I   =  0  1*P/8  2*P/8  3*P/8  4*P/8  5*P/8  6*P/8  7*P/8  8*P/8
520! RAND=  0  5*P/8  2*P/8  7*P/8  4*P/8  1*P/8  6*P/8  3*P/8  0
521! Thus the 4*P/8 = 2097152 random number is 2097152/2**22.
522!
523!      Several multipliers have been subjected to the spectral test
524! (see Knuth, p. 82).  Four suitable multipliers roughly in order of
525! goodness according to the spectral test are
526!    3146757 = 1536*2048 + 1029 = 2**21 + 2**20 + 2**10 + 5
527!    2098181 = 1024*2048 + 1029 = 2**21 + 2**10 + 5
528!    3146245 = 1536*2048 +  517 = 2**21 + 2**20 + 2**9 + 5
529!    2776669 = 1355*2048 + 1629 = 5**9 + 7**7 + 1
530!
531!      In the table below LOG10(NU(I)) gives roughly the number of
532! random decimal digits in the random numbers considered I at a time.
533! C is the primary measure of goodness.  In both cases bigger is better.
534!
535!                   LOG10 NU(I)              C(I)
536!       A       I=2  I=3  I=4  I=5    I=2  I=3  I=4  I=5
537!
538!    3146757    3.3  2.0  1.6  1.3    3.1  1.3  4.6  2.6
539!    2098181    3.3  2.0  1.6  1.2    3.2  1.3  4.6  1.7
540!    3146245    3.3  2.2  1.5  1.1    3.2  4.2  1.1  0.4
541!    2776669    3.3  2.1  1.6  1.3    2.5  2.0  1.9  2.6
542!   Best
543!    Possible   3.3  2.3  1.7  1.4    3.6  5.9  9.7  14.9
544!
545!             Input Argument --
546! R      If R=0., the next random number of the sequence is generated.
547!        If R .LT. 0., the last generated number will be returned for
548!          possible use in a restart procedure.
549!        If R .GT. 0., the sequence of random numbers will start with
550!          the seed R mod 1.  This seed is also returned as the value of
551!          RAND provided the arithmetic is done exactly.
552!
553!             Output Value --
554! RAND   a pseudo-random number between 0. and 1.
555!
556!***REFERENCES  (NONE)
557!***ROUTINES CALLED  (NONE)
558!***REVISION HISTORY  (YYMMDD)
559!   770401  DATE WRITTEN
560!   890531  Changed all specific intrinsics to generic.  (WRB)
561!   890531  REVISION DATE from Version 3.2
562!   891214  Prologue converted to Version 4.0 format.  (BAB)
563!***END PROLOGUE  RAND
564      SAVE IA1, IA0, IA1MA0, IC, IX1, IX0
565      DATA IA1, IA0, IA1MA0 /1536, 1029, 507/
566      DATA IC /1731/
567      DATA IX1, IX0 /0, 0/
568!***FIRST EXECUTABLE STATEMENT  RAND
569!
570!           A*X = 2**22*IA1*IX1 + 2**11*(IA1*IX1 + (IA1-IA0)*(IX0-IX1)
571!                   + IA0*IX0) + IA0*IX0
572!
573      IF (R.EQ.0.) THEN
574       IY0 = IA0*IX0
575       IY1 = IA1*IX1 + IA1MA0*(IX0-IX1) + IY0
576       IY0 = IY0 + IC
577       IX0 = MOD (IY0, 2048)
578       IY1 = IY1 + (IY0-IX0)/2048
579       IX1 = MOD (IY1, 2048)
580      ENDIF
581
582      IF (R.GT.0.) THEN
583       IX1 = MOD(R,1.)*4194304. + 0.5
584       IX0 = MOD (IX1, 2048)
585       IX1 = (IX1-IX0)/2048
586      ENDIF
587
588      ALEAS = IX1*2048 + IX0
589      ALEAS = ALEAS / 4194304.
590      RETURN
591
592      END
593
594
Note: See TracBrowser for help on using the repository browser.