source: LMDZ6/branches/Amaury_dev/libf/phylmd/hines_gwd.F90 @ 5201

Last change on this file since 5201 was 5158, checked in by abarral, 7 weeks ago

Add missing klon on strataer_emiss_mod.F90
Correct various missing explicit declarations
Replace tabs by spaces (tabs are not part of the fortran charset)
Continue cleaning modules
Removed unused arguments and variables

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