source: LMDZ6/trunk/libf/phylmd/hines_gwd.f90 @ 5308

Last change on this file since 5308 was 5285, checked in by abarral, 4 days ago

As discussed internally, remove generic ONLY: ... for new _mod_h modules

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 62.1 KB
Line 
1
2! $Id: hines_gwd.f90 5285 2024-10-28 13:33:29Z fairhead $
3
4SUBROUTINE hines_gwd(nlon, nlev, dtime, paphm1x, papm1x, rlat, tx, ux, vx, &
5    zustrhi, zvstrhi, d_t_hin, d_u_hin, d_v_hin)
6
7  ! ########################################################################
8  ! Parametrization of the momentum flux deposition due to a broad band
9  ! spectrum of gravity waves, following Hines (1997a,b), as coded by
10  ! McLANDRESS (1995). Modified by McFARLANE and MANZINI (1995-1997)
11  ! MAECHAM model stand alone version
12  ! ########################################################################
13
14
15  USE dimphy
16  USE yomcst_mod_h
17IMPLICIT NONE
18
19  include "YOEGWD.h"
20
21
22  INTEGER nazmth
23  PARAMETER (nazmth=8)
24
25  ! INPUT ARGUMENTS.
26  ! ----- ----------
27
28  ! - 2D
29  ! PAPHM1   : HALF LEVEL PRESSURE (T-DT)
30  ! PAPM1    : FULL LEVEL PRESSURE (T-DT)
31  ! PTM1     : TEMPERATURE (T-DT)
32  ! PUM1     : ZONAL WIND (T-DT)
33  ! PVM1     : MERIDIONAL WIND (T-DT)
34
35
36  ! REFERENCE.
37  ! ----------
38  ! SEE MODEL DOCUMENTATION
39
40  ! AUTHOR.
41  ! -------
42
43  ! N. MCFARLANE   DKRZ-HAMBURG   MAY 1995
44  ! STAND ALONE E. MANZINI MPI-HAMBURG FEBRUARY 1997
45
46  ! BASED ON A COMBINATION OF THE OROGRAPHIC SCHEME BY N.MCFARLANE 1987
47  ! AND THE HINES SCHEME AS CODED BY C. MCLANDRESS 1995.
48
49
50
51  ! ym      INTEGER KLEVM1
52
53  REAL paphm1(klon, klev+1), papm1(klon, klev)
54  REAL ptm1(klon, klev), pum1(klon, klev), pvm1(klon, klev)
55  REAL prflux(klon)
56  ! 1
57  ! 1
58  ! 1
59  REAL rlat(klon), coslat(klon)
60
61  REAL th(klon, klev), utendgw(klon, klev), vtendgw(klon, klev), &
62    pressg(klon), uhs(klon, klev), vhs(klon, klev), zpr(klon)
63
64  ! * VERTICAL POSITIONING ARRAYS.
65
66  REAL sgj(klon, klev), shj(klon, klev), shxkj(klon, klev), dsgj(klon, klev)
67
68  ! * LOGICAL SWITCHES TO CONTROL ROOF DRAG, ENVELOP GW DRAG AND
69  ! * HINES' DOPPLER SPREADING EXTROWAVE GW DRAG.
70  ! * LOZPR IS TRUE FOR ZPR ENHANCEMENT
71
72
73  ! * WORK ARRAYS.
74
75  REAL m_alpha(klon, klev, nazmth), v_alpha(klon, klev, nazmth), &
76    sigma_alpha(klon, klev, nazmth), sigsqh_alpha(klon, klev, nazmth), &
77    drag_u(klon, klev), drag_v(klon, klev), flux_u(klon, klev), &
78    flux_v(klon, klev), heat(klon, klev), diffco(klon, klev), &
79    bvfreq(klon, klev), density(klon, klev), sigma_t(klon, klev), &
80    visc_mol(klon, klev), alt(klon, klev), sigsqmcw(klon, klev, nazmth), &
81    sigmatm(klon, klev), ak_alpha(klon, nazmth), k_alpha(klon, nazmth), &
82    mmin_alpha(klon, nazmth), i_alpha(klon, nazmth), rmswind(klon), &
83    bvfbot(klon), densbot(klon)
84  REAL smoothr1(klon, klev), smoothr2(klon, klev)
85  REAL sigalpmc(klon, klev, nazmth)
86  REAL f2mod(klon, klev)
87
88  ! * THES ARE THE INPUT PARAMETERS FOR HINES ROUTINE AND
89  ! * ARE SPECIFIED IN ROUTINE HINES_SETUP. SINCE THIS IS CALLED
90  ! * ONLY AT FIRST CALL TO THIS ROUTINE THESE VARIABLES MUST BE SAVED
91  ! * FOR USE AT SUBSEQUENT CALLS. THIS CAN BE AVOIDED BY CALLING
92  ! * HINES_SETUP IN MAIN PROGRAM AND PASSING THE PARAMETERS AS
93  ! * SUBROUTINE ARGUEMENTS.
94
95
96  REAL rmscon
97  INTEGER nmessg, iprint, ilrms
98  INTEGER ifl
99
100  INTEGER naz, icutoff, nsmax, iheatcal
101  REAL slope, f1, f2, f3, f5, f6, kstar(klon), alt_cutoff, smco
102
103  ! PROVIDED AS INPUT
104
105  INTEGER nlon, nlev
106
107  REAL dtime
108  REAL paphm1x(nlon, nlev+1), papm1x(nlon, nlev)
109  REAL ux(nlon, nlev), vx(nlon, nlev), tx(nlon, nlev)
110
111  ! VARIABLES FOR OUTPUT
112
113
114  REAL d_t_hin(nlon, nlev), d_u_hin(nlon, nlev), d_v_hin(nlon, nlev)
115  REAL zustrhi(nlon), zvstrhi(nlon)
116
117
118  ! * LOGICAL SWITCHES TO CONTROL PRECIP ENHANCEMENT AND
119  ! * HINES' DOPPLER SPREADING EXTROWAVE GW DRAG.
120  ! * LOZPR IS TRUE FOR ZPR ENHANCEMENT
121
122  LOGICAL lozpr, lorms(klon)
123
124  ! LOCAL PARAMETERS TO MAKE THINGS WORK (TEMPORARY VARIABLE)
125
126  REAL rhoh2o, zpcons, rgocp, zlat, dttdsf, ratio, hscal
127  INTEGER i, j, l, jl, jk, le, lref, lrefp, levbot
128
129  ! DATA PARAMETERS NEEDED, EXPLAINED LATER
130
131  REAL v0, vmin, dmpscal, taufac, hmin, apibt, cpart, fcrit
132  REAL pcrit, pcons
133  INTEGER iplev, ierror
134
135
136
137  ! PRINT *,' IT IS STARTED HINES GOING ON...'
138
139
140
141
142  ! *    COMPUTATIONAL CONSTANTS.
143  ! ------------- ----------
144
145
146  d_t_hin(:, :) = 0.
147
148  rhoh2o = 1000.
149  zpcons = (1000.*86400.)/rhoh2o
150  ! ym      KLEVM1=KLEV-1
151
152
153  DO jl = kidia, kfdia
154    paphm1(jl, 1) = paphm1x(jl, klev+1)
155    DO jk = 1, klev
156      le = klev + 1 - jk
157      paphm1(jl, jk+1) = paphm1x(jl, le)
158      papm1(jl, jk) = papm1x(jl, le)
159      ptm1(jl, jk) = tx(jl, le)
160      pum1(jl, jk) = ux(jl, le)
161      pvm1(jl, jk) = vx(jl, le)
162    END DO
163  END DO
164
165  ! Define constants and arrays needed for the ccc/mam gwd scheme
166  ! *Constants:
167
168  rgocp = rd/rcpd
169  lrefp = klev - 1
170  lref = klev - 2
171  ! 1
172  ! 1    *Arrays
173  ! 1
174  DO jk = 1, klev
175    DO jl = kidia, kfdia
176      shj(jl, jk) = papm1(jl, jk)/paphm1(jl, klev+1)
177      sgj(jl, jk) = papm1(jl, jk)/paphm1(jl, klev+1)
178      dsgj(jl, jk) = (paphm1(jl,jk+1)-paphm1(jl,jk))/paphm1(jl, klev+1)
179      shxkj(jl, jk) = (papm1(jl,jk)/paphm1(jl,klev+1))**rgocp
180      th(jl, jk) = ptm1(jl, jk)
181    END DO
182  END DO
183
184  ! C
185  DO jl = kidia, kfdia
186    pressg(jl) = paphm1(jl, klev+1)
187  END DO
188
189
190  DO jl = kidia, kfdia
191    prflux(jl) = 0.0
192    zpr(jl) = zpcons*prflux(jl)
193    zlat = (rlat(jl)/180.)*rpi
194    coslat(jl) = cos(zlat)
195  END DO
196
197  ! /#########################################################################
198  ! /
199  ! /
200
201  ! * AUG. 14/95 - C. MCLANDRESS.
202  ! * SEP.    95   N. MCFARLANE.
203
204  ! * THIS ROUTINE CALCULATES THE HORIZONTAL WIND TENDENCIES
205  ! * DUE TO MCFARLANE'S OROGRAPHIC GW DRAG SCHEME, HINES'
206  ! * DOPPLER SPREAD SCHEME FOR "EXTROWAVES" AND ADDS ON
207  ! * ROOF DRAG. IT IS BASED ON THE ROUTINE GWDFLX8.
208
209  ! * LREFP IS THE INDEX OF THE MODEL LEVEL BELOW THE REFERENCE LEVEL
210  ! * I/O ARRAYS PASSED FROM MAIN.
211  ! * (PRESSG = SURFACE PRESSURE)
212
213
214
215
216  ! * CONSTANTS VALUES DEFINED IN DATA STATEMENT ARE :
217  ! * VMIN     = MIMINUM WIND IN THE DIRECTION OF REFERENCE LEVEL
218  ! *            WIND BEFORE WE CONSIDER BREAKING TO HAVE OCCURED.
219  ! * DMPSCAL  = DAMPING TIME FOR GW DRAG IN SECONDS.
220  ! * TAUFAC   = 1/(LENGTH SCALE).
221  ! * HMIN     = MIMINUM ENVELOPE HEIGHT REQUIRED TO PRODUCE GW DRAG.
222  ! * V0       = VALUE OF WIND THAT APPROXIMATES ZERO.
223
224
225  DATA vmin/5.0/, v0/1.E-10/, taufac/5.E-6/, hmin/40000./, dmpscal/6.5E+6/, &
226    apibt/1.5708/, cpart/0.7/, fcrit/1./
227
228  ! * HINES EXTROWAVE GWD CONSTANTS DEFINED IN DATA STATEMENT ARE:
229  ! * RMSCON = ROOT MEAN SQUARE GRAVITY WAVE WIND AT LOWEST LEVEL (M/S).
230  ! * NMESSG  = UNIT NUMBER FOR PRINTED MESSAGES.
231  ! * IPRINT  = 1 TO DO PRINT OUT SOME HINES ARRAYS.
232  ! * IFL     = FIRST CALL FLAG TO HINES_SETUP ("SAVE" IT)
233  ! * PCRIT = CRITICAL VALUE OF ZPR (MM/D)
234  ! * IPLEV = LEVEL OF APPLICATION OF PRCIT
235  ! * PCONS = FACTOR OF ZPR ENHANCEMENT
236
237
238  DATA pcrit/5./, pcons/4.75/
239
240  iplev = lrefp - 1
241
242  DATA rmscon/1.00/iprint/2/, nmessg/6/
243  DATA ifl/0/
244
245  lozpr = .FALSE.
246
247  ! -----------------------------------------------------------------------
248
249
250
251  ! * SET ERROR FLAG
252
253  ierror = 0
254
255  ! * SPECIFY VARIOUS PARAMETERS FOR HINES ROUTINE AT VERY FIRST CALL.
256  ! * (NOTE THAT ARRAY K_ALPHA IS SPECIFIED SO MAKE SURE THAT
257  ! * IT IS NOT OVERWRITTEN LATER ON).
258
259  CALL hines_setup(naz, slope, f1, f2, f3, f5, f6, kstar, icutoff, &
260    alt_cutoff, smco, nsmax, iheatcal, k_alpha, ierror, nmessg, klon, nazmth, &
261    coslat)
262  IF (ierror/=0) GO TO 999
263
264  ! * START GWD CALCULATIONS.
265
266  lref = lrefp - 1
267
268
269  DO j = 1, nazmth
270    DO l = 1, klev
271      DO i = kidia, klon
272        sigsqmcw(i, l, j) = 0.
273      END DO
274    END DO
275  END DO
276
277
278
279  ! * INITIALIZE NECESSARY ARRAYS.
280
281  DO l = 1, klev
282    DO i = kidia, kfdia
283      utendgw(i, l) = 0.
284      vtendgw(i, l) = 0.
285
286      uhs(i, l) = 0.
287      vhs(i, l) = 0.
288
289    END DO
290  END DO
291
292  ! * IF USING HINES SCHEME THEN CALCULATE B V FREQUENCY AT ALL POINTS
293  ! * AND SMOOTH BVFREQ.
294
295  DO l = 2, klev
296    DO i = kidia, kfdia
297      dttdsf = (th(i,l)/shxkj(i,l)-th(i,l-1)/shxkj(i,l-1))/ &
298        (shj(i,l)-shj(i,l-1))
299      dttdsf = min(dttdsf, -5./sgj(i,l))
300      bvfreq(i, l) = sqrt(-dttdsf*sgj(i,l)*(sgj(i,l)**rgocp)/rd)*rg/ptm1(i, l &
301        )
302    END DO
303  END DO
304  DO l = 1, klev
305    DO i = kidia, kfdia
306      IF (l==1) THEN
307        bvfreq(i, l) = bvfreq(i, l+1)
308      END IF
309      IF (l>1) THEN
310        ratio = 5.*log(sgj(i,l)/sgj(i,l-1))
311        bvfreq(i, l) = (bvfreq(i,l-1)+ratio*bvfreq(i,l))/(1.+ratio)
312      END IF
313    END DO
314  END DO
315
316  ! * CALCULATE GW DRAG DUE TO HINES' EXTROWAVES
317  ! * SET MOLECULAR VISCOSITY TO A VERY SMALL VALUE.
318  ! * IF THE MODEL TOP IS GREATER THAN 100 KM THEN THE ACTUAL
319  ! * VISCOSITY COEFFICIENT COULD BE SPECIFIED HERE.
320
321  DO l = 1, klev
322    DO i = kidia, kfdia
323      visc_mol(i, l) = 1.5E-5
324      drag_u(i, l) = 0.
325      drag_v(i, l) = 0.
326      flux_u(i, l) = 0.
327      flux_v(i, l) = 0.
328      heat(i, l) = 0.
329      diffco(i, l) = 0.
330    END DO
331  END DO
332
333  ! * ALTITUDE AND DENSITY AT BOTTOM.
334
335  DO i = kidia, kfdia
336    hscal = rd*ptm1(i, klev)/rg
337    density(i, klev) = sgj(i, klev)*pressg(i)/(rg*hscal)
338    alt(i, klev) = 0.
339  END DO
340
341  ! * ALTITUDE AND DENSITY AT REMAINING LEVELS.
342
343  DO l = klev - 1, 1, -1
344    DO i = kidia, kfdia
345      hscal = rd*ptm1(i, l)/rg
346      alt(i, l) = alt(i, l+1) + hscal*dsgj(i, l)/sgj(i, l)
347      density(i, l) = sgj(i, l)*pressg(i)/(rg*hscal)
348    END DO
349  END DO
350
351
352  ! * INITIALIZE SWITCHES FOR HINES GWD CALCULATION
353
354  ilrms = 0
355
356  DO i = kidia, kfdia
357    lorms(i) = .FALSE.
358  END DO
359
360
361  ! * DEFILE BOTTOM LAUNCH LEVEL
362
363  levbot = iplev
364
365  ! * BACKGROUND WIND MINUS VALUE AT BOTTOM LAUNCH LEVEL.
366
367  DO l = 1, levbot
368    DO i = kidia, kfdia
369      uhs(i, l) = pum1(i, l) - pum1(i, levbot)
370      vhs(i, l) = pvm1(i, l) - pvm1(i, levbot)
371    END DO
372  END DO
373
374  ! * SPECIFY ROOT MEAN SQUARE WIND AT BOTTOM LAUNCH LEVEL.
375
376  DO i = kidia, kfdia
377    rmswind(i) = rmscon
378  END DO
379
380  IF (lozpr) THEN
381    DO i = kidia, kfdia
382      IF (zpr(i)>pcrit) THEN
383        rmswind(i) = rmscon + ((zpr(i)-pcrit)/zpr(i))*pcons
384      END IF
385    END DO
386  END IF
387
388  DO i = kidia, kfdia
389    IF (rmswind(i)>0.0) THEN
390      ilrms = ilrms + 1
391      lorms(i) = .TRUE.
392    END IF
393  END DO
394
395  ! * CALCULATE GWD (NOTE THAT DIFFUSION COEFFICIENT AND
396  ! * HEATING RATE ONLY CALCULATED IF IHEATCAL = 1).
397
398  IF (ilrms>0) THEN
399
400    CALL hines_extro0(drag_u, drag_v, heat, diffco, flux_u, flux_v, uhs, vhs, &
401      bvfreq, density, visc_mol, alt, rmswind, k_alpha, m_alpha, v_alpha, &
402      sigma_alpha, sigsqh_alpha, ak_alpha, mmin_alpha, i_alpha, sigma_t, &
403      densbot, bvfbot, 1, iheatcal, icutoff, iprint, nsmax, smco, alt_cutoff, &
404      kstar, slope, f1, f2, f3, f5, f6, naz, sigsqmcw, sigmatm, kidia, klon, &
405      1, levbot, klon, klev, nazmth, lorms, smoothr1, smoothr2, sigalpmc, &
406      f2mod)
407
408    ! * ADD ON HINES' GWD TENDENCIES TO OROGRAPHIC TENDENCIES AND
409    ! * APPLY HINES' GW DRAG ON (UROW,VROW) WORK ARRAYS.
410
411    DO l = 1, klev
412      DO i = kidia, kfdia
413        utendgw(i, l) = utendgw(i, l) + drag_u(i, l)
414        vtendgw(i, l) = vtendgw(i, l) + drag_v(i, l)
415      END DO
416    END DO
417
418
419    ! * END OF HINES CALCULATIONS.
420
421  END IF
422
423  ! -----------------------------------------------------------------------
424
425  DO jl = kidia, kfdia
426    zustrhi(jl) = flux_u(jl, 1)
427    zvstrhi(jl) = flux_v(jl, 1)
428    DO jk = 1, klev
429      le = klev - jk + 1
430      d_u_hin(jl, jk) = utendgw(jl, le)*dtime
431      d_v_hin(jl, jk) = vtendgw(jl, le)*dtime
432    END DO
433  END DO
434
435  ! PRINT *,'UTENDGW:',UTENDGW
436
437  ! PRINT *,' HINES HAS BEEN COMPLETED (LONG ISNT IT...)'
438
439  RETURN
440999 CONTINUE
441
442  ! * IF ERROR DETECTED THEN ABORT.
443
444  WRITE (nmessg, 6000)
445  WRITE (nmessg, 6010) ierror
4466000 FORMAT (/' EXECUTION ABORTED IN GWDOREXV')
4476010 FORMAT ('     ERROR FLAG =', I4)
448
449
450  RETURN
451END SUBROUTINE hines_gwd
452! /
453! /
454
455
456SUBROUTINE hines_extro0(drag_u, drag_v, heat, diffco, flux_u, flux_v, vel_u, &
457    vel_v, bvfreq, density, visc_mol, alt, rmswind, k_alpha, m_alpha, &
458    v_alpha, sigma_alpha, sigsqh_alpha, ak_alpha, mmin_alpha, i_alpha, &
459    sigma_t, densb, bvfb, iorder, iheatcal, icutoff, iprint, nsmax, smco, &
460    alt_cutoff, kstar, slope, f1, f2, f3, f5, f6, naz, sigsqmcw, sigmatm, &
461    il1, il2, lev1, lev2, nlons, nlevs, nazmth, lorms, smoothr1, smoothr2, &
462    sigalpmc, f2mod)
463
464  IMPLICIT NONE
465
466  ! Main routine for Hines' "extrowave" gravity wave parameterization based
467  ! on Hines' Doppler spread theory. This routine calculates zonal
468  ! and meridional components of gravity wave drag, heating rates
469  ! and diffusion coefficient on a longitude by altitude grid.
470  ! No "mythical" lower boundary region calculation is made so it
471  ! is assumed that lowest level winds are weak (i.e, approximately zero).
472
473  ! Aug. 13/95 - C. McLandress
474  ! SEPT. /95  - N.McFarlane
475
476  ! Modifications:
477
478  ! Output arguements:
479
480  ! * DRAG_U = zonal component of gravity wave drag (m/s^2).
481  ! * DRAG_V = meridional component of gravity wave drag (m/s^2).
482  ! * HEAT   = gravity wave heating (K/sec).
483  ! * DIFFCO = diffusion coefficient (m^2/sec)
484  ! * FLUX_U = zonal component of vertical momentum flux (Pascals)
485  ! * FLUX_V = meridional component of vertical momentum flux (Pascals)
486
487  ! Input arguements:
488
489  ! * VEL_U      = background zonal wind component (m/s).
490  ! * VEL_V      = background meridional wind component (m/s).
491  ! * BVFREQ     = background Brunt Vassala frequency (radians/sec).
492  ! * DENSITY    = background density (kg/m^3)
493  ! * VISC_MOL   = molecular viscosity (m^2/s)
494  ! * ALT        = altitude of momentum, density, buoyancy levels (m)
495  ! *              (NOTE: levels ordered so that ALT(I,1) > ALT(I,2), etc.)
496  ! * RMSWIND   = root mean square gravity wave wind at lowest level (m/s).
497  ! * K_ALPHA    = horizontal wavenumber of each azimuth (1/m).
498  ! * IORDER       = 1 means vertical levels are indexed from top down
499  ! *              (i.e., highest level indexed 1 and lowest level NLEVS);
500  ! *           .NE. 1 highest level is index NLEVS.
501  ! * IHEATCAL   = 1 to calculate heating rates and diffusion coefficient.
502  ! * IPRINT     = 1 to print out various arrays.
503  ! * ICUTOFF    = 1 to exponentially damp GWD, heating and diffusion
504  ! *              arrays above ALT_CUTOFF; otherwise arrays not modified.
505  ! * ALT_CUTOFF = altitude in meters above which exponential decay applied.
506  ! * SMCO       = smoothing factor used to smooth cutoff vertical
507  ! *              wavenumbers and total rms winds in vertical direction
508  ! *              before calculating drag or heating
509  ! *              (SMCO >= 1 ==> 1:SMCO:1 stencil used).
510  ! * NSMAX      = number of times smoother applied ( >= 1),
511  ! *            = 0 means no smoothing performed.
512  ! * KSTAR      = typical gravity wave horizontal wavenumber (1/m).
513  ! * SLOPE      = slope of incident vertical wavenumber spectrum
514  ! *              (SLOPE must equal 1., 1.5 or 2.).
515  ! * F1 to F6   = Hines's fudge factors (F4 not needed since used for
516  ! *              vertical flux of vertical momentum).
517  ! * NAZ        = actual number of horizontal azimuths used.
518  ! * IL1        = first longitudinal index to use (IL1 >= 1).
519  ! * IL2        = last longitudinal index to use (IL1 <= IL2 <= NLONS).
520  ! * LEV1       = index of first level for drag calculation.
521  ! * LEV2       = index of last level for drag calculation
522  ! *              (i.e., LEV1 < LEV2 <= NLEVS).
523  ! * NLONS      = number of longitudes.
524  ! * NLEVS      = number of vertical levels.
525  ! * NAZMTH     = azimuthal array dimension (NAZMTH >= NAZ).
526
527  ! Work arrays.
528
529  ! * M_ALPHA      = cutoff vertical wavenumber (1/m).
530  ! * V_ALPHA      = wind component at each azimuth (m/s) and if IHEATCAL=1
531  ! *                holds vertical derivative of cutoff wavenumber.
532  ! * SIGMA_ALPHA  = total rms wind in each azimuth (m/s).
533  ! * SIGSQH_ALPHA = portion of wind variance from waves having wave
534  ! *                normals in the alpha azimuth (m/s).
535  ! * SIGMA_T      = total rms horizontal wind (m/s).
536  ! * AK_ALPHA     = spectral amplitude factor at each azimuth
537  ! *                (i.e.,{AjKj}) in m^4/s^2.
538  ! * I_ALPHA      = Hines' integral.
539  ! * MMIN_ALPHA   = minimum value of cutoff wavenumber.
540  ! * DENSB        = background density at bottom level.
541  ! * BVFB         = buoyancy frequency at bottom level and
542  ! *                work array for ICUTOFF = 1.
543
544  ! * LORMS       = .TRUE. for drag computation
545
546  INTEGER naz, nlons, nlevs, nazmth, il1, il2, lev1, lev2
547  INTEGER icutoff, nsmax, iorder, iheatcal, iprint
548  REAL kstar(nlons), f1, f2, f3, f5, f6, slope
549  REAL alt_cutoff, smco
550  REAL drag_u(nlons, nlevs), drag_v(nlons, nlevs)
551  REAL heat(nlons, nlevs), diffco(nlons, nlevs)
552  REAL flux_u(nlons, nlevs), flux_v(nlons, nlevs)
553  REAL vel_u(nlons, nlevs), vel_v(nlons, nlevs)
554  REAL bvfreq(nlons, nlevs), density(nlons, nlevs)
555  REAL visc_mol(nlons, nlevs), alt(nlons, nlevs)
556  REAL rmswind(nlons), bvfb(nlons), densb(nlons)
557  REAL sigma_t(nlons, nlevs), sigsqmcw(nlons, nlevs, nazmth)
558  REAL sigma_alpha(nlons, nlevs, nazmth), sigmatm(nlons, nlevs)
559  REAL sigsqh_alpha(nlons, nlevs, nazmth)
560  REAL m_alpha(nlons, nlevs, nazmth), v_alpha(nlons, nlevs, nazmth)
561  REAL ak_alpha(nlons, nazmth), k_alpha(nlons, nazmth)
562  REAL mmin_alpha(nlons, nazmth), i_alpha(nlons, nazmth)
563  REAL smoothr1(nlons, nlevs), smoothr2(nlons, nlevs)
564  REAL sigalpmc(nlons, nlevs, nazmth)
565  REAL f2mod(nlons, nlevs)
566
567  LOGICAL lorms(nlons)
568
569  ! Internal variables.
570
571  INTEGER levbot, levtop, i, n, l, lev1p, lev2m
572  INTEGER ilprt1, ilprt2
573  ! -----------------------------------------------------------------------
574
575  ! PRINT *,' IN HINES_EXTRO0'
576  lev1p = lev1 + 1
577  lev2m = lev2 - 1
578
579  ! Index of lowest altitude level (bottom of drag calculation).
580
581  levbot = lev2
582  levtop = lev1
583  IF (iorder/=1) THEN
584    WRITE (6, 1)
5851   FORMAT (2X, ' error: IORDER NOT ONE! ')
586  END IF
587
588  ! Buoyancy and density at bottom level.
589
590  DO i = il1, il2
591    bvfb(i) = bvfreq(i, levbot)
592    densb(i) = density(i, levbot)
593  END DO
594
595  ! initialize some variables
596
597  DO n = 1, naz
598    DO l = lev1, lev2
599      DO i = il1, il2
600        m_alpha(i, l, n) = 0.0
601      END DO
602    END DO
603  END DO
604  DO l = lev1, lev2
605    DO i = il1, il2
606      sigma_t(i, l) = 0.0
607    END DO
608  END DO
609  DO n = 1, naz
610    DO i = il1, il2
611      i_alpha(i, n) = 0.0
612    END DO
613  END DO
614
615  ! Compute azimuthal wind components from zonal and meridional winds.
616
617  CALL hines_wind(v_alpha, vel_u, vel_v, naz, il1, il2, lev1, lev2, nlons, &
618    nlevs, nazmth)
619
620  ! Calculate cutoff vertical wavenumber and velocity variances.
621
622  CALL hines_wavnum(m_alpha, sigma_alpha, sigsqh_alpha, sigma_t, ak_alpha, &
623    v_alpha, visc_mol, density, densb, bvfreq, bvfb, rmswind, i_alpha, &
624    mmin_alpha, kstar, slope, f1, f2, f3, naz, levbot, levtop, il1, il2, &
625    nlons, nlevs, nazmth, sigsqmcw, sigmatm, lorms, sigalpmc, f2mod)
626  ! Smooth cutoff wavenumbers and total rms velocity in the vertical
627  ! direction NSMAX times, using FLUX_U as temporary work array.
628
629  IF (nsmax>0) THEN
630    DO n = 1, naz
631      DO l = lev1, lev2
632        DO i = il1, il2
633          smoothr1(i, l) = m_alpha(i, l, n)
634        END DO
635      END DO
636      CALL vert_smooth(smoothr1, smoothr2, smco, nsmax, il1, il2, lev1, lev2, &
637        nlons, nlevs)
638      DO l = lev1, lev2
639        DO i = il1, il2
640          m_alpha(i, l, n) = smoothr1(i, l)
641        END DO
642      END DO
643    END DO
644    CALL vert_smooth(sigma_t, smoothr2, smco, nsmax, il1, il2, lev1, lev2, &
645      nlons, nlevs)
646  END IF
647
648  ! Calculate zonal and meridional components of the
649  ! momentum flux and drag.
650
651  CALL hines_flux(flux_u, flux_v, drag_u, drag_v, alt, density, densb, &
652    m_alpha, ak_alpha, k_alpha, slope, naz, il1, il2, lev1, lev2, nlons, &
653    nlevs, nazmth, lorms)
654
655  ! Cutoff drag above ALT_CUTOFF, using BVFB as temporary work array.
656
657  IF (icutoff==1) THEN
658    CALL hines_exp(drag_u, bvfb, alt, alt_cutoff, iorder, il1, il2, lev1, &
659      lev2, nlons, nlevs)
660    CALL hines_exp(drag_v, bvfb, alt, alt_cutoff, iorder, il1, il2, lev1, &
661      lev2, nlons, nlevs)
662  END IF
663
664  ! Print out various arrays for diagnostic purposes.
665
666  IF (iprint==1) THEN
667    ilprt1 = 15
668    ilprt2 = 16
669    CALL hines_print(flux_u, flux_v, drag_u, drag_v, alt, sigma_t, &
670      sigma_alpha, v_alpha, m_alpha, 1, 1, 6, ilprt1, ilprt2, lev1, lev2, &
671      naz, nlons, nlevs, nazmth)
672  END IF
673
674  ! If not calculating heating rate and diffusion coefficient then finished.
675
676  IF (iheatcal/=1) RETURN
677
678  ! Calculate vertical derivative of cutoff wavenumber (store
679  ! in array V_ALPHA) using centered differences at interior gridpoints
680  ! and one-sided differences at first and last levels.
681
682  DO n = 1, naz
683    DO l = lev1p, lev2m
684      DO i = il1, il2
685        v_alpha(i, l, n) = (m_alpha(i,l+1,n)-m_alpha(i,l-1,n))/ &
686          (alt(i,l+1)-alt(i,l-1))
687      END DO
688    END DO
689    DO i = il1, il2
690      v_alpha(i, lev1, n) = (m_alpha(i,lev1p,n)-m_alpha(i,lev1,n))/ &
691        (alt(i,lev1p)-alt(i,lev1))
692    END DO
693    DO i = il1, il2
694      v_alpha(i, lev2, n) = (m_alpha(i,lev2,n)-m_alpha(i,lev2m,n))/ &
695        (alt(i,lev2)-alt(i,lev2m))
696    END DO
697  END DO
698
699  ! Heating rate and diffusion coefficient.
700
701  CALL hines_heat(heat, diffco, m_alpha, v_alpha, ak_alpha, k_alpha, bvfreq, &
702    density, densb, sigma_t, visc_mol, kstar, slope, f2, f3, f5, f6, naz, &
703    il1, il2, lev1, lev2, nlons, nlevs, nazmth)
704
705  ! Finished.
706
707  RETURN
708  ! -----------------------------------------------------------------------
709END SUBROUTINE hines_extro0
710
711SUBROUTINE hines_wavnum(m_alpha, sigma_alpha, sigsqh_alpha, sigma_t, &
712    ak_alpha, v_alpha, visc_mol, density, densb, bvfreq, bvfb, rms_wind, &
713    i_alpha, mmin_alpha, kstar, slope, f1, f2, f3, naz, levbot, levtop, il1, &
714    il2, nlons, nlevs, nazmth, sigsqmcw, sigmatm, lorms, sigalpmc, f2mod)
715  IMPLICIT NONE
716  ! This routine calculates the cutoff vertical wavenumber and velocity
717  ! variances on a longitude by altitude grid for the Hines' Doppler
718  ! spread gravity wave drag parameterization scheme.
719  ! NOTE: (1) only values of four or eight can be used for # azimuths (NAZ).
720  ! (2) only values of 1.0, 1.5 or 2.0 can be used for slope (SLOPE).
721
722  ! Aug. 10/95 - C. McLandress
723
724  ! Output arguements:
725
726  ! * M_ALPHA      = cutoff wavenumber at each azimuth (1/m).
727  ! * SIGMA_ALPHA  = total rms wind in each azimuth (m/s).
728  ! * SIGSQH_ALPHA = portion of wind variance from waves having wave
729  ! *                normals in the alpha azimuth (m/s).
730  ! * SIGMA_T      = total rms horizontal wind (m/s).
731  ! * AK_ALPHA     = spectral amplitude factor at each azimuth
732  ! *                (i.e.,{AjKj}) in m^4/s^2.
733
734  ! Input arguements:
735
736  ! * V_ALPHA  = wind component at each azimuth (m/s).
737  ! * VISC_MOL = molecular viscosity (m^2/s)
738  ! * DENSITY  = background density (kg/m^3).
739  ! * DENSB    = background density at model bottom (kg/m^3).
740  ! * BVFREQ   = background Brunt Vassala frequency (radians/sec).
741  ! * BVFB     = background Brunt Vassala frequency at model bottom.
742  ! * RMS_WIND = root mean square gravity wave wind at lowest level (m/s).
743  ! * KSTAR    = typical gravity wave horizontal wavenumber (1/m).
744  ! * SLOPE    = slope of incident vertical wavenumber spectrum
745  ! *            (SLOPE = 1., 1.5 or 2.).
746  ! * F1,F2,F3 = Hines's fudge factors.
747  ! * NAZ      = actual number of horizontal azimuths used (4 or 8).
748  ! * LEVBOT   = index of lowest vertical level.
749  ! * LEVTOP   = index of highest vertical level
750  ! *            (NOTE: if LEVTOP < LEVBOT then level index
751  ! *             increases from top down).
752  ! * IL1      = first longitudinal index to use (IL1 >= 1).
753  ! * IL2      = last longitudinal index to use (IL1 <= IL2 <= NLONS).
754  ! * NLONS    = number of longitudes.
755  ! * NLEVS    = number of vertical levels.
756  ! * NAZMTH   = azimuthal array dimension (NAZMTH >= NAZ).
757
758  ! * LORMS       = .TRUE. for drag computation
759
760  ! Input work arrays:
761
762  ! * I_ALPHA    = Hines' integral at a single level.
763  ! * MMIN_ALPHA = minimum value of cutoff wavenumber.
764
765  INTEGER naz, levbot, levtop, il1, il2, nlons, nlevs, nazmth
766  REAL slope, kstar(nlons), f1, f2, f3, f2mfac
767  REAL m_alpha(nlons, nlevs, nazmth)
768  REAL sigma_alpha(nlons, nlevs, nazmth)
769  REAL sigalpmc(nlons, nlevs, nazmth)
770  REAL sigsqh_alpha(nlons, nlevs, nazmth)
771  REAL sigsqmcw(nlons, nlevs, nazmth)
772  REAL sigma_t(nlons, nlevs)
773  REAL sigmatm(nlons, nlevs)
774  REAL ak_alpha(nlons, nazmth)
775  REAL v_alpha(nlons, nlevs, nazmth)
776  REAL visc_mol(nlons, nlevs)
777  REAL f2mod(nlons, nlevs)
778  REAL density(nlons, nlevs), densb(nlons)
779  REAL bvfreq(nlons, nlevs), bvfb(nlons), rms_wind(nlons)
780  REAL i_alpha(nlons, nazmth), mmin_alpha(nlons, nazmth)
781
782  LOGICAL lorms(nlons)
783
784  ! Internal variables.
785
786  INTEGER i, l, n, lstart, lend, lincr, lbelow
787  REAL m_sub_m_turb, m_sub_m_mol, m_trial
788  REAL visc, visc_min, azfac, sp1
789
790  ! c      REAL  N_OVER_M(1000), SIGFAC(1000)
791
792  REAL n_over_m(nlons), sigfac(nlons)
793  DATA visc_min/1.E-10/
794  ! -----------------------------------------------------------------------
795
796
797  ! PRINT *,'IN HINES_WAVNUM'
798  sp1 = slope + 1.
799
800  ! Indices of levels to process.
801
802  IF (levbot>levtop) THEN
803    lstart = levbot - 1
804    lend = levtop
805    lincr = -1
806  ELSE
807    WRITE (6, 1)
8081   FORMAT (2X, ' error: IORDER NOT ONE! ')
809  END IF
810
811  ! Use horizontal isotropy to calculate azimuthal variances at bottom level.
812
813  azfac = 1./real(naz)
814  DO n = 1, naz
815    DO i = il1, il2
816      sigsqh_alpha(i, levbot, n) = azfac*rms_wind(i)**2
817    END DO
818  END DO
819
820  ! Velocity variances at bottom level.
821
822  CALL hines_sigma(sigma_t, sigma_alpha, sigsqh_alpha, naz, levbot, il1, il2, &
823    nlons, nlevs, nazmth)
824
825  CALL hines_sigma(sigmatm, sigalpmc, sigsqmcw, naz, levbot, il1, il2, nlons, &
826    nlevs, nazmth)
827
828  ! Calculate cutoff wavenumber and spectral amplitude factor
829  ! at bottom level where it is assumed that background winds vanish
830  ! and also initialize minimum value of cutoff wavnumber.
831
832  DO n = 1, naz
833    DO i = il1, il2
834      IF (lorms(i)) THEN
835        m_alpha(i, levbot, n) = bvfb(i)/(f1*sigma_alpha(i,levbot,n)+f2* &
836          sigma_t(i,levbot))
837        ak_alpha(i, n) = sigsqh_alpha(i, levbot, n)/ &
838          (m_alpha(i,levbot,n)**sp1/sp1)
839        mmin_alpha(i, n) = m_alpha(i, levbot, n)
840      END IF
841    END DO
842  END DO
843
844  ! Calculate quantities from the bottom upwards,
845  ! starting one level above bottom.
846
847  DO l = lstart, lend, lincr
848
849    ! Level beneath present level.
850
851    lbelow = l - lincr
852
853    ! Calculate N/m_M where m_M is maximum permissible value of the vertical
854    ! wavenumber (i.e., m > m_M are obliterated) and N is buoyancy frequency.
855    ! m_M is taken as the smaller of the instability-induced
856    ! wavenumber (M_SUB_M_TURB) and that imposed by molecular viscosity
857    ! (M_SUB_M_MOL). Since variance at this level is not yet known
858    ! use value at level below.
859
860    DO i = il1, il2
861      IF (lorms(i)) THEN
862
863        f2mfac = sigmatm(i, lbelow)**2
864        f2mod(i, lbelow) = 1. + 2.*f2mfac/(f2mfac+sigma_t(i,lbelow)**2)
865
866        visc = amax1(visc_mol(i,l), visc_min)
867        m_sub_m_turb = bvfreq(i, l)/(f2*f2mod(i,lbelow)*sigma_t(i,lbelow))
868        m_sub_m_mol = (bvfreq(i,l)*kstar(i)/visc)**0.33333333/f3
869        IF (m_sub_m_turb<m_sub_m_mol) THEN
870          n_over_m(i) = f2*f2mod(i, lbelow)*sigma_t(i, lbelow)
871        ELSE
872          n_over_m(i) = bvfreq(i, l)/m_sub_m_mol
873        END IF
874      END IF
875    END DO
876
877    ! Calculate cutoff wavenumber at this level.
878
879    DO n = 1, naz
880      DO i = il1, il2
881        IF (lorms(i)) THEN
882
883          ! Calculate trial value (since variance at this level is not yet
884          ! known
885          ! use value at level below). If trial value is negative or if it
886          ! exceeds
887          ! minimum value (not permitted) then set it to minimum value.
888
889          m_trial = bvfb(i)/(f1*(sigma_alpha(i,lbelow,n)+sigalpmc(i,lbelow, &
890            n))+n_over_m(i)+v_alpha(i,l,n))
891          IF (m_trial<=0. .OR. m_trial>mmin_alpha(i,n)) THEN
892            m_trial = mmin_alpha(i, n)
893          END IF
894          m_alpha(i, l, n) = m_trial
895
896          ! Reset minimum value of cutoff wavenumber if necessary.
897
898          IF (m_alpha(i,l,n)<mmin_alpha(i,n)) THEN
899            mmin_alpha(i, n) = m_alpha(i, l, n)
900          END IF
901
902        END IF
903      END DO
904    END DO
905
906    ! Calculate the Hines integral at this level.
907
908    CALL hines_intgrl(i_alpha, v_alpha, m_alpha, bvfb, slope, naz, l, il1, &
909      il2, nlons, nlevs, nazmth, lorms)
910
911
912    ! Calculate the velocity variances at this level.
913
914    DO i = il1, il2
915      sigfac(i) = densb(i)/density(i, l)*bvfreq(i, l)/bvfb(i)
916    END DO
917    DO n = 1, naz
918      DO i = il1, il2
919        sigsqh_alpha(i, l, n) = sigfac(i)*ak_alpha(i, n)*i_alpha(i, n)
920      END DO
921    END DO
922    CALL hines_sigma(sigma_t, sigma_alpha, sigsqh_alpha, naz, l, il1, il2, &
923      nlons, nlevs, nazmth)
924
925    CALL hines_sigma(sigmatm, sigalpmc, sigsqmcw, naz, l, il1, il2, nlons, &
926      nlevs, nazmth)
927
928    ! End of level loop.
929
930  END DO
931
932  RETURN
933  ! -----------------------------------------------------------------------
934END SUBROUTINE hines_wavnum
935
936SUBROUTINE hines_wind(v_alpha, vel_u, vel_v, naz, il1, il2, lev1, lev2, &
937    nlons, nlevs, nazmth)
938    IMPLICIT NONE
939  ! This routine calculates the azimuthal horizontal background wind
940  ! components
941  ! on a longitude by altitude grid for the case of 4 or 8 azimuths for
942  ! the Hines' Doppler spread GWD parameterization scheme.
943
944  ! Aug. 7/95 - C. McLandress
945
946  ! Output arguement:
947
948  ! * V_ALPHA   = background wind component at each azimuth (m/s).
949  ! *             (note: first azimuth is in eastward direction
950  ! *              and rotate in counterclockwise direction.)
951
952  ! Input arguements:
953
954  ! * VEL_U     = background zonal wind component (m/s).
955  ! * VEL_V     = background meridional wind component (m/s).
956  ! * NAZ       = actual number of horizontal azimuths used (must be 4 or 8).
957  ! * IL1       = first longitudinal index to use (IL1 >= 1).
958  ! * IL2       = last longitudinal index to use (IL1 <= IL2 <= NLONS).
959  ! * LEV1      = first altitude level to use (LEV1 >=1).
960  ! * LEV2      = last altitude level to use (LEV1 < LEV2 <= NLEVS).
961  ! * NLONS     = number of longitudes.
962  ! * NLEVS     = number of vertical levels.
963  ! * NAZMTH    = azimuthal array dimension (NAZMTH >= NAZ).
964
965  ! Constants in DATA statements.
966
967  ! * COS45 = cosine of 45 degrees.
968  ! * UMIN  = minimum allowable value for zonal or meridional
969  ! *         wind component (m/s).
970
971  ! Subroutine arguements.
972
973  INTEGER naz, il1, il2, lev1, lev2
974  INTEGER nlons, nlevs, nazmth
975  REAL v_alpha(nlons, nlevs, nazmth)
976  REAL vel_u(nlons, nlevs), vel_v(nlons, nlevs)
977
978  ! Internal variables.
979
980  INTEGER i, l
981  REAL u, v, cos45, umin
982
983  DATA cos45/0.7071068/
984  DATA umin/0.001/
985  ! -----------------------------------------------------------------------
986
987  ! Case with 4 azimuths.
988
989
990  ! PRINT *,'IN HINES_WIND'
991  IF (naz==4) THEN
992    DO l = lev1, lev2
993      DO i = il1, il2
994        u = vel_u(i, l)
995        v = vel_v(i, l)
996        IF (abs(u)<umin) u = umin
997        IF (abs(v)<umin) v = umin
998        v_alpha(i, l, 1) = u
999        v_alpha(i, l, 2) = v
1000        v_alpha(i, l, 3) = -u
1001        v_alpha(i, l, 4) = -v
1002      END DO
1003    END DO
1004  END IF
1005
1006  ! Case with 8 azimuths.
1007
1008  IF (naz==8) THEN
1009    DO l = lev1, lev2
1010      DO i = il1, il2
1011        u = vel_u(i, l)
1012        v = vel_v(i, l)
1013        IF (abs(u)<umin) u = umin
1014        IF (abs(v)<umin) v = umin
1015        v_alpha(i, l, 1) = u
1016        v_alpha(i, l, 2) = cos45*(v+u)
1017        v_alpha(i, l, 3) = v
1018        v_alpha(i, l, 4) = cos45*(v-u)
1019        v_alpha(i, l, 5) = -u
1020        v_alpha(i, l, 6) = -v_alpha(i, l, 2)
1021        v_alpha(i, l, 7) = -v
1022        v_alpha(i, l, 8) = -v_alpha(i, l, 4)
1023      END DO
1024    END DO
1025  END IF
1026
1027  RETURN
1028  ! -----------------------------------------------------------------------
1029END SUBROUTINE hines_wind
1030
1031SUBROUTINE hines_flux(flux_u, flux_v, drag_u, drag_v, alt, density, densb, &
1032    m_alpha, ak_alpha, k_alpha, slope, naz, il1, il2, lev1, lev2, nlons, &
1033    nlevs, nazmth, lorms)
1034    IMPLICIT NONE
1035  ! Calculate zonal and meridional components of the vertical flux
1036  ! of horizontal momentum and corresponding wave drag (force per unit mass)
1037  ! on a longitude by altitude grid for the Hines' Doppler spread
1038  ! GWD parameterization scheme.
1039  ! NOTE: only 4 or 8 azimuths can be used.
1040
1041  ! Aug. 6/95 - C. McLandress
1042
1043  ! Output arguements:
1044
1045  ! * FLUX_U = zonal component of vertical momentum flux (Pascals)
1046  ! * FLUX_V = meridional component of vertical momentum flux (Pascals)
1047  ! * DRAG_U = zonal component of drag (m/s^2).
1048  ! * DRAG_V = meridional component of drag (m/s^2).
1049
1050  ! Input arguements:
1051
1052  ! * ALT       = altitudes (m).
1053  ! * DENSITY   = background density (kg/m^3).
1054  ! * DENSB     = background density at bottom level (kg/m^3).
1055  ! * M_ALPHA   = cutoff vertical wavenumber (1/m).
1056  ! * AK_ALPHA  = spectral amplitude factor (i.e., {AjKj} in m^4/s^2).
1057  ! * K_ALPHA   = horizontal wavenumber (1/m).
1058  ! * SLOPE     = slope of incident vertical wavenumber spectrum.
1059  ! * NAZ       = actual number of horizontal azimuths used (must be 4 or 8).
1060  ! * IL1       = first longitudinal index to use (IL1 >= 1).
1061  ! * IL2       = last longitudinal index to use (IL1 <= IL2 <= NLONS).
1062  ! * LEV1      = first altitude level to use (LEV1 >=1).
1063  ! * LEV2      = last altitude level to use (LEV1 < LEV2 <= NLEVS).
1064  ! * NLONS     = number of longitudes.
1065  ! * NLEVS     = number of vertical levels.
1066  ! * NAZMTH    = azimuthal array dimension (NAZMTH >= NAZ).
1067
1068  ! * LORMS       = .TRUE. for drag computation
1069
1070  ! Constant in DATA statement.
1071
1072  ! * COS45 = cosine of 45 degrees.
1073
1074  ! Subroutine arguements.
1075
1076  INTEGER naz, il1, il2, lev1, lev2
1077  INTEGER nlons, nlevs, nazmth
1078  REAL slope
1079  REAL flux_u(nlons, nlevs), flux_v(nlons, nlevs)
1080  REAL drag_u(nlons, nlevs), drag_v(nlons, nlevs)
1081  REAL alt(nlons, nlevs), density(nlons, nlevs), densb(nlons)
1082  REAL m_alpha(nlons, nlevs, nazmth)
1083  REAL ak_alpha(nlons, nazmth), k_alpha(nlons, nazmth)
1084
1085  LOGICAL lorms(nlons)
1086
1087  ! Internal variables.
1088
1089  INTEGER i, l, lev1p, lev2m, lev2p
1090  REAL cos45, prod2, prod4, prod6, prod8, dendz, dendz2
1091  DATA cos45/0.7071068/
1092  ! -----------------------------------------------------------------------
1093
1094  lev1p = lev1 + 1
1095  lev2m = lev2 - 1
1096  lev2p = lev2 + 1
1097
1098  ! Sum over azimuths for case where SLOPE = 1.
1099
1100  IF (slope==1.) THEN
1101
1102    ! Case with 4 azimuths.
1103
1104    IF (naz==4) THEN
1105      DO l = lev1, lev2
1106        DO i = il1, il2
1107          flux_u(i, l) = ak_alpha(i, 1)*k_alpha(i, 1)*m_alpha(i, l, 1) - &
1108            ak_alpha(i, 3)*k_alpha(i, 3)*m_alpha(i, l, 3)
1109          flux_v(i, l) = ak_alpha(i, 2)*k_alpha(i, 2)*m_alpha(i, l, 2) - &
1110            ak_alpha(i, 4)*k_alpha(i, 4)*m_alpha(i, l, 4)
1111        END DO
1112      END DO
1113    END IF
1114
1115    ! Case with 8 azimuths.
1116
1117    IF (naz==8) THEN
1118      DO l = lev1, lev2
1119        DO i = il1, il2
1120          prod2 = ak_alpha(i, 2)*k_alpha(i, 2)*m_alpha(i, l, 2)
1121          prod4 = ak_alpha(i, 4)*k_alpha(i, 4)*m_alpha(i, l, 4)
1122          prod6 = ak_alpha(i, 6)*k_alpha(i, 6)*m_alpha(i, l, 6)
1123          prod8 = ak_alpha(i, 8)*k_alpha(i, 8)*m_alpha(i, l, 8)
1124          flux_u(i, l) = ak_alpha(i, 1)*k_alpha(i, 1)*m_alpha(i, l, 1) - &
1125            ak_alpha(i, 5)*k_alpha(i, 5)*m_alpha(i, l, 5) + &
1126            cos45*(prod2-prod4-prod6+prod8)
1127          flux_v(i, l) = ak_alpha(i, 3)*k_alpha(i, 3)*m_alpha(i, l, 3) - &
1128            ak_alpha(i, 7)*k_alpha(i, 7)*m_alpha(i, l, 7) + &
1129            cos45*(prod2+prod4-prod6-prod8)
1130        END DO
1131      END DO
1132    END IF
1133
1134  END IF
1135
1136  ! Sum over azimuths for case where SLOPE not equal to 1.
1137
1138  IF (slope/=1.) THEN
1139
1140    ! Case with 4 azimuths.
1141
1142    IF (naz==4) THEN
1143      DO l = lev1, lev2
1144        DO i = il1, il2
1145          flux_u(i, l) = ak_alpha(i, 1)*k_alpha(i, 1)* &
1146            m_alpha(i, l, 1)**slope - ak_alpha(i, 3)*k_alpha(i, 3)*m_alpha(i, &
1147            l, 3)**slope
1148          flux_v(i, l) = ak_alpha(i, 2)*k_alpha(i, 2)* &
1149            m_alpha(i, l, 2)**slope - ak_alpha(i, 4)*k_alpha(i, 4)*m_alpha(i, &
1150            l, 4)**slope
1151        END DO
1152      END DO
1153    END IF
1154
1155    ! Case with 8 azimuths.
1156
1157    IF (naz==8) THEN
1158      DO l = lev1, lev2
1159        DO i = il1, il2
1160          prod2 = ak_alpha(i, 2)*k_alpha(i, 2)*m_alpha(i, l, 2)**slope
1161          prod4 = ak_alpha(i, 4)*k_alpha(i, 4)*m_alpha(i, l, 4)**slope
1162          prod6 = ak_alpha(i, 6)*k_alpha(i, 6)*m_alpha(i, l, 6)**slope
1163          prod8 = ak_alpha(i, 8)*k_alpha(i, 8)*m_alpha(i, l, 8)**slope
1164          flux_u(i, l) = ak_alpha(i, 1)*k_alpha(i, 1)* &
1165            m_alpha(i, l, 1)**slope - ak_alpha(i, 5)*k_alpha(i, 5)*m_alpha(i, &
1166            l, 5)**slope + cos45*(prod2-prod4-prod6+prod8)
1167          flux_v(i, l) = ak_alpha(i, 3)*k_alpha(i, 3)* &
1168            m_alpha(i, l, 3)**slope - ak_alpha(i, 7)*k_alpha(i, 7)*m_alpha(i, &
1169            l, 7)**slope + cos45*(prod2+prod4-prod6-prod8)
1170        END DO
1171      END DO
1172    END IF
1173
1174  END IF
1175
1176  ! Calculate flux from sum.
1177
1178  DO l = lev1, lev2
1179    DO i = il1, il2
1180      flux_u(i, l) = flux_u(i, l)*densb(i)/slope
1181      flux_v(i, l) = flux_v(i, l)*densb(i)/slope
1182    END DO
1183  END DO
1184
1185  ! Calculate drag at intermediate levels using centered differences
1186
1187  DO l = lev1p, lev2m
1188    DO i = il1, il2
1189      IF (lorms(i)) THEN
1190        ! cc       DENDZ2 = DENSITY(I,L) * ( ALT(I,L+1) - ALT(I,L-1) )
1191        dendz2 = density(i, l)*(alt(i,l-1)-alt(i,l))
1192        ! cc       DRAG_U(I,L) = - ( FLUX_U(I,L+1) - FLUX_U(I,L-1) ) / DENDZ2
1193        drag_u(i, l) = -(flux_u(i,l-1)-flux_u(i,l))/dendz2
1194        ! cc       DRAG_V(I,L) = - ( FLUX_V(I,L+1) - FLUX_V(I,L-1) ) / DENDZ2
1195        drag_v(i, l) = -(flux_v(i,l-1)-flux_v(i,l))/dendz2
1196
1197      END IF
1198    END DO
1199  END DO
1200
1201  ! Drag at first and last levels using one-side differences.
1202
1203  DO i = il1, il2
1204    IF (lorms(i)) THEN
1205      dendz = density(i, lev1)*(alt(i,lev1)-alt(i,lev1p))
1206      drag_u(i, lev1) = flux_u(i, lev1)/dendz
1207      drag_v(i, lev1) = flux_v(i, lev1)/dendz
1208    END IF
1209  END DO
1210  DO i = il1, il2
1211    IF (lorms(i)) THEN
1212      dendz = density(i, lev2)*(alt(i,lev2m)-alt(i,lev2))
1213      drag_u(i, lev2) = -(flux_u(i,lev2m)-flux_u(i,lev2))/dendz
1214      drag_v(i, lev2) = -(flux_v(i,lev2m)-flux_v(i,lev2))/dendz
1215    END IF
1216  END DO
1217  IF (nlevs>lev2) THEN
1218    DO i = il1, il2
1219      IF (lorms(i)) THEN
1220        dendz = density(i, lev2p)*(alt(i,lev2)-alt(i,lev2p))
1221        drag_u(i, lev2p) = -flux_u(i, lev2)/dendz
1222        drag_v(i, lev2p) = -flux_v(i, lev2)/dendz
1223      END IF
1224    END DO
1225  END IF
1226
1227  RETURN
1228  ! -----------------------------------------------------------------------
1229END SUBROUTINE hines_flux
1230
1231SUBROUTINE hines_heat(heat, diffco, m_alpha, dmdz_alpha, ak_alpha, k_alpha, &
1232    bvfreq, density, densb, sigma_t, visc_mol, kstar, slope, f2, f3, f5, f6, &
1233    naz, il1, il2, lev1, lev2, nlons, nlevs, nazmth)
1234  IMPLICIT NONE
1235  ! This routine calculates the gravity wave induced heating and
1236  ! diffusion coefficient on a longitude by altitude grid for
1237  ! the Hines' Doppler spread gravity wave drag parameterization scheme.
1238
1239  ! Aug. 6/95 - C. McLandress
1240
1241  ! Output arguements:
1242
1243  ! * HEAT   = gravity wave heating (K/sec).
1244  ! * DIFFCO = diffusion coefficient (m^2/sec)
1245
1246  ! Input arguements:
1247
1248  ! * M_ALPHA     = cutoff vertical wavenumber (1/m).
1249  ! * DMDZ_ALPHA  = vertical derivative of cutoff wavenumber.
1250  ! * AK_ALPHA    = spectral amplitude factor of each azimuth
1251  ! (i.e., {AjKj} in m^4/s^2).
1252  ! * K_ALPHA     = horizontal wavenumber of each azimuth (1/m).
1253  ! * BVFREQ      = background Brunt Vassala frequency (rad/sec).
1254  ! * DENSITY     = background density (kg/m^3).
1255  ! * DENSB       = background density at bottom level (kg/m^3).
1256  ! * SIGMA_T     = total rms horizontal wind (m/s).
1257  ! * VISC_MOL    = molecular viscosity (m^2/s).
1258  ! * KSTAR       = typical gravity wave horizontal wavenumber (1/m).
1259  ! * SLOPE       = slope of incident vertical wavenumber spectrum.
1260  ! * F2,F3,F5,F6 = Hines's fudge factors.
1261  ! * NAZ         = actual number of horizontal azimuths used.
1262  ! * IL1         = first longitudinal index to use (IL1 >= 1).
1263  ! * IL2         = last longitudinal index to use (IL1 <= IL2 <= NLONS).
1264  ! * LEV1        = first altitude level to use (LEV1 >=1).
1265  ! * LEV2        = last altitude level to use (LEV1 < LEV2 <= NLEVS).
1266  ! * NLONS       = number of longitudes.
1267  ! * NLEVS       = number of vertical levels.
1268  ! * NAZMTH      = azimuthal array dimension (NAZMTH >= NAZ).
1269
1270  INTEGER naz, il1, il2, lev1, lev2, nlons, nlevs, nazmth
1271  REAL kstar(nlons), slope, f2, f3, f5, f6
1272  REAL heat(nlons, nlevs), diffco(nlons, nlevs)
1273  REAL m_alpha(nlons, nlevs, nazmth), dmdz_alpha(nlons, nlevs, nazmth)
1274  REAL ak_alpha(nlons, nazmth), k_alpha(nlons, nazmth)
1275  REAL bvfreq(nlons, nlevs), density(nlons, nlevs), densb(nlons)
1276  REAL sigma_t(nlons, nlevs), visc_mol(nlons, nlevs)
1277
1278  ! Internal variables.
1279
1280  INTEGER i, l, n
1281  REAL m_sub_m_turb, m_sub_m_mol, m_sub_m, heatng
1282  REAL visc, visc_min, cpgas, sm1
1283
1284  ! specific heat at constant pressure
1285
1286  DATA cpgas/1004./
1287
1288  ! minimum permissible viscosity
1289
1290  DATA visc_min/1.E-10/
1291  ! -----------------------------------------------------------------------
1292
1293  ! Initialize heating array.
1294
1295  DO l = 1, nlevs
1296    DO i = 1, nlons
1297      heat(i, l) = 0.
1298    END DO
1299  END DO
1300
1301  ! Perform sum over azimuths for case where SLOPE = 1.
1302
1303  IF (slope==1.) THEN
1304    DO n = 1, naz
1305      DO l = lev1, lev2
1306        DO i = il1, il2
1307          heat(i, l) = heat(i, l) + ak_alpha(i, n)*k_alpha(i, n)*dmdz_alpha(i &
1308            , l, n)
1309        END DO
1310      END DO
1311    END DO
1312  END IF
1313
1314  ! Perform sum over azimuths for case where SLOPE not 1.
1315
1316  IF (slope/=1.) THEN
1317    sm1 = slope - 1.
1318    DO n = 1, naz
1319      DO l = lev1, lev2
1320        DO i = il1, il2
1321          heat(i, l) = heat(i, l) + ak_alpha(i, n)*k_alpha(i, n)*m_alpha(i, l &
1322            , n)**sm1*dmdz_alpha(i, l, n)
1323        END DO
1324      END DO
1325    END DO
1326  END IF
1327
1328  ! Heating and diffusion.
1329
1330  DO l = lev1, lev2
1331    DO i = il1, il2
1332
1333      ! Maximum permissible value of cutoff wavenumber is the smaller
1334      ! of the instability-induced wavenumber (M_SUB_M_TURB) and
1335      ! that imposed by molecular viscosity (M_SUB_M_MOL).
1336
1337      visc = amax1(visc_mol(i,l), visc_min)
1338      m_sub_m_turb = bvfreq(i, l)/(f2*sigma_t(i,l))
1339      m_sub_m_mol = (bvfreq(i,l)*kstar(i)/visc)**0.33333333/f3
1340      m_sub_m = amin1(m_sub_m_turb, m_sub_m_mol)
1341
1342      heatng = -heat(i, l)*f5*bvfreq(i, l)/m_sub_m*densb(i)/density(i, l)
1343      diffco(i, l) = f6*heatng**0.33333333/m_sub_m**1.33333333
1344      heat(i, l) = heatng/cpgas
1345
1346    END DO
1347  END DO
1348
1349  RETURN
1350  ! -----------------------------------------------------------------------
1351END SUBROUTINE hines_heat
1352
1353SUBROUTINE hines_sigma(sigma_t, sigma_alpha, sigsqh_alpha, naz, lev, il1, &
1354    il2, nlons, nlevs, nazmth)
1355  IMPLICIT NONE
1356  ! This routine calculates the total rms and azimuthal rms horizontal
1357  ! velocities at a given level on a longitude by altitude grid for
1358  ! the Hines' Doppler spread GWD parameterization scheme.
1359  ! NOTE: only four or eight azimuths can be used.
1360
1361  ! Aug. 7/95 - C. McLandress
1362
1363  ! Output arguements:
1364
1365  ! * SIGMA_T      = total rms horizontal wind (m/s).
1366  ! * SIGMA_ALPHA  = total rms wind in each azimuth (m/s).
1367
1368  ! Input arguements:
1369
1370  ! * SIGSQH_ALPHA = portion of wind variance from waves having wave
1371  ! *                normals in the alpha azimuth (m/s).
1372  ! * NAZ       = actual number of horizontal azimuths used (must be 4 or 8).
1373  ! * LEV       = altitude level to process.
1374  ! * IL1       = first longitudinal index to use (IL1 >= 1).
1375  ! * IL2       = last longitudinal index to use (IL1 <= IL2 <= NLONS).
1376  ! * NLONS     = number of longitudes.
1377  ! * NLEVS     = number of vertical levels.
1378  ! * NAZMTH    = azimuthal array dimension (NAZMTH >= NAZ).
1379
1380  ! Subroutine arguements.
1381
1382  INTEGER lev, naz, il1, il2
1383  INTEGER nlons, nlevs, nazmth
1384  REAL sigma_t(nlons, nlevs)
1385  REAL sigma_alpha(nlons, nlevs, nazmth)
1386  REAL sigsqh_alpha(nlons, nlevs, nazmth)
1387
1388  ! Internal variables.
1389
1390  INTEGER i, n
1391  REAL sum_even, sum_odd
1392  ! -----------------------------------------------------------------------
1393
1394  ! Calculate azimuthal rms velocity for the 4 azimuth case.
1395
1396  IF (naz==4) THEN
1397    DO i = il1, il2
1398      sigma_alpha(i, lev, 1) = sqrt(sigsqh_alpha(i,lev,1)+sigsqh_alpha(i,lev, &
1399        3))
1400      sigma_alpha(i, lev, 2) = sqrt(sigsqh_alpha(i,lev,2)+sigsqh_alpha(i,lev, &
1401        4))
1402      sigma_alpha(i, lev, 3) = sigma_alpha(i, lev, 1)
1403      sigma_alpha(i, lev, 4) = sigma_alpha(i, lev, 2)
1404    END DO
1405  END IF
1406
1407  ! Calculate azimuthal rms velocity for the 8 azimuth case.
1408
1409  IF (naz==8) THEN
1410    DO i = il1, il2
1411      sum_odd = (sigsqh_alpha(i,lev,1)+sigsqh_alpha(i,lev,3)+ &
1412        sigsqh_alpha(i,lev,5)+sigsqh_alpha(i,lev,7))/2.
1413      sum_even = (sigsqh_alpha(i,lev,2)+sigsqh_alpha(i,lev,4)+ &
1414        sigsqh_alpha(i,lev,6)+sigsqh_alpha(i,lev,8))/2.
1415      sigma_alpha(i, lev, 1) = sqrt(sigsqh_alpha(i,lev,1)+sigsqh_alpha(i,lev, &
1416        5)+sum_even)
1417      sigma_alpha(i, lev, 2) = sqrt(sigsqh_alpha(i,lev,2)+sigsqh_alpha(i,lev, &
1418        6)+sum_odd)
1419      sigma_alpha(i, lev, 3) = sqrt(sigsqh_alpha(i,lev,3)+sigsqh_alpha(i,lev, &
1420        7)+sum_even)
1421      sigma_alpha(i, lev, 4) = sqrt(sigsqh_alpha(i,lev,4)+sigsqh_alpha(i,lev, &
1422        8)+sum_odd)
1423      sigma_alpha(i, lev, 5) = sigma_alpha(i, lev, 1)
1424      sigma_alpha(i, lev, 6) = sigma_alpha(i, lev, 2)
1425      sigma_alpha(i, lev, 7) = sigma_alpha(i, lev, 3)
1426      sigma_alpha(i, lev, 8) = sigma_alpha(i, lev, 4)
1427    END DO
1428  END IF
1429
1430  ! Calculate total rms velocity.
1431
1432  DO i = il1, il2
1433    sigma_t(i, lev) = 0.
1434  END DO
1435  DO n = 1, naz
1436    DO i = il1, il2
1437      sigma_t(i, lev) = sigma_t(i, lev) + sigsqh_alpha(i, lev, n)
1438    END DO
1439  END DO
1440  DO i = il1, il2
1441    sigma_t(i, lev) = sqrt(sigma_t(i,lev))
1442  END DO
1443
1444  RETURN
1445  ! -----------------------------------------------------------------------
1446END SUBROUTINE hines_sigma
1447
1448SUBROUTINE hines_intgrl(i_alpha, v_alpha, m_alpha, bvfb, slope, naz, lev, &
1449    il1, il2, nlons, nlevs, nazmth, lorms)
1450  IMPLICIT NONE
1451  ! This routine calculates the vertical wavenumber integral
1452  ! for a single vertical level at each azimuth on a longitude grid
1453  ! for the Hines' Doppler spread GWD parameterization scheme.
1454  ! NOTE: (1) only spectral slopes of 1, 1.5 or 2 are permitted.
1455  ! (2) the integral is written in terms of the product QM
1456  ! which by construction is always less than 1. Series
1457  ! solutions are used for small |QM| and analytical solutions
1458  ! for remaining values.
1459
1460  ! Aug. 8/95 - C. McLandress
1461
1462  ! Output arguement:
1463
1464  ! * I_ALPHA = Hines' integral.
1465
1466  ! Input arguements:
1467
1468  ! * V_ALPHA = azimuthal wind component (m/s).
1469  ! * M_ALPHA = azimuthal cutoff vertical wavenumber (1/m).
1470  ! * BVFB    = background Brunt Vassala frequency at model bottom.
1471  ! * SLOPE   = slope of initial vertical wavenumber spectrum
1472  ! *           (must use SLOPE = 1., 1.5 or 2.)
1473  ! * NAZ     = actual number of horizontal azimuths used.
1474  ! * LEV     = altitude level to process.
1475  ! * IL1     = first longitudinal index to use (IL1 >= 1).
1476  ! * IL2     = last longitudinal index to use (IL1 <= IL2 <= NLONS).
1477  ! * NLONS   = number of longitudes.
1478  ! * NLEVS   = number of vertical levels.
1479  ! * NAZMTH  = azimuthal array dimension (NAZMTH >= NAZ).
1480
1481  ! * LORMS       = .TRUE. for drag computation
1482
1483  ! Constants in DATA statements:
1484
1485  ! * QMIN = minimum value of Q_ALPHA (avoids indeterminant form of integral)
1486  ! * QM_MIN = minimum value of Q_ALPHA * M_ALPHA (used to avoid numerical
1487  ! *          problems).
1488
1489  INTEGER lev, naz, il1, il2, nlons, nlevs, nazmth
1490  REAL i_alpha(nlons, nazmth)
1491  REAL v_alpha(nlons, nlevs, nazmth)
1492  REAL m_alpha(nlons, nlevs, nazmth)
1493  REAL bvfb(nlons), slope
1494
1495  LOGICAL lorms(nlons)
1496
1497  ! Internal variables.
1498
1499  INTEGER i, n
1500  REAL q_alpha, qm, sqrtqm, q_min, qm_min
1501
1502  DATA q_min/1.0/, qm_min/0.01/
1503  ! -----------------------------------------------------------------------
1504
1505  ! For integer value SLOPE = 1.
1506
1507  IF (slope==1.) THEN
1508
1509    DO n = 1, naz
1510      DO i = il1, il2
1511        IF (lorms(i)) THEN
1512
1513          q_alpha = v_alpha(i, lev, n)/bvfb(i)
1514          qm = q_alpha*m_alpha(i, lev, n)
1515
1516          ! If |QM| is small then use first 4 terms series of Taylor series
1517          ! expansion of integral in order to avoid indeterminate form of
1518          ! integral,
1519          ! otherwise use analytical form of integral.
1520
1521          IF (abs(q_alpha)<q_min .OR. abs(qm)<qm_min) THEN
1522            IF (q_alpha==0.) THEN
1523              i_alpha(i, n) = m_alpha(i, lev, n)**2/2.
1524            ELSE
1525              i_alpha(i, n) = (qm**2/2.+qm**3/3.+qm**4/4.+qm**5/5.)/ &
1526                q_alpha**2
1527            END IF
1528          ELSE
1529            i_alpha(i, n) = -(alog(1.-qm)+qm)/q_alpha**2
1530          END IF
1531
1532        END IF
1533      END DO
1534    END DO
1535
1536  END IF
1537
1538  ! For integer value SLOPE = 2.
1539
1540  IF (slope==2.) THEN
1541
1542    DO n = 1, naz
1543      DO i = il1, il2
1544        IF (lorms(i)) THEN
1545
1546          q_alpha = v_alpha(i, lev, n)/bvfb(i)
1547          qm = q_alpha*m_alpha(i, lev, n)
1548
1549          ! If |QM| is small then use first 4 terms series of Taylor series
1550          ! expansion of integral in order to avoid indeterminate form of
1551          ! integral,
1552          ! otherwise use analytical form of integral.
1553
1554          IF (abs(q_alpha)<q_min .OR. abs(qm)<qm_min) THEN
1555            IF (q_alpha==0.) THEN
1556              i_alpha(i, n) = m_alpha(i, lev, n)**3/3.
1557            ELSE
1558              i_alpha(i, n) = (qm**3/3.+qm**4/4.+qm**5/5.+qm**6/6.)/ &
1559                q_alpha**3
1560            END IF
1561          ELSE
1562            i_alpha(i, n) = -(alog(1.-qm)+qm+qm**2/2.)/q_alpha**3
1563          END IF
1564
1565        END IF
1566      END DO
1567    END DO
1568
1569  END IF
1570
1571  ! For real value SLOPE = 1.5
1572
1573  IF (slope==1.5) THEN
1574
1575    DO n = 1, naz
1576      DO i = il1, il2
1577        IF (lorms(i)) THEN
1578
1579          q_alpha = v_alpha(i, lev, n)/bvfb(i)
1580          qm = q_alpha*m_alpha(i, lev, n)
1581
1582          ! If |QM| is small then use first 4 terms series of Taylor series
1583          ! expansion of integral in order to avoid indeterminate form of
1584          ! integral,
1585          ! otherwise use analytical form of integral.
1586
1587          IF (abs(q_alpha)<q_min .OR. abs(qm)<qm_min) THEN
1588            IF (q_alpha==0.) THEN
1589              i_alpha(i, n) = m_alpha(i, lev, n)**2.5/2.5
1590            ELSE
1591              i_alpha(i, n) = (qm/2.5+qm**2/3.5+qm**3/4.5+qm**4/5.5)* &
1592                m_alpha(i, lev, n)**1.5/q_alpha
1593            END IF
1594          ELSE
1595            qm = abs(qm)
1596            sqrtqm = sqrt(qm)
1597            IF (q_alpha>=0.) THEN
1598              i_alpha(i, n) = (alog((1.+sqrtqm)/(1.-sqrtqm))-2.*sqrtqm*(1.+qm &
1599                /3.))/q_alpha**2.5
1600            ELSE
1601              i_alpha(i, n) = 2.*(atan(sqrtqm)+sqrtqm*(qm/3.-1.))/ &
1602                abs(q_alpha)**2.5
1603            END IF
1604          END IF
1605
1606        END IF
1607      END DO
1608    END DO
1609
1610  END IF
1611
1612  ! If integral is negative (which in principal should not happen) then
1613  ! print a message and some info since execution will abort when calculating
1614  ! the variances.
1615
1616  ! DO 80 N = 1,NAZ
1617  ! DO 70 I = IL1,IL2
1618  ! IF (I_ALPHA(I,N).LT.0.)  THEN
1619  ! WRITE (6,*)
1620  ! WRITE (6,*) '******************************'
1621  ! WRITE (6,*) 'Hines integral I_ALPHA < 0 '
1622  ! WRITE (6,*) '  longitude I=',I
1623  ! WRITE (6,*) '  azimuth   N=',N
1624  ! WRITE (6,*) '  level   LEV=',LEV
1625  ! WRITE (6,*) '  I_ALPHA =',I_ALPHA(I,N)
1626  ! WRITE (6,*) '  V_ALPHA =',V_ALPHA(I,LEV,N)
1627  ! WRITE (6,*) '  M_ALPHA =',M_ALPHA(I,LEV,N)
1628  ! WRITE (6,*) '  Q_ALPHA =',V_ALPHA(I,LEV,N) / BVFB(I)
1629  ! WRITE (6,*) '  QM      =',V_ALPHA(I,LEV,N) / BVFB(I)
1630  ! ^                                * M_ALPHA(I,LEV,N)
1631  ! WRITE (6,*) '******************************'
1632  ! END IF
1633  ! 70     CONTINUE
1634  ! 80   CONTINUE
1635
1636  RETURN
1637  ! -----------------------------------------------------------------------
1638END SUBROUTINE hines_intgrl
1639
1640SUBROUTINE hines_setup(naz, slope, f1, f2, f3, f5, f6, kstar, icutoff, &
1641    alt_cutoff, smco, nsmax, iheatcal, k_alpha, ierror, nmessg, nlons, &
1642    nazmth, coslat)
1643  IMPLICIT NONE
1644  ! This routine specifies various parameters needed for the
1645  ! the Hines' Doppler spread gravity wave drag parameterization scheme.
1646
1647  ! Aug. 8/95 - C. McLandress
1648
1649  ! Output arguements:
1650
1651  ! * NAZ        = actual number of horizontal azimuths used
1652  ! *              (code set up presently for only NAZ = 4 or 8).
1653  ! * SLOPE      = slope of incident vertical wavenumber spectrum
1654  ! *              (code set up presently for SLOPE 1., 1.5 or 2.).
1655  ! * F1         = "fudge factor" used in calculation of trial value of
1656  ! *              azimuthal cutoff wavenumber M_ALPHA (1.2 <= F1 <= 1.9).
1657  ! * F2         = "fudge factor" used in calculation of maximum
1658  ! *              permissible instabiliy-induced cutoff wavenumber
1659  ! *              M_SUB_M_TURB (0.1 <= F2 <= 1.4).
1660  ! * F3         = "fudge factor" used in calculation of maximum
1661  ! *              permissible molecular viscosity-induced cutoff wavenumber
1662  ! *              M_SUB_M_MOL (0.1 <= F2 <= 1.4).
1663  ! * F5         = "fudge factor" used in calculation of heating rate
1664  ! *              (1 <= F5 <= 3).
1665  ! * F6         = "fudge factor" used in calculation of turbulent
1666  ! *              diffusivity coefficient.
1667  ! * KSTAR      = typical gravity wave horizontal wavenumber (1/m)
1668  ! *              used in calculation of M_SUB_M_TURB.
1669  ! * ICUTOFF    = 1 to exponentially damp off GWD, heating and diffusion
1670  ! *              arrays above ALT_CUTOFF; otherwise arrays not modified.
1671  ! * ALT_CUTOFF = altitude in meters above which exponential decay applied.
1672  ! * SMCO       = smoother used to smooth cutoff vertical wavenumbers
1673  ! *              and total rms winds before calculating drag or heating.
1674  ! *              (==> a 1:SMCO:1 stencil used; SMCO >= 1.).
1675  ! * NSMAX      = number of times smoother applied ( >= 1),
1676  ! *            = 0 means no smoothing performed.
1677  ! * IHEATCAL   = 1 to calculate heating rates and diffusion coefficient.
1678  ! *            = 0 means only drag and flux calculated.
1679  ! * K_ALPHA    = horizontal wavenumber of each azimuth (1/m) which
1680  ! *              is set here to KSTAR.
1681  ! * IERROR     = error flag.
1682  ! *            = 0 no errors.
1683  ! *            = 10 ==> NAZ > NAZMTH
1684  ! *            = 20 ==> invalid number of azimuths (NAZ must be 4 or 8).
1685  ! *            = 30 ==> invalid slope (SLOPE must be 1., 1.5 or 2.).
1686  ! *            = 40 ==> invalid smoother (SMCO must be >= 1.)
1687
1688  ! Input arguements:
1689
1690  ! * NMESSG  = output unit number where messages to be printed.
1691  ! * NLONS   = number of longitudes.
1692  ! * NAZMTH  = azimuthal array dimension (NAZMTH >= NAZ).
1693
1694  INTEGER naz, nlons, nazmth, iheatcal, icutoff
1695  INTEGER nmessg, nsmax, ierror
1696  REAL kstar(nlons), slope, f1, f2, f3, f5, f6, alt_cutoff, smco
1697  REAL k_alpha(nlons, nazmth), coslat(nlons)
1698  REAL ksmin, ksmax
1699
1700  ! Internal variables.
1701
1702  INTEGER i, n
1703  ! -----------------------------------------------------------------------
1704
1705  ! Specify constants.
1706
1707  naz = 8
1708  slope = 1.
1709  f1 = 1.5
1710  f2 = 0.3
1711  f3 = 1.0
1712  f5 = 3.0
1713  f6 = 1.0
1714  ksmin = 1.E-5
1715  ksmax = 1.E-4
1716  DO i = 1, nlons
1717    kstar(i) = ksmin/(coslat(i)+(ksmin/ksmax))
1718  END DO
1719  icutoff = 1
1720  alt_cutoff = 105.E3
1721  smco = 2.0
1722  ! SMCO       = 1.0
1723  nsmax = 5
1724  ! NSMAX      = 2
1725  iheatcal = 0
1726
1727  ! Print information to output file.
1728
1729  ! WRITE (NMESSG,6000)
1730  ! 6000 FORMAT (/' Subroutine HINES_SETUP:')
1731  ! WRITE (NMESSG,*)  '  SLOPE = ', SLOPE
1732  ! WRITE (NMESSG,*)  '  NAZ = ', NAZ
1733  ! WRITE (NMESSG,*)  '  F1,F2,F3  = ', F1, F2, F3
1734  ! WRITE (NMESSG,*)  '  F5,F6     = ', F5, F6
1735  ! WRITE (NMESSG,*)  '  KSTAR     = ', KSTAR
1736  ! >           ,'  COSLAT     = ', COSLAT
1737  ! IF (ICUTOFF .EQ. 1)  THEN
1738  ! WRITE (NMESSG,*) '  Drag exponentially damped above ',
1739  ! &                       ALT_CUTOFF/1.E3
1740  ! END IF
1741  ! IF (NSMAX.LT.1 )  THEN
1742  ! WRITE (NMESSG,*) '  No smoothing of cutoff wavenumbers, etc'
1743  ! ELSE
1744  ! WRITE (NMESSG,*) '  Cutoff wavenumbers and sig_t smoothed:'
1745  ! WRITE (NMESSG,*) '    SMCO  =', SMCO
1746  ! WRITE (NMESSG,*) '    NSMAX =', NSMAX
1747  ! END IF
1748
1749  ! Check that things are setup correctly and log error if not
1750
1751  ierror = 0
1752  IF (naz>nazmth) ierror = 10
1753  IF (naz/=4 .AND. naz/=8) ierror = 20
1754  IF (slope/=1. .AND. slope/=1.5 .AND. slope/=2.) ierror = 30
1755  IF (smco<1.) ierror = 40
1756
1757  ! Use single value for azimuthal-dependent horizontal wavenumber.
1758
1759  DO n = 1, naz
1760    DO i = 1, nlons
1761      k_alpha(i, n) = kstar(i)
1762    END DO
1763  END DO
1764
1765  RETURN
1766  ! -----------------------------------------------------------------------
1767END SUBROUTINE hines_setup
1768
1769SUBROUTINE hines_print(flux_u, flux_v, drag_u, drag_v, alt, sigma_t, &
1770    sigma_alpha, v_alpha, m_alpha, iu_print, iv_print, nmessg, ilprt1, &
1771    ilprt2, levprt1, levprt2, naz, nlons, nlevs, nazmth)
1772  IMPLICIT NONE
1773  ! Print out altitude profiles of various quantities from
1774  ! Hines' Doppler spread gravity wave drag parameterization scheme.
1775  ! (NOTE: only for NAZ = 4 or 8).
1776
1777  ! Aug. 8/95 - C. McLandress
1778
1779  ! Input arguements:
1780
1781  ! * IU_PRINT = 1 to print out values in east-west direction.
1782  ! * IV_PRINT = 1 to print out values in north-south direction.
1783  ! * NMESSG   = unit number for printed output.
1784  ! * ILPRT1   = first longitudinal index to print.
1785  ! * ILPRT2   = last longitudinal index to print.
1786  ! * LEVPRT1  = first altitude level to print.
1787  ! * LEVPRT2  = last altitude level to print.
1788
1789  INTEGER naz, ilprt1, ilprt2, levprt1, levprt2
1790  INTEGER nlons, nlevs, nazmth
1791  INTEGER iu_print, iv_print, nmessg
1792  REAL flux_u(nlons, nlevs), flux_v(nlons, nlevs)
1793  REAL drag_u(nlons, nlevs), drag_v(nlons, nlevs)
1794  REAL alt(nlons, nlevs), sigma_t(nlons, nlevs)
1795  REAL sigma_alpha(nlons, nlevs, nazmth)
1796  REAL v_alpha(nlons, nlevs, nazmth), m_alpha(nlons, nlevs, nazmth)
1797
1798  ! Internal variables.
1799
1800  INTEGER n_east, n_west, n_north, n_south
1801  INTEGER i, l
1802  ! -----------------------------------------------------------------------
1803
1804  ! Azimuthal indices of cardinal directions.
1805
1806  n_east = 1
1807  IF (naz==4) THEN
1808    n_west = 3
1809    n_north = 2
1810    n_south = 4
1811  ELSE IF (naz==8) THEN
1812    n_west = 5
1813    n_north = 3
1814    n_south = 7
1815  END IF
1816
1817  ! Print out values for range of longitudes.
1818
1819  DO i = ilprt1, ilprt2
1820
1821    ! Print east-west wind, sigmas, cutoff wavenumbers, flux and drag.
1822
1823    IF (iu_print==1) THEN
1824      WRITE (nmessg, *)
1825      WRITE (nmessg, 6001) i
1826      WRITE (nmessg, 6005)
18276001  FORMAT ('Hines GW (east-west) at longitude I =', I3)
18286005  FORMAT (15X, ' U ', 2X, 'sig_E', 2X, 'sig_T', 3X, 'm_E', 4X, 'm_W', 4X, &
1829        'fluxU', 5X, 'gwdU')
1830      DO l = levprt1, levprt2
1831        WRITE (nmessg, 6701) alt(i, l)/1.E3, v_alpha(i, l, n_east), &
1832          sigma_alpha(i, l, n_east), sigma_t(i, l), &
1833          m_alpha(i, l, n_east)*1.E3, m_alpha(i, l, n_west)*1.E3, &
1834          flux_u(i, l)*1.E5, drag_u(i, l)*24.*3600.
1835      END DO
18366701  FORMAT (' z=', F7.2, 1X, 3F7.1, 2F7.3, F9.4, F9.3)
1837    END IF
1838
1839    ! Print north-south winds, sigmas, cutoff wavenumbers, flux and drag.
1840
1841    IF (iv_print==1) THEN
1842      WRITE (nmessg, *)
1843      WRITE (nmessg, 6002) i
18446002  FORMAT ('Hines GW (north-south) at longitude I =', I3)
1845      WRITE (nmessg, 6006)
18466006  FORMAT (15X, ' V ', 2X, 'sig_N', 2X, 'sig_T', 3X, 'm_N', 4X, 'm_S', 4X, &
1847        'fluxV', 5X, 'gwdV')
1848      DO l = levprt1, levprt2
1849        WRITE (nmessg, 6701) alt(i, l)/1.E3, v_alpha(i, l, n_north), &
1850          sigma_alpha(i, l, n_north), sigma_t(i, l), &
1851          m_alpha(i, l, n_north)*1.E3, m_alpha(i, l, n_south)*1.E3, &
1852          flux_v(i, l)*1.E5, drag_v(i, l)*24.*3600.
1853      END DO
1854    END IF
1855
1856  END DO
1857
1858  RETURN
1859  ! -----------------------------------------------------------------------
1860END SUBROUTINE hines_print
1861
1862SUBROUTINE hines_exp(data, data_zmax, alt, alt_exp, iorder, il1, il2, lev1, &
1863    lev2, nlons, nlevs)
1864  IMPLICIT NONE
1865  ! This routine exponentially damps a longitude by altitude array
1866  ! of data above a specified altitude.
1867
1868  ! Aug. 13/95 - C. McLandress
1869
1870  ! Output arguements:
1871
1872  ! * DATA = modified data array.
1873
1874  ! Input arguements:
1875
1876  ! * DATA    = original data array.
1877  ! * ALT     = altitudes.
1878  ! * ALT_EXP = altitude above which exponential decay applied.
1879  ! * IORDER    = 1 means vertical levels are indexed from top down
1880  ! *           (i.e., highest level indexed 1 and lowest level NLEVS);
1881  ! *           .NE. 1 highest level is index NLEVS.
1882  ! * IL1     = first longitudinal index to use (IL1 >= 1).
1883  ! * IL2     = last longitudinal index to use (IL1 <= IL2 <= NLONS).
1884  ! * LEV1    = first altitude level to use (LEV1 >=1).
1885  ! * LEV2    = last altitude level to use (LEV1 < LEV2 <= NLEVS).
1886  ! * NLONS   = number of longitudes.
1887  ! * NLEVS   = number of vertical
1888
1889  ! Input work arrays:
1890
1891  ! * DATA_ZMAX = data values just above altitude ALT_EXP.
1892
1893  INTEGER iorder, il1, il2, lev1, lev2, nlons, nlevs
1894  REAL alt_exp
1895  REAL data(nlons, nlevs), data_zmax(nlons), alt(nlons, nlevs)
1896
1897  ! Internal variables.
1898
1899  INTEGER levbot, levtop, lincr, i, l
1900  REAL hscale
1901  DATA hscale/5.E3/
1902  ! -----------------------------------------------------------------------
1903
1904  ! Index of lowest altitude level (bottom of drag calculation).
1905
1906  levbot = lev2
1907  levtop = lev1
1908  lincr = 1
1909  IF (iorder/=1) THEN
1910    levbot = lev1
1911    levtop = lev2
1912    lincr = -1
1913  END IF
1914
1915  ! Data values at first level above ALT_EXP.
1916
1917  DO i = il1, il2
1918    DO l = levtop, levbot, lincr
1919      IF (alt(i,l)>=alt_exp) THEN
1920        data_zmax(i) = data(i, l)
1921      END IF
1922    END DO
1923  END DO
1924
1925  ! Exponentially damp field above ALT_EXP to model top at L=1.
1926
1927  DO l = 1, lev2
1928    DO i = il1, il2
1929      IF (alt(i,l)>=alt_exp) THEN
1930        data(i, l) = data_zmax(i)*exp((alt_exp-alt(i,l))/hscale)
1931      END IF
1932    END DO
1933  END DO
1934
1935  RETURN
1936  ! -----------------------------------------------------------------------
1937END SUBROUTINE hines_exp
1938
1939SUBROUTINE vert_smooth(data, work, coeff, nsmooth, il1, il2, lev1, lev2, &
1940    nlons, nlevs)
1941  IMPLICIT NONE
1942  ! Smooth a longitude by altitude array in the vertical over a
1943  ! specified number of levels using a three point smoother.
1944
1945  ! NOTE: input array DATA is modified on output!
1946
1947  ! Aug. 3/95 - C. McLandress
1948
1949  ! Output arguement:
1950
1951  ! * DATA    = smoothed array (on output).
1952
1953  ! Input arguements:
1954
1955  ! * DATA    = unsmoothed array of data (on input).
1956  ! * WORK    = work array of same dimension as DATA.
1957  ! * COEFF   = smoothing coefficient for a 1:COEFF:1 stencil.
1958  ! *           (e.g., COEFF = 2 will result in a smoother which
1959  ! *           weights the level L gridpoint by two and the two
1960  ! *           adjecent levels (L+1 and L-1) by one).
1961  ! * NSMOOTH = number of times to smooth in vertical.
1962  ! *           (e.g., NSMOOTH=1 means smoothed only once,
1963  ! *           NSMOOTH=2 means smoothing repeated twice, etc.)
1964  ! * IL1     = first longitudinal index to use (IL1 >= 1).
1965  ! * IL2     = last longitudinal index to use (IL1 <= IL2 <= NLONS).
1966  ! * LEV1    = first altitude level to use (LEV1 >=1).
1967  ! * LEV2    = last altitude level to use (LEV1 < LEV2 <= NLEVS).
1968  ! * NLONS   = number of longitudes.
1969  ! * NLEVS   = number of vertical levels.
1970
1971  ! Subroutine arguements.
1972
1973  INTEGER nsmooth, il1, il2, lev1, lev2, nlons, nlevs
1974  REAL coeff
1975  REAL data(nlons, nlevs), work(nlons, nlevs)
1976
1977  ! Internal variables.
1978
1979  INTEGER i, l, ns, lev1p, lev2m
1980  REAL sum_wts
1981  ! -----------------------------------------------------------------------
1982
1983  ! Calculate sum of weights.
1984
1985  sum_wts = coeff + 2.
1986
1987  lev1p = lev1 + 1
1988  lev2m = lev2 - 1
1989
1990  ! Smooth NSMOOTH times
1991
1992  DO ns = 1, nsmooth
1993
1994    ! Copy data into work array.
1995
1996    DO l = lev1, lev2
1997      DO i = il1, il2
1998        work(i, l) = data(i, l)
1999      END DO
2000    END DO
2001
2002    ! Smooth array WORK in vertical direction and put into DATA.
2003
2004    DO l = lev1p, lev2m
2005      DO i = il1, il2
2006        data(i, l) = (work(i,l+1)+coeff*work(i,l)+work(i,l-1))/sum_wts
2007      END DO
2008    END DO
2009
2010  END DO
2011
2012  RETURN
2013  ! -----------------------------------------------------------------------
2014END SUBROUTINE vert_smooth
2015
2016
2017
2018
2019
2020
Note: See TracBrowser for help on using the repository browser.