GFS Physics Documentation
gwdps.f
Go to the documentation of this file.
1 
4 
115 
337  SUBROUTINE gwdps(IM,IX,IY,KM,A,B,C,U1,V1,T1,Q1,KPBL, &
338  & PRSI,DEL,PRSL,PRSLK,PHII, PHIL,DELTIM,KDT, &
339  & HPRIME,OC,OA4,CLX4,THETA,SIGMA,GAMMA,ELVMAX, &
340  & DUSFC,DVSFC,G, CP, RD, RV, IMX, &
341  & nmtvr, cdmbgwd, me, lprnt, ipr, rdxzb)
342 !
343 ! ********************************************************************
344 ! -----> I M P L E M E N T A T I O N V E R S I O N <----------
345 !
346 ! --- Not in this code -- History of GWDP at NCEP----
347 ! ---------------- -----------------------
348 ! VERSION 3 MODIFIED FOR GRAVITY WAVES, LOCATION: .FR30(V3GWD) *J*
349 !--- 3.1 INCLUDES VARIABLE SATURATION FLUX PROFILE CF ISIGST
350 !--- 3.G INCLUDES PS COMBINED W/ PH (GLAS AND GFDL)
351 !----- ALSO INCLUDED IS RI SMOOTH OVER A THICK LOWER LAYER
352 !----- ALSO INCLUDED IS DECREASE IN DE-ACC AT TOP BY 1/2
353 !----- THE NMC GWD INCORPORATING BOTH GLAS(P&S) AND GFDL(MIGWD)
354 !----- MOUNTAIN INDUCED GRAVITY WAVE DRAG
355 !----- CODE FROM .FR30(V3MONNX) FOR MONIN3
356 !----- THIS VERSION (06 MAR 1987)
357 !----- THIS VERSION (26 APR 1987) 3.G
358 !----- THIS VERSION (01 MAY 1987) 3.9
359 !----- CHANGE TO FORTRAN 77 (FEB 1989) --- HANN-MING HENRY JUANG
360 !----- 20070601 ELVMAX bug fix (*j*)
361 !
362 ! VERSION 4
363 ! ----- This code -----
364 !
365 !----- MODIFIED TO IMPLEMENT THE ENHANCED LOW TROPOSPHERIC GRAVITY
366 !----- WAVE DRAG DEVELOPED BY KIM AND ARAKAWA(JAS, 1995).
367 ! Orographic Std Dev (hprime), Convexity (OC), Asymmetry (OA4)
368 ! and Lx (CLX4) are input topographic statistics needed.
369 !
370 !----- PROGRAMMED AND DEBUGGED BY HONG, ALPERT AND KIM --- JAN 1996.
371 !----- debugged again - moorthi and iredell --- may 1998.
372 !-----
373 ! Further Cleanup, optimization and modification
374 ! - S. Moorthi May 98, March 99.
375 !----- modified for usgs orography data (ncep office note 424)
376 ! and with several bugs fixed - moorthi and hong --- july 1999.
377 !
378 !----- Modified & implemented into NRL NOGAPS
379 ! - Young-Joon Kim, July 2000
380 !-----
381 ! VERSION lm MB (6): oz fix 8/2003
382 ! ----- This code -----
383 !
384 !------ Changed to include the Lott and Miller Mtn Blocking
385 ! with some modifications by (*j*) 4/02
386 ! From a Principal Coordinate calculation using the
387 ! Hi Res 8 minute orography, the Angle of the
388 ! mtn with that to the East (x) axis is THETA, the slope
389 ! parameter SIGMA. The anisotropy is in GAMMA - all are input
390 ! topographic statistics needed. These are calculated off-line
391 ! as a function of model resolution in the fortran code ml01rg2.f,
392 ! with script mlb2.sh. (*j*)
393 !----- gwdps_mb.f version (following lmi) elvmax < hncrit (*j*)
394 ! MB3a expt to enhance elvmax mtn hgt see sigfac & hncrit
395 ! gwdps_GWDFIX_v6.f FIXGWD GF6.0 20070608 sigfac=4.
396 !-----
397 !----------------------------------------------------------------------C
398 ! USE
399 ! ROUTINE IS CALLED FROM GBPHYS (AFTER CALL TO MONNIN)
400 !
401 ! PURPOSE
402 ! USING THE GWD PARAMETERIZATIONS OF PS-GLAS AND PH-
403 ! GFDL TECHNIQUE. THE TIME TENDENCIES OF U V
404 ! ARE ALTERED TO INCLUDE THE EFFECT OF MOUNTAIN INDUCED
405 ! GRAVITY WAVE DRAG FROM SUB-GRID SCALE OROGRAPHY INCLUDING
406 ! CONVECTIVE BREAKING, SHEAR BREAKING AND THE PRESENCE OF
407 ! CRITICAL LEVELS
408 !
409 ! INPUT
410 ! A(IY,KM) NON-LIN TENDENCY FOR V WIND COMPONENT
411 ! B(IY,KM) NON-LIN TENDENCY FOR U WIND COMPONENT
412 ! C(IY,KM) NON-LIN TENDENCY FOR TEMPERATURE
413 ! U1(IX,KM) ZONAL WIND M/SEC AT T0-DT
414 ! V1(IX,KM) MERIDIONAL WIND M/SEC AT T0-DT
415 ! T1(IX,KM) TEMPERATURE DEG K AT T0-DT
416 ! Q1(IX,KM) SPECIFIC HUMIDITY AT T0-DT
417 !
418 ! DELTIM TIME STEP SECS
419 ! SI(N) P/PSFC AT BASE OF LAYER N
420 ! SL(N) P/PSFC AT MIDDLE OF LAYER N
421 ! DEL(N) POSITIVE INCREMENT OF P/PSFC ACROSS LAYER N
422 ! KPBL(IM) is the index of the top layer of the PBL
423 ! ipr & lprnt for diagnostics
424 !
425 ! OUTPUT
426 ! A, B AS AUGMENTED BY TENDENCY DUE TO GWDPS
427 ! OTHER INPUT VARIABLES UNMODIFIED.
428 ! revision log:
429 ! May 2013 J. Wang change cleff back to opn setting
430 ! Jan 2014 J. Wang merge Henry and Fangin's dissipation heat in gfs to nems
431 !
432 !
433 ! ********************************************************************
434  USE machine , ONLY : kind_phys
435  implicit none
436  integer im, iy, ix, km, imx, kdt, ipr, me
437  integer KPBL(im) ! Index for the PBL top layer!
438  real(kind=kind_phys) deltim, G, CP, RD, RV, cdmbgwd(2)
439  real(kind=kind_phys) A(iy,km), B(iy,km), C(iy,km), &
440  & U1(IX,KM), V1(IX,KM), T1(IX,KM), &
441  & Q1(IX,KM), PRSI(IX,KM+1), DEL(IX,KM), &
442  & PRSL(IX,KM), PRSLK(IX,KM), PHIL(IX,KM), &
443  & PHII(IX,KM+1),RDXZB(IY)
444  real(kind=kind_phys) OC(im), OA4(iy,4), CLX4(iy,4) &
445  &, HPRIME(IM)
446 ! for lm mtn blocking
447  real(kind=kind_phys) ELVMAX(im),THETA(im),SIGMA(im),GAMMA(im)
448  real(kind=kind_phys) wk(im)
449  real(kind=kind_phys) bnv2lm(im,km),PE(im),EK(im),ZBK(im),UP(im)
450  real(kind=kind_phys) DB(im,km),ANG(im,km),UDS(im,km)
451  real(kind=kind_phys) ZLEN, DBTMP, R, PHIANG, CDmb, DBIM
452  real(kind=kind_phys) ENG0, ENG1
453 !
454 ! Some constants
455 !
456  real(kind=kind_phys) pi, dw2min, rimin, ric, bnv2min, efmin
457  &, efmax,hpmax,hpmin, rad_to_deg, deg_to_rad
458  parameter(pi=3.1415926535897931)
459  parameter(rad_to_deg=180.0/pi, deg_to_rad=pi/180.0)
460  parameter(dw2min=1., rimin=-100., ric=0.25, bnv2min=1.0e-5)
461 ! PARAMETER (EFMIN=0.0, EFMAX=10.0, hpmax=200.0)
462  parameter(efmin=0.0, efmax=10.0, hpmax=2400.0, hpmin=1.0)
463 !
464  real(kind=kind_phys) FRC, CE, CEOFRC, frmax, CG, GMAX
465  &, veleps, factop, rlolev, rdi
466 ! &, CRITAC, VELEPS, FACTOP, RLOLEV, RDI
467  parameter(frc=1.0, ce=0.8, ceofrc=ce/frc, frmax=100., cg=0.5)
468  parameter(gmax=1.0, veleps=1.0, factop=0.5)
469 ! parameter (GMAX=1.0, CRITAC=5.0E-4, VELEPS=1.0, FACTOP=0.5)
470  parameter(rlolev=50000.0)
471 ! parameter (RLOLEV=500.0)
472 ! parameter (RLOLEV=0.5)
473 !
474  real(kind=kind_phys) dpmin,hminmt,hncrit,minwnd,sigfac
475 ! --- for lm mtn blocking
476 ! parameter (cdmb = 1.0) ! non-dim sub grid mtn drag Amp (*j*)
477  parameter(hncrit=8000.) ! Max value in meters for ELVMAX (*j*)
478 ! hncrit set to 8000m and sigfac added to enhance elvmax mtn hgt
479  parameter(sigfac=4.0) ! MB3a expt test for ELVMAX factor (*j*)
480  parameter(hminmt=50.) ! min mtn height (*j*)
481  parameter(minwnd=0.1) ! min wind component (*j*)
482 
483 ! parameter (dpmin=00.0) ! Minimum thickness of the reference layer
484 !! parameter (dpmin=05.0) ! Minimum thickness of the reference layer
485 ! parameter (dpmin=20.0) ! Minimum thickness of the reference layer
486  ! in centibars
487  parameter(dpmin=5000.0) ! Minimum thickness of the reference layer
488  ! in Pa
489 !
490  real(kind=kind_phys) FDIR
491  integer mdir
492  parameter(mdir=8, fdir=mdir/(pi+pi))
493  integer nwdir(mdir)
494  data nwdir/6,7,5,8,2,3,1,4/
495  save nwdir
496 !
497  LOGICAL ICRILV(im)
498 !
499 !---- MOUNTAIN INDUCED GRAVITY WAVE DRAG
500 !
501  real(kind=kind_phys) TAUB(im), XN(im), YN(im), UBAR(im) &
502  &, VBAR(IM), ULOW(IM), OA(IM), CLX(IM) &
503  &, ROLL(IM), ULOI(IM), DUSFC(IM), DVSFC(IM) &
504  &, DTFAC(IM), XLINV(IM), DELKS(IM), DELKS1(IM)
505 !
506  real(kind=kind_phys) BNV2(im,km), TAUP(im,km+1), ri_n(im,km) &
507  &, TAUD(IM,KM), RO(IM,KM), VTK(IM,KM) &
508  &, VTJ(IM,KM), SCOR(IM), VELCO(IM,KM-1) &
509  &, bnv2bar(im)
510 !
511 ! real(kind=kind_phys) VELKO(KM-1)
512  Integer kref(im), kint(im), iwk(im), ipt(im)
513 ! for lm mtn blocking
514  Integer kreflm(im), iwklm(im)
515  Integer idxzb(im), ktrial, klevm1, nmtvr
516 !
517  real(kind=kind_phys) gor, gocp, fv, gr2, bnv, fr &
518  &, brvf, cleff, tem, tem1, tem2, temc, temv &
519  &, wdir, ti, rdz, dw2, shr2, bvf2 &
520  &, rdelks, efact, coefm, gfobnv &
521  &, scork, rscor, hd, fro, rim, sira &
522  &, dtaux, dtauy, pkp1log, pklog
523  integer kmm1, kmm2, lcap, lcapp1, kbps, kbpsp1,kbpsm1 &
524  &, kmps, idir, nwd, i, j, k, klcap, kp1, kmpbl, npt, npr &
525  &, kmll
526 ! &, kmll,kmds,ihit,jhit
527  logical lprnt
528 !
529 ! parameter (cdmb = 1.0) ! non-dim sub grid mtn drag Amp (*j*)
530 ! non-dim sub grid mtn drag Amp (*j*)
531 ! cdmb = 1.0/float(IMX/192)
532 ! cdmb = 192.0/float(IMX)
533  cdmb = 4.0 * 192.0/float(imx)
534  if (cdmbgwd(1) >= 0.0) cdmb = cdmb * cdmbgwd(1)
535 !
536  npr = 0
537  DO i = 1, im
538  dusfc(i) = 0.
539  dvsfc(i) = 0.
540  ENDDO
541 !
542  DO k = 1, km
543  DO i = 1, im
544  db(i,k) = 0.
545  ang(i,k) = 0.
546  uds(i,k) = 0.
547  ENDDO
548  ENDDO
549 !
550  rdi = 1.0 / rd
551  gor = g/rd
552  gr2 = g*gor
553  gocp = g/cp
554  fv = rv/rd - 1
555 !
556 ! NCNT = 0
557  kmm1 = km - 1
558  kmm2 = km - 2
559  lcap = km
560  lcapp1 = lcap + 1
561 !
562 !
563  IF ( nmtvr .eq. 14) then
564 ! ---- for lm and gwd calculation points
565  rdxzb(:) = 0
566  ipt = 0
567  npt = 0
568  DO i = 1,im
569  IF ( (elvmax(i) .GT. hminmt)
570  & .and. (hprime(i) .GT. hpmin) ) then
571  npt = npt + 1
572  ipt(npt) = i
573  if (ipr .eq. i) npr = npt
574  ENDIF
575  ENDDO
576  IF (npt .eq. 0) RETURN ! No gwd/mb calculation done!
577 !
578 ! if (lprnt) print *,' npt=',npt,' npr=',npr,' ipr=',ipr,' im=',im
579 ! &,' ipt(npt)=',ipt(npt)
580 !
581 ! --- iwklm is the level above the height of the of the mountain.
582 ! --- idxzb is the level of the dividing streamline.
583 ! INITIALIZE DIVIDING STREAMLINE (DS) CONTROL VECTOR
584 !
585  do i=1,npt
586  iwklm(i) = 2
587  idxzb(i) = 0
588  kreflm(i) = 0
589  enddo
590 ! if (lprnt)
591 ! & print *,' in gwdps_lm.f npt,IM,IX,IY,km,me=',npt,IM,IX,IY,km,me
592 !
593 !
595 !
596 !..............................
597 !..............................
598 !
599 ! (*j*) 11/03: test upper limit on KMLL=km - 1
600 ! then do not need hncrit -- test with large hncrit first.
601 ! KMLL = km / 2 ! maximum mtnlm height : # of vertical levels / 2
602  kmll = kmm1
603 ! --- No mtn should be as high as KMLL (so we do not have to start at
604 ! --- the top of the model but could do calc for all levels).
605 !
606  DO i = 1, npt
607  j = ipt(i)
608  elvmax(j) = min(elvmax(j) + sigfac * hprime(j), hncrit)
609  ENDDO
610 !
611  DO k = 1,kmll
612  DO i = 1, npt
613  j = ipt(i)
614 ! --- interpolate to max mtn height for index, iwklm(I) wk[gz]
615 ! --- ELVMAX is limited to hncrit because to hi res topo30 orog.
616  pkp1log = phil(j,k+1) / g
617  pklog = phil(j,k) / g
618 !!!------- ELVMAX(J) = min (ELVMAX(J) + sigfac * hprime(j), hncrit)
619  if ( ( elvmax(j) .le. pkp1log ) .and.
620  & ( elvmax(j) .ge. pklog ) ) THEN
621 ! print *,' in gwdps_lm.f 1 =',k,ELVMAX(j),pklog,pkp1log,me
622 ! --- wk for diags but can be saved and reused.
623  wk(i) = g * elvmax(j) / ( phil(j,k+1) - phil(j,k) )
624  iwklm(i) = max(iwklm(i), k+1 )
625 ! print *,' in gwdps_lm.f 2 npt=',npt,i,j,wk(i),iwklm(i),me
626  endif
627 !
628 ! --- find at prsl levels large scale environment variables
629 ! --- these cover all possible mtn max heights
630  vtj(i,k) = t1(j,k) * (1.+fv*q1(j,k))
631  vtk(i,k) = vtj(i,k) / prslk(j,k)
632  ro(i,k) = rdi * prsl(j,k) / vtj(i,k) ! DENSITY Kg/M**3
633  ENDDO
634  ENDDO
635 !
636 ! testing for highest model level of mountain top
637 !
638 ! ihit = 2
639 ! jhit = 0
640 ! do i = 1, npt
641 ! j=ipt(i)
642 ! if ( iwklm(i) .gt. ihit ) then
643 ! ihit = iwklm(i)
644 ! jhit = j
645 ! endif
646 ! enddo
647 ! print *, ' mb: kdt,max(iwklm),jhit,phil,me=',
648 ! & kdt,ihit,jhit,phil(jhit,ihit),me
649 
650  klevm1 = kmll - 1
651  DO k = 1, klevm1
652  DO i = 1, npt
653  j = ipt(i)
654  rdz = g / ( phil(j,k+1) - phil(j,k) )
655 ! --- Brunt-Vaisala Frequency
657  bnv2lm(i,k) = (g+g) * rdz * ( vtk(i,k+1)-vtk(i,k) )
658  & / ( vtk(i,k+1)+vtk(i,k) )
659  bnv2lm(i,k) = max( bnv2lm(i,k), bnv2min )
660  ENDDO
661  ENDDO
662 ! print *,' in gwdps_lm.f 3 npt=',npt,j,RDZ,me
663 !
664  DO i = 1, npt
665  j = ipt(i)
666  delks(i) = 1.0 / (prsi(j,1) - prsi(j,iwklm(i)))
667  delks1(i) = 1.0 / (prsl(j,1) - prsl(j,iwklm(i)))
668  ubar(i) = 0.0
669  vbar(i) = 0.0
670  roll(i) = 0.0
671  pe(i) = 0.0
672  ek(i) = 0.0
673  bnv2bar(i) = (prsl(j,1)-prsl(j,2)) * delks1(i) * bnv2lm(i,1)
674  ENDDO
675 
676 ! --- find the dividing stream line height
677 ! --- starting from the level above the max mtn downward
678 ! --- iwklm(i) is the k-index of mtn elvmax elevation
681  DO ktrial = kmll, 1, -1
682  DO i = 1, npt
683  IF ( ktrial .LT. iwklm(i) .and. kreflm(i) .eq. 0 ) then
684  kreflm(i) = ktrial
685  ENDIF
686  ENDDO
687  ENDDO
688 ! print *,' in gwdps_lm.f 4 npt=',npt,kreflm(npt),me
689 !
690 ! --- in the layer kreflm(I) to 1 find PE (which needs N, ELVMAX)
691 ! --- make averages, guess dividing stream (DS) line layer.
692 ! --- This is not used in the first cut except for testing and
693 ! --- is the vert ave of quantities from the surface to mtn top.
694 !
695  DO i = 1, npt
696  DO k = 1, kreflm(i)
697  j = ipt(i)
698  rdelks = del(j,k) * delks(i)
699  ubar(i) = ubar(i) + rdelks * u1(j,k) ! trial Mean U below
700  vbar(i) = vbar(i) + rdelks * v1(j,k) ! trial Mean V below
701  roll(i) = roll(i) + rdelks * ro(i,k) ! trial Mean RO below
702  rdelks = (prsl(j,k)-prsl(j,k+1)) * delks1(i)
703  bnv2bar(i) = bnv2bar(i) + bnv2lm(i,k) * rdelks
704 ! --- these vert ave are for diags, testing and GWD to follow (*j*).
705  ENDDO
706  ENDDO
707 ! print *,' in gwdps_lm.f 5 =',i,kreflm(npt),BNV2bar(npt),me
708 !
709 ! --- integrate to get PE in the trial layer.
710 ! --- Need the first layer where PE>EK - as soon as
711 ! --- IDXZB is not 0 we have a hit and Zb is found.
712 !
713  DO i = 1, npt
714  j = ipt(i)
715  DO k = iwklm(i), 1, -1
716  phiang = atan2(v1(j,k),u1(j,k))*rad_to_deg
717  ang(i,k) = ( theta(j) - phiang )
718  if ( ang(i,k) .gt. 90. ) ang(i,k) = ang(i,k) - 180.
719  if ( ang(i,k) .lt. -90. ) ang(i,k) = ang(i,k) + 180.
720  ang(i,k) = ang(i,k) * deg_to_rad
721 !
728  uds(i,k) =
729  & max(sqrt(u1(j,k)*u1(j,k) + v1(j,k)*v1(j,k)), minwnd)
730 ! --- Test to see if we found Zb previously
731  IF (idxzb(i) .eq. 0 ) then
732  pe(i) = pe(i) + bnv2lm(i,k) *
733  & ( g * elvmax(j) - phil(j,k) ) *
734  & ( phii(j,k+1) - phii(j,k) ) / (g*g)
735 ! --- KE
736 ! --- Wind projected on the line perpendicular to mtn range, U(Zb(K)).
737 ! --- kenetic energy is at the layer Zb
738 ! --- THETA ranges from -+90deg |_ to the mtn "largest topo variations"
739  up(i) = uds(i,k) * cos(ang(i,k))
740  ek(i) = 0.5 * up(i) * up(i)
741 
742 ! --- Dividing Stream lime is found when PE =exceeds EK.
743  IF ( pe(i) .ge. ek(i) ) THEN
744  idxzb(i) = k
745  rdxzb(j) = real(k,kind=kind_phys)
746  ENDIF
747 ! --- Then mtn blocked flow is between Zb=k(IDXZB(I)) and surface
748 !
762  ENDIF
763  ENDDO
764  ENDDO
765 !
766 ! print *,' in gwdps_lm.f 6 =',phiang,THETA(ipt(npt)),me
767 ! print *,' in gwdps_lm.f 7 =',IDXZB(npt),PE(npt)
768 !
769 ! if (lprnt .and. npr .gt. 0) then
770 ! print *,' BNV2bar,BNV2lm=',bnv2bar(npr),BNV2lm(npr,1:klevm1)
771 ! print *,' npr,IDXZB,UDS=',npr,IDXZB(npr),UDS(npr,:)
772 ! print *,' PE,UP,EK=',PE(npr),UP(npr),EK(npr)
773 ! endif
774 !
775  DO i = 1, npt
776  j = ipt(i)
777 ! --- Calc if N constant in layers (Zb guess) - a diagnostic only.
778  zbk(i) = elvmax(j)
779  & - sqrt(ubar(i)*ubar(i) + vbar(i)*vbar(i))/bnv2bar(i)
780  ENDDO
781 !
782 ! if (lprnt .and. npr .gt. 0) then
783 ! print *,' iwklm,ZBK=',iwklm(npr),ZBK(npr),IDXZB(npr)
784 ! print *,' Zb=',PHIL(ipr),IDXZB(npr))/G
785 ! print *,' in gwdps_lm.f 8 npt =',npt,ZBK(npt),UP(npt),me
786 ! endif
787 !
788 ! --- The drag for mtn blocked flow
789 !
790  DO i = 1, npt
791  j = ipt(i)
792  zlen = 0.
793 ! print *,' in gwdps_lm.f 9 =',i,j,IDXZB(i),me
794  IF ( idxzb(i) .gt. 0 ) then
795  DO k = idxzb(i), 1, -1
796  IF ( phil(j,idxzb(i)) .gt. phil(j,k) ) then
797 
805  zlen = sqrt( ( phil(j,idxzb(i)) - phil(j,k) ) /
806  & ( phil(j,k ) + g * hprime(j) ) )
807 ! --- lm eq 14:
817  r = (cos(ang(i,k))**2 + gamma(j) * sin(ang(i,k))**2) /
818  & (gamma(j) * cos(ang(i,k))**2 + sin(ang(i,k))**2)
819 ! --- (negitive of DB -- see sign at tendency)
829 
830  dbtmp = 0.25 * cdmb *
831  & max( 2. - 1. / r, 0. ) * sigma(j) *
832  & max(cos(ang(i,k)), gamma(j)*sin(ang(i,k))) *
833  & zlen / hprime(j)
834  db(i,k) = dbtmp * uds(i,k)
835 !
836 ! if(lprnt .and. i .eq. npr) then
837 ! print *,' in gwdps_lmi.f 10 npt=',npt,i,j,idxzb(i)
838 ! &, DBTMP,R' ang=',ang(i,k),' gamma=',gamma(j),' K=',K
839 ! print *,' in gwdps_lmi.f 11 K=',k,ZLEN,cos(ANG(I,K))
840 ! print *,' in gwdps_lmi.f 12 DB=',DB(i,k),sin(ANG(I,K))
841 ! endif
842  endif
843  ENDDO
844 ! if(lprnt) print *,' @K=1,ZLEN,DBTMP=',K,ZLEN,DBTMP
845  endif
846  ENDDO
847 !
848 !.............................
849 !.............................
850 ! end mtn blocking section
851 !
852  ELSEIF ( nmtvr .ne. 14) then
853 ! ---- for mb not present and gwd (nmtvr .ne .14)
854  ipt = 0
855  npt = 0
856  DO i = 1,im
857  IF ( hprime(i) .GT. hpmin ) then
858  npt = npt + 1
859  ipt(npt) = i
860  if (ipr .eq. i) npr = npt
861  ENDIF
862  ENDDO
863  IF (npt .eq. 0) RETURN ! No gwd/mb calculation done!
864 !
865 ! if (lprnt) print *,' NPR=',npr,' npt=',npt,' IPR=',IPR
866 ! &,' ipt(npt)=',ipt(npt)
867 !
868  do i=1,npt
869  idxzb(i) = 0
870  rdxzb(i) = 0.
871  enddo
872  ENDIF
873 !
874 !.............................
875 !.............................
876 !
878  kmpbl = km / 2 ! maximum pbl height : # of vertical levels / 2
879 !
880 ! Scale cleff between IM=384*2 and 192*2 for T126/T170 and T62
881 !
882  if (imx .gt. 0) then
883 ! cleff = 1.0E-5 * SQRT(FLOAT(IMX)/384.0) ! this is inverse of CLEFF!
884 ! cleff = 1.0E-5 * SQRT(FLOAT(IMX)/192.0) ! this is inverse of CLEFF!
885 ! cleff = 0.5E-5 * SQRT(FLOAT(IMX)/192.0) ! this is inverse of CLEFF!
886 ! cleff = 1.0E-5 * SQRT(FLOAT(IMX)/192)/float(IMX/192)
887 ! cleff = 1.0E-5 / SQRT(FLOAT(IMX)/192.0) ! this is inverse of CLEFF!
888  cleff = 0.5e-5 / sqrt(float(imx)/192.0) ! this is inverse of CLEFF!
889 ! hmhj for ndsl
890 ! jw cleff = 0.1E-5 / SQRT(FLOAT(IMX)/192.0) ! this is inverse of CLEFF!
891 ! cleff = 2.0E-5 * SQRT(FLOAT(IMX)/192.0) ! this is inverse of CLEFF!
892 ! cleff = 2.5E-5 * SQRT(FLOAT(IMX)/192.0) ! this is inverse of CLEFF!
893  endif
894  if (cdmbgwd(2) >= 0.0) cleff = cleff * cdmbgwd(2)
895 !
896  DO k = 1,km
897  DO i =1,npt
898  j = ipt(i)
899  vtj(i,k) = t1(j,k) * (1.+fv*q1(j,k))
900  vtk(i,k) = vtj(i,k) / prslk(j,k)
901  ro(i,k) = rdi * prsl(j,k) / vtj(i,k) ! DENSITY TONS/M**3
902  taup(i,k) = 0.0
903  ENDDO
904  ENDDO
905  DO k = 1,kmm1
906  DO i =1,npt
907  j = ipt(i)
908  ti = 2.0 / (t1(j,k)+t1(j,k+1))
909  tem = ti / (prsl(j,k)-prsl(j,k+1))
910  rdz = g / (phil(j,k+1) - phil(j,k))
911  tem1 = u1(j,k) - u1(j,k+1)
912  tem2 = v1(j,k) - v1(j,k+1)
913  dw2 = tem1*tem1 + tem2*tem2
914  shr2 = max(dw2,dw2min) * rdz * rdz
915  bvf2 = g*(gocp+rdz*(vtj(i,k+1)-vtj(i,k))) * ti
916  ri_n(i,k) = max(bvf2/shr2,rimin) ! Richardson number
917 ! Brunt-Vaisala Frequency
918 ! TEM = GR2 * (PRSL(J,K)+PRSL(J,K+1)) * TEM
919 ! BNV2(I,K) = TEM * (VTK(I,K+1)-VTK(I,K))/(VTK(I,K+1)+VTK(I,K))
920  bnv2(i,k) = (g+g) * rdz * (vtk(i,k+1)-vtk(i,k))
921  & / (vtk(i,k+1)+vtk(i,k))
922  bnv2(i,k) = max( bnv2(i,k), bnv2min )
923  ENDDO
924  ENDDO
925 ! print *,' in gwdps_lm.f GWD:14 =',npt,kmm1,bnv2(npt,kmm1)
926 !
927 ! Apply 3 point smoothing on BNV2
928 !
929 ! do k=1,km
930 ! do i=1,im
931 ! vtk(i,k) = bnv2(i,k)
932 ! enddo
933 ! enddo
934 ! do k=2,kmm1
935 ! do i=1,im
936 ! bnv2(i,k) = 0.25*(vtk(i,k-1)+vtk(i,k+1)) + 0.5*vtk(i,k)
937 ! enddo
938 ! enddo
939 !
940 ! Finding the first interface index above 50 hPa level
941 !
942  do i=1,npt
943  iwk(i) = 2
944  enddo
945  DO k=3,kmpbl
946  DO i=1,npt
947  j = ipt(i)
948  tem = (prsi(j,1) - prsi(j,k))
949  if (tem .lt. dpmin) iwk(i) = k
950  enddo
951  enddo
952 !
955  kbps = 1
956  kmps = km
957  DO i=1,npt
958  j = ipt(i)
959  kref(i) = max(iwk(i), kpbl(j)+1 ) ! reference level
960  delks(i) = 1.0 / (prsi(j,1) - prsi(j,kref(i)))
961  delks1(i) = 1.0 / (prsl(j,1) - prsl(j,kref(i)))
962  ubar(i) = 0.0
963  vbar(i) = 0.0
964  roll(i) = 0.0
965  kbps = max(kbps, kref(i))
966  kmps = min(kmps, kref(i))
967 !
968  bnv2bar(i) = (prsl(j,1)-prsl(j,2)) * delks1(i) * bnv2(i,1)
969  ENDDO
970 ! print *,' in gwdps_lm.f GWD:15 =',KBPS,KMPS
971  kbpsp1 = kbps + 1
972  kbpsm1 = kbps - 1
973  DO k = 1,kbps
974  DO i = 1,npt
975  IF (k .LT. kref(i)) THEN
976  j = ipt(i)
977  rdelks = del(j,k) * delks(i)
978  ubar(i) = ubar(i) + rdelks * u1(j,k) ! Mean U below kref
979  vbar(i) = vbar(i) + rdelks * v1(j,k) ! Mean V below kref
980 !
981  roll(i) = roll(i) + rdelks * ro(i,k) ! Mean RO below kref
982  rdelks = (prsl(j,k)-prsl(j,k+1)) * delks1(i)
983  bnv2bar(i) = bnv2bar(i) + bnv2(i,k) * rdelks
984  ENDIF
985  ENDDO
986  ENDDO
987 ! print *,' in gwdps_lm.f GWD:15B =',bnv2bar(npt)
988 !
989 ! FIGURE OUT LOW-LEVEL HORIZONTAL WIND DIRECTION AND FIND 'OA'
990 !
991 ! NWD 1 2 3 4 5 6 7 8
992 ! WD W S SW NW E N NE SE
993 !
996  DO i = 1,npt
997  j = ipt(i)
998  wdir = atan2(ubar(i),vbar(i)) + pi
999  idir = mod(nint(fdir*wdir),mdir) + 1
1000  nwd = nwdir(idir)
1001  oa(i) = (1-2*int( (nwd-1)/4 )) * oa4(j,mod(nwd-1,4)+1)
1002  clx(i) = clx4(j,mod(nwd-1,4)+1)
1003  ENDDO
1004 !
1005 !-----XN,YN "LOW-LEVEL" WIND PROJECTIONS IN ZONAL
1006 ! & MERIDIONAL DIRECTIONS
1007 !-----ULOW "LOW-LEVEL" WIND MAGNITUDE - (= U)
1008 !-----BNV2 BNV2 = N**2
1009 !-----TAUB BASE MOMENTUM FLUX
1010 !-----= -(RO * U**3/(N*XL)*GF(FR) FOR N**2 > 0
1011 !-----= 0. FOR N**2 < 0
1012 !-----FR FROUDE = N*HPRIME / U
1013 !-----G GMAX*FR**2/(FR**2+CG/OC)
1014 !
1015 !-----INITIALIZE SOME ARRAYS
1016 !
1017  DO i = 1,npt
1018  xn(i) = 0.0
1019  yn(i) = 0.0
1020  taub(i) = 0.0
1021  ulow(i) = 0.0
1022  dtfac(i) = 1.0
1023  icrilv(i) = .false. ! INITIALIZE CRITICAL LEVEL CONTROL VECTOR
1024 
1025 !
1026 !----COMPUTE THE "LOW LEVEL" WIND MAGNITUDE (M/S)
1027 !
1028  ulow(i) = max(sqrt(ubar(i)*ubar(i) + vbar(i)*vbar(i)), 1.0)
1029  uloi(i) = 1.0 / ulow(i)
1030  ENDDO
1031 !
1032  DO k = 1,kmm1
1033  DO i = 1,npt
1034  j = ipt(i)
1035  velco(i,k) = 0.5 * ((u1(j,k)+u1(j,k+1))*ubar(i)
1036  & + (v1(j,k)+v1(j,k+1))*vbar(i))
1037  velco(i,k) = velco(i,k) * uloi(i)
1038 ! IF ((VELCO(I,K).LT.VELEPS) .AND. (VELCO(I,K).GT.0.)) THEN
1039 ! VELCO(I,K) = VELEPS
1040 ! ENDIF
1041  ENDDO
1042  ENDDO
1043 !
1044 !
1045 ! find the interface level of the projected wind where
1046 ! low levels & upper levels meet above pbl
1047 !
1048 ! do i=1,npt
1049 ! kint(i) = km
1050 ! enddo
1051 ! do k = 1,kmm1
1052 ! do i = 1,npt
1053 ! IF (K .GT. kref(I)) THEN
1054 ! if(velco(i,k) .lt. veleps .and. kint(i) .eq. km) then
1055 ! kint(i) = k+1
1056 ! endif
1057 ! endif
1058 ! enddo
1059 ! enddo
1060 ! WARNING KINT = KREF !!!!!!!!!
1061  do i=1,npt
1062  kint(i) = kref(i)
1063  enddo
1064 !
1065 ! if(lprnt) print *,' ubar=',ubar
1066 ! &,' vbar=',vbar,' ulow=',ulow,' veleps=',veleps
1067 !
1068  DO i = 1,npt
1069  j = ipt(i)
1070  bnv = sqrt( bnv2bar(i) )
1071  fr = bnv * uloi(i) * min(hprime(j),hpmax)
1072  fr = min(fr, frmax)
1073  xn(i) = ubar(i) * uloi(i)
1074  yn(i) = vbar(i) * uloi(i)
1075 !
1076 ! Compute the base level stress and store it in TAUB
1077 ! CALCULATE ENHANCEMENT FACTOR, NUMBER OF MOUNTAINS & ASPECT
1078 ! RATIO CONST. USE SIMPLIFIED RELATIONSHIP BETWEEN STANDARD
1079 ! DEVIATION & CRITICAL HGT
1080 !
1109 
1118 
1119  efact = (oa(i) + 2.) ** (ceofrc*fr)
1120  efact = min( max(efact,efmin), efmax )
1121 !
1122  coefm = (1. + clx(i)) ** (oa(i)+1.)
1123 !
1124  xlinv(i) = coefm * cleff
1125 !
1126  tem = fr * fr * oc(j)
1127  gfobnv = gmax * tem / ((tem + cg)*bnv) ! G/N0
1128 !
1129  taub(i) = xlinv(i) * roll(i) * ulow(i) * ulow(i)
1130  & * ulow(i) * gfobnv * efact ! BASE FLUX Tau0
1131 !
1132 ! tem = min(HPRIME(I),hpmax)
1133 ! TAUB(I) = XLINV(I) * ROLL(I) * ULOW(I) * BNV * tem * tem
1134 !
1135  k = max(1, kref(i)-1)
1136  tem = max(velco(i,k)*velco(i,k), 0.1)
1137  scor(i) = bnv2(i,k) / tem ! Scorer parameter below ref level
1138  ENDDO
1139 ! if(lprnt) print *,' taub=',taub
1140 !
1141 !----SET UP BOTTOM VALUES OF STRESS
1142 !
1143  DO k = 1, kbps
1144  DO i = 1,npt
1145  IF (k .LE. kref(i)) taup(i,k) = taub(i)
1146  ENDDO
1147  ENDDO
1148 !
1149 ! Now compute vertical structure of the stress.
1150 !
1151  DO k = kmps, kmm1 ! Vertical Level K Loop!
1152  kp1 = k + 1
1153  DO i = 1, npt
1154 !
1155 !-----UNSTABLE LAYER IF RI < RIC
1156 !-----UNSTABLE LAYER IF UPPER AIR VEL COMP ALONG SURF VEL <=0 (CRIT LAY)
1157 !---- AT (U-C)=0. CRIT LAYER EXISTS AND BIT VECTOR SHOULD BE SET (.LE.)
1158 !
1159  IF (k .GE. kref(i)) THEN
1160  icrilv(i) = icrilv(i) .OR. ( ri_n(i,k) .LT. ric)
1161  & .OR. (velco(i,k) .LE. 0.0)
1162  ENDIF
1163  ENDDO
1164 !
1177  DO i = 1,npt
1178  IF (k .GE. kref(i)) THEN
1179  IF (.NOT.icrilv(i) .AND. taup(i,k) .GT. 0.0 ) THEN
1180  temv = 1.0 / max(velco(i,k), 0.01)
1181 ! IF (OA(I) .GT. 0. .AND. PRSI(ipt(i),KP1).GT.RLOLEV) THEN
1182  IF (oa(i).GT.0. .AND. kp1 .lt. kint(i)) THEN
1183  scork = bnv2(i,k) * temv * temv
1184  rscor = min(1.0, scork / scor(i))
1185  scor(i) = scork
1186  ELSE
1187  rscor = 1.
1188  ENDIF
1189 !
1201 
1202  brvf = sqrt(bnv2(i,k)) ! Brunt-Vaisala Frequency
1203 ! TEM1 = XLINV(I)*(RO(I,KP1)+RO(I,K))*BRVF*VELCO(I,K)*0.5
1204  tem1 = xlinv(i)*(ro(i,kp1)+ro(i,k))*brvf*0.5
1205  & * max(velco(i,k),0.01)
1206  hd = sqrt(taup(i,k) / tem1)
1207  fro = brvf * hd * temv
1208 !
1209 ! RIM is the MINIMUM-RICHARDSON NUMBER BY SHUTTS (1985)
1210 !
1219 
1220  tem2 = sqrt(ri_n(i,k))
1221  tem = 1. + tem2 * fro
1222  rim = ri_n(i,k) * (1.-fro) / (tem * tem)
1223 !
1224 ! CHECK STABILITY TO EMPLOY THE 'SATURATION HYPOTHESIS'
1225 ! OF LINDZEN (1981) EXCEPT AT TROPOSPHERIC DOWNSTREAM REGIONS
1226 !
1242 ! ----------------------
1243  IF (rim .LE. ric .AND.
1244 ! & (OA(I) .LE. 0. .OR. PRSI(ipt(I),KP1).LE.RLOLEV )) THEN
1245  & (oa(i) .LE. 0. .OR. kp1 .ge. kint(i) )) THEN
1246  temc = 2.0 + 1.0 / tem2
1247  hd = velco(i,k) * (2.*sqrt(temc)-temc) / brvf
1248  taup(i,kp1) = tem1 * hd * hd
1249  ELSE
1250  taup(i,kp1) = taup(i,k) * rscor
1251  ENDIF
1252  taup(i,kp1) = min(taup(i,kp1), taup(i,k))
1253  ENDIF
1254  ENDIF
1255  ENDDO
1256  ENDDO
1257 !
1258 ! DO I=1,IM
1259 ! taup(i,km+1) = taup(i,km)
1260 ! ENDDO
1261 !
1262  IF(lcap .LE. km) THEN
1263  DO klcap = lcapp1, km+1
1264  DO i = 1,npt
1265  sira = prsi(ipt(i),klcap) / prsi(ipt(i),lcap)
1266  taup(i,klcap) = sira * taup(i,lcap)
1267  ENDDO
1268  ENDDO
1269  ENDIF
1270 !
1271 ! Calculate - (g/p*)*d(tau)/d(sigma) and Decel terms DTAUX, DTAUY
1272 !
1273  DO k = 1,km
1274  DO i = 1,npt
1275  taud(i,k) = g * (taup(i,k+1) - taup(i,k)) / del(ipt(i),k)
1276  ENDDO
1277  ENDDO
1278 !
1279 !------LIMIT DE-ACCELERATION (MOMENTUM DEPOSITION ) AT TOP TO 1/2 VALUE
1280 !------THE IDEA IS SOME STUFF MUST GO OUT THE 'TOP'
1281 !
1282  DO klcap = lcap, km
1283  DO i = 1,npt
1284  taud(i,klcap) = taud(i,klcap) * factop
1285  ENDDO
1286  ENDDO
1287 !
1288 !------IF THE GRAVITY WAVE DRAG WOULD FORCE A CRITICAL LINE IN THE
1289 !------LAYERS BELOW SIGMA=RLOLEV DURING THE NEXT DELTIM TIMESTEP,
1290 !------THEN ONLY APPLY DRAG UNTIL THAT CRITICAL LINE IS REACHED.
1291 !
1292  DO k = 1,kmm1
1293  DO i = 1,npt
1294  IF (k .GT. kref(i) .and. prsi(ipt(i),k) .GE. rlolev) THEN
1295  IF(taud(i,k).NE.0.) THEN
1296  tem = deltim * taud(i,k)
1297  dtfac(i) = min(dtfac(i),abs(velco(i,k)/tem))
1298  ENDIF
1299  ENDIF
1300  ENDDO
1301  ENDDO
1302 !
1303 ! if(lprnt .and. npr .gt. 0) then
1304 ! print *,' before A=',A(npr,:)
1305 ! print *,' before B=',B(npr,:)
1306 ! endif
1307 
1312  DO k = 1,km
1313  DO i = 1,npt
1314  j = ipt(i)
1315  taud(i,k) = taud(i,k) * dtfac(i)
1316  dtaux = taud(i,k) * xn(i)
1317  dtauy = taud(i,k) * yn(i)
1318  eng0 = 0.5*(u1(j,k)*u1(j,k)+v1(j,k)*v1(j,k))
1319 ! --- lm mb (*j*) changes overwrite GWD
1320  if ( k .lt. idxzb(i) .AND. idxzb(i) .ne. 0 ) then
1321  dbim = db(i,k) / (1.+db(i,k)*deltim)
1322  a(j,k) = - dbim * v1(j,k) + a(j,k)
1323  b(j,k) = - dbim * u1(j,k) + b(j,k)
1324  eng1 = eng0*(1.0-dbim*deltim)*(1.0-dbim*deltim)
1325 ! if ( ABS(DBIM * U1(J,K)) .gt. .01 )
1326 ! & print *,' in gwdps_lmi.f KDT=',KDT,I,K,DB(I,K),
1327 ! & dbim,idxzb(I),U1(J,K),V1(J,K),me
1328  dusfc(j) = dusfc(j) - dbim * u1(j,k) * del(j,k)
1329  dvsfc(j) = dvsfc(j) - dbim * v1(j,k) * del(j,k)
1330  else
1331 !
1332  a(j,k) = dtauy + a(j,k)
1333  b(j,k) = dtaux + b(j,k)
1334  eng1 = 0.5*(
1335  & (u1(j,k)+dtaux*deltim)*(u1(j,k)+dtaux*deltim)
1336  & + (v1(j,k)+dtauy*deltim)*(v1(j,k)+dtauy*deltim))
1337  dusfc(j) = dusfc(j) + dtaux * del(j,k)
1338  dvsfc(j) = dvsfc(j) + dtauy * del(j,k)
1339  endif
1340  c(j,k) = c(j,k) + max(eng0-eng1,0.)/cp/deltim
1341  ENDDO
1342  ENDDO
1343 ! if (lprnt) then
1344 ! print *,' in gwdps_lm.f after A=',A(ipr,:)
1345 ! print *,' in gwdps_lm.f after B=',B(ipr,:)
1346 ! print *,' DB=',DB(ipr,:)
1347 ! endif
1348  tem = -1.0/g
1349  DO i = 1,npt
1350  j = ipt(i)
1351 ! TEM = (-1.E3/G)
1352  dusfc(j) = tem * dusfc(j)
1353  dvsfc(j) = tem * dvsfc(j)
1354  ENDDO
1355 !
1356 ! MONITOR FOR EXCESSIVE GRAVITY WAVE DRAG TENDENCIES IF NCNT>0
1357 !
1358 ! IF(NCNT.GT.0) THEN
1359 ! IF(LAT.GE.38.AND.LAT.LE.42) THEN
1360 !CMIC$ GUARD 37
1361 ! DO 92 I = 1,IM
1362 ! IF(IKOUNT.GT.NCNT) GO TO 92
1363 ! IF(I.LT.319.OR.I.GT.320) GO TO 92
1364 ! DO 91 K = 1,KM
1365 ! IF(ABS(TAUD(I,K)) .GT. CRITAC) THEN
1366 ! IF(I.LE.IM) THEN
1367 ! IKOUNT = IKOUNT+1
1368 ! PRINT 123,I,LAT,KDT
1369 ! PRINT 124,TAUB(I),BNV(I),ULOW(I),
1370 ! 1 GF(I),FR(I),ROLL(I),HPRIME(I),XN(I),YN(I)
1371 ! PRINT 124,(TAUD(I,KK),KK = 1,KM)
1372 ! PRINT 124,(TAUP(I,KK),KK = 1,KM+1)
1373 ! PRINT 124,(ri_n(I,KK),KK = 1,KM)
1374 ! DO 93 KK = 1,KMM1
1375 ! VELKO(KK) =
1376 ! 1 0.5*((U1(I,KK)+U1(I,KK+1))*UBAR(I)+
1377 ! 2 (V1(I,KK)+V1(I,KK+1))*VBAR(I))*ULOI(I)
1378 !93 CONTINUE
1379 ! PRINT 124,(VELKO(KK),KK = 1,KMM1)
1380 ! PRINT 124,(A (I,KK),KK = 1,KM)
1381 ! PRINT 124,(DTAUY(I,KK),KK = 1,KM)
1382 ! PRINT 124,(B (I,KK),KK = 1,KM)
1383 ! PRINT 124,(DTAUX(I,KK),KK = 1,KM)
1384 ! GO TO 92
1385 ! ENDIF
1386 ! ENDIF
1387 !91 CONTINUE
1388 !92 CONTINUE
1389 !CMIC$ END GUARD 37
1390 !123 FORMAT(' *** MIGWD PRINT *** I=',I3,' LAT=',I3,' KDT=',I3)
1391 !124 FORMAT(2X, 10E13.6)
1392 ! ENDIF
1393 ! ENDIF
1394 !
1395 ! print *,' in gwdps_lm.f 18 =',A(ipt(1),idxzb(1))
1396 ! &, B(ipt(1),idxzb(1)),me
1397  RETURN
1398  END
1399 
subroutine gwdps(IM, IX, IY, KM, A, B, C, U1, V1, T1, Q1, KPBL, PRSI, DEL, PRSL, PRSLK, PHII, PHIL, DELTIM, KDT, HPRIME, OC, OA4, CLX4, THETA, SIGMA, GAMMA, ELVMAX, DUSFC, DVSFC, G, CP, RD, RV, IMX, nmtvr, cdmbgwd, me, lprnt, ipr, rdxzb)
Definition: gwdps.f:342