Line data Source code
1 : module beljaars_drag
2 :
3 : implicit none
4 : private
5 : save
6 :
7 : public init_blj ! Initialization
8 : public compute_blj ! Full routine
9 :
10 : ! ------------ !
11 : ! Private data !
12 : ! ------------ !
13 :
14 : integer, parameter :: r8 = selected_real_kind(12) ! 8 byte real
15 :
16 : real(r8), parameter :: horomin= 1._r8 ! Minimum value of subgrid orographic height for mountain stress [ m ]
17 : real(r8), parameter :: z0max = 100._r8 ! Maximum value of z_0 for orography [ m ]
18 : real(r8), parameter :: dv2min = 0.01_r8 ! Minimum shear squared [ m2/s2 ]
19 : real(r8) :: orocnst ! Converts from standard deviation to height [ no unit ]
20 : real(r8) :: z0fac ! Factor determining z_0 from orographic standard deviation [ no unit ]
21 : real(r8) :: karman ! von Karman constant
22 : real(r8) :: gravit ! Acceleration due to gravity
23 : real(r8) :: rair ! Gas constant for dry air
24 :
25 : contains
26 :
27 : !============================================================================ !
28 : ! !
29 : !============================================================================ !
30 :
31 1536 : subroutine init_blj( kind, gravit_in, rair_in , errstring )
32 :
33 : integer, intent(in) :: kind
34 : real(r8), intent(in) :: gravit_in, rair_in
35 :
36 : character(len=*), intent(out) :: errstring
37 :
38 1536 : errstring = ' '
39 :
40 1536 : if ( kind /= r8 ) then
41 0 : errstring = 'inconsistent KIND of reals passed to init_blj'
42 0 : return
43 : endif
44 :
45 1536 : gravit = gravit_in
46 1536 : rair = rair_in
47 :
48 1536 : end subroutine init_blj
49 :
50 : !============================================================================ !
51 : ! !
52 : !============================================================================ !
53 :
54 1489176 : subroutine compute_blj( pcols , pver , ncol , &
55 1489176 : u , v , t , pmid , delp , &
56 1489176 : zm , sgh , drag , taux , tauy , &
57 : landfrac )
58 :
59 : !------------------------------------------------------------------------------ !
60 : ! Beljaars Sub-Grid Orographic (SGO) Form drag parameterization !
61 : ! !
62 : ! Returns drag profile and integrated stress associated with subgrid mountains !
63 : ! with horizontal length scales nominally below 3km. Similar to TMS but !
64 : ! drag is distributed in the vertical (Beljaars et al., 2003, QJRMS). !
65 : ! !
66 : ! First cut follows TMS. J. Bacmeister, March 2016 !
67 : !------------------------------------------------------------------------------ !
68 :
69 : ! ---------------------- !
70 : ! Input-Output Arguments !
71 : ! ---------------------- !
72 :
73 : integer, intent(in) :: pcols ! Number of columns dimensioned
74 : integer, intent(in) :: pver ! Number of model layers
75 : integer, intent(in) :: ncol ! Number of columns actually used
76 :
77 : real(r8), intent(in) :: u(pcols,pver) ! Layer mid-point zonal wind [ m/s ]
78 : real(r8), intent(in) :: v(pcols,pver) ! Layer mid-point meridional wind [ m/s ]
79 : real(r8), intent(in) :: t(pcols,pver) ! Layer mid-point temperature [ K ]
80 : real(r8), intent(in) :: pmid(pcols,pver) ! Layer mid-point pressure [ Pa ]
81 : real(r8), intent(in) :: delp(pcols,pver) ! Layer thickness [ Pa ]
82 : real(r8), intent(in) :: zm(pcols,pver) ! Layer mid-point height [ m ]
83 : real(r8), intent(in) :: sgh(pcols) ! Standard deviation of orography [ m ]
84 : real(r8), intent(in) :: landfrac(pcols) ! Land fraction [ fraction ]
85 :
86 : real(r8), intent(out) :: drag(pcols,pver) ! SGO drag profile [ kg/s/m2 ]
87 : real(r8), intent(out) :: taux(pcols) ! Surface zonal wind stress [ N/m2 ]
88 : real(r8), intent(out) :: tauy(pcols) ! Surface meridional wind stress [ N/m2 ]
89 :
90 : ! --------------- !
91 : ! Local Variables !
92 : ! --------------- !
93 :
94 : integer :: i,k ! Loop indices
95 : integer :: kb, kt ! Bottom and top of source region
96 :
97 : real(r8) :: vmag ! Velocity magnitude [ m /s ]
98 :
99 : real(r8) :: alpha,beta,Cmd,Ccorr,n1,n2,k1,kflt,k2,IH
100 2978352 : real(r8) :: a1(pcols),a2(pcols)
101 :
102 1489176 : alpha = 12._r8
103 1489176 : beta = 1._r8
104 1489176 : n1 = -1.9_r8
105 1489176 : n2 = -2.8_r8
106 :
107 1489176 : Cmd = 0.005_r8
108 1489176 : Ccorr = 0.6_r8 * 5._r8
109 :
110 1489176 : kflt = 0.00035_r8 ! m-1
111 1489176 : k1 = 0.003_r8 ! m-1
112 1489176 : IH = 0.00102_r8 ! m-1
113 :
114 24865776 : a1(1:ncol) = (sgh(1:ncol)*sgh(1:ncol)) / ( IH* (kflt**n1) )
115 24865776 : a2(1:ncol) = a1(1:ncol) * k1**(n1-n2)
116 :
117 :
118 : ! ----------------------- !
119 : ! Main Computation Begins !
120 : ! ----------------------- !
121 :
122 139982544 : do k = 1, pver
123 2314006344 : do i = 1, ncol
124 2174023800 : Vmag = SQRT( u(i,k)**2 + v(i,k)**2)
125 : drag(i,k) = -alpha * beta * Cmd * Ccorr * Vmag * 2.109_r8 * &
126 : EXP ( -(zm(i,k)/1500._r8 )*SQRT(zm(i,k)/1500._r8) ) * ( zm(i,k)**(-1.2_r8) ) &
127 2312517168 : * a2(i)
128 : end do
129 : end do
130 :
131 :
132 : !---------------------------------!
133 : ! Diagnose effective surface drag !
134 : ! in X and Y by integrating in !
135 : ! the vertical !
136 : !---------------------------------!
137 : ! FIXME: uses 'state' u and v.
138 : ! Should updated u and v's be used?
139 :
140 25315992 : taux=0._r8
141 25315992 : tauy=0._r8
142 139982544 : do k = 1, pver
143 2314006344 : do i = 1, ncol
144 2174023800 : taux(i) = taux(i) + drag(i,k)*u(i,k)*delp(i,k)/gravit
145 2312517168 : tauy(i) = tauy(i) + drag(i,k)*v(i,k)*delp(i,k)/gravit
146 : end do
147 : end do
148 :
149 1489176 : return
150 : end subroutine compute_blj
151 :
152 : end module beljaars_drag
|