1 |
c swd -- bulkf formula used in bulkf and ice pkgs |
2 |
|
3 |
c taken from exf package |
4 |
#include "CPP_OPTIONS.h" |
5 |
subroutine bulkf_formula_lanl( |
6 |
I uw, vw, us, Ta, Qa, nc, tsf_in, |
7 |
I flwupa, flha, fsha, df0dT, |
8 |
I ust, vst, evp, ssq, iceornot, windread |
9 |
& ) |
10 |
|
11 |
IMPLICIT NONE |
12 |
c |
13 |
#include "EEPARAMS.h" |
14 |
#include "SIZE.h" |
15 |
#include "PARAMS.h" |
16 |
#include "DYNVARS.h" |
17 |
#include "GRID.h" |
18 |
#include "FFIELDS.h" |
19 |
#ifdef ALLOW_BULK_FORCE |
20 |
#include "BULKF_ICE_CONSTANTS.h" |
21 |
#include "BULKF.h" |
22 |
#endif |
23 |
|
24 |
c |
25 |
c calculate bulk formula fluxes over open ocean or seaice |
26 |
c |
27 |
c input |
28 |
_RL us ! wind speed |
29 |
_RL uw ! zonal wind speed (at grid center) |
30 |
_RL vw ! meridional wind speed (at grid center) |
31 |
_RL Ta ! air temperature at ht |
32 |
_RL Qa ! specific humidity at ht |
33 |
_RL tsf_in ! surface temperature (either ice or open water) |
34 |
_RL nc ! fraction cloud cover |
35 |
integer iceornot ! 0=open water, 1=ice cover |
36 |
logical windread ! |
37 |
c output |
38 |
_RL flwupa ! upward long wave radiation |
39 |
_RL flha ! latent heat flux |
40 |
_RL fsha ! sensible heat flux |
41 |
_RL df0dT ! derivative of heat flux with respect to Tsf |
42 |
_RL ust ! zonal wind stress (at grid center) |
43 |
_RL vst ! meridional wind stress (at grid center) |
44 |
_RL evp ! evaporation rate (over open water) |
45 |
_RL ssq ! surface specific humidity (kg/kg) |
46 |
c |
47 |
c local variables |
48 |
_RL tsf ! surface temperature in K |
49 |
_RL ht ! reference height (m) |
50 |
_RL hq ! reference height for humidity (m) |
51 |
_RL hu ! reference height for wind speed (m) |
52 |
_RL zref ! reference height |
53 |
_RL czol |
54 |
_RL usm ! wind speed limited |
55 |
_RL t0 ! virtual temperature (K) |
56 |
_RL deltap ! potential temperature diff (K) |
57 |
_RL delq ! specific humidity difference |
58 |
_RL stable |
59 |
_RL rdn ,ren, rhn |
60 |
_RL ustar |
61 |
_RL tstar |
62 |
_RL qstar |
63 |
_RL huol |
64 |
_RL htol |
65 |
_RL hqol |
66 |
_RL xsq |
67 |
_RL x |
68 |
_RL re |
69 |
_RL rh |
70 |
_RL tau |
71 |
_RL psimh |
72 |
_RL psixh |
73 |
_RL rd |
74 |
_RL uzn |
75 |
_RL shn |
76 |
_RL aln |
77 |
_RL cdalton |
78 |
_RL dflhdT ! derivative of latent heat with respect to T |
79 |
_RL dfshdT ! derivative of sensible heat with respect to T |
80 |
_RL dflwupdT ! derivative of long wave with respect to T |
81 |
_RL mixratio |
82 |
_RL ea |
83 |
_RL psim_fac |
84 |
_RL umin |
85 |
_RL lath |
86 |
_RL csha |
87 |
_RL clha |
88 |
_RL zice |
89 |
integer niter_bulk, iter |
90 |
|
91 |
#ifdef ALLOW_BULK_FORCE |
92 |
|
93 |
c --- external functions --- |
94 |
_RL exf_BulkCdn |
95 |
external exf_BulkCdn |
96 |
_RL exf_BulkqSat |
97 |
external exf_BulkqSat |
98 |
_RL exf_BulkRhn |
99 |
external exf_BulkRhn |
100 |
|
101 |
cQQQQQQQQQQQQQQQQQQQQq |
102 |
c -- compute turbulent surface fluxes |
103 |
ht = 2.d0 |
104 |
hq = 2.d0 |
105 |
hu = 10.d0 |
106 |
zref = 10.d0 |
107 |
zice = 0.0005 |
108 |
aln = log(ht/zref) |
109 |
niter_bulk = 5 |
110 |
cdalton = 0.0346000 |
111 |
czol = zref*xkar*gravity |
112 |
psim_fac=5.d0 |
113 |
umin=1.d0 |
114 |
c |
115 |
lath=Lvap |
116 |
if (iceornot.gt.0) lath=Lvap+Lfresh |
117 |
Tsf=Tsf_in+Tf0kel |
118 |
c Wind speed |
119 |
if (us.eq.0) then |
120 |
us = sqrt(uw*uw + vw*vw) |
121 |
endif |
122 |
usm = max(us,umin) |
123 |
cQQQ try to limit drag |
124 |
cQQ usm = min(usm,5.d0) |
125 |
c |
126 |
t0 = Ta*(1.d0 + humid_fac*Qa) |
127 |
cQQ ssq = 0.622*6.11*exp(22.47*(1.d0-Tf0kel/tsf))/p0 |
128 |
ssq =3.797915*exp(lath*(7.93252d-6-2.166847d-3/Tsf))/p0 |
129 |
cBB debugging |
130 |
cBB print*,'ice, ssq', ssq |
131 |
c |
132 |
deltap = ta - tsf + gamma_blk*ht |
133 |
delq = Qa - ssq |
134 |
c |
135 |
c initialize estimate exchange coefficients |
136 |
rdn=xkar/(log(zref/zice)) |
137 |
rhn=rdn |
138 |
ren=rdn |
139 |
c calculate turbulent scales |
140 |
ustar=rdn*usm |
141 |
tstar=rhn*deltap |
142 |
qstar=ren*delq |
143 |
c |
144 |
c interation with psi-functions to find transfer coefficients |
145 |
do iter=1,niter_bulk |
146 |
huol = czol/ustar**2 *(tstar/t0 + |
147 |
& qstar/(1.d0/humid_fac+Qa)) |
148 |
huol = sign( min(abs(huol),10.d0), huol) |
149 |
stable = 5. _d -1 + sign(5. _d -1 , huol) |
150 |
xsq = max(sqrt(abs(1.d0 - 16.d0*huol)),1.d0) |
151 |
x = sqrt(xsq) |
152 |
psimh = -5.d0*huol*stable + (1.d0-stable)* |
153 |
& (2.d0*log(5.e-1*(1.d0+x)) + |
154 |
& 2.d0*log(5.e-1*(1.d0+xsq)) - |
155 |
& 2.d0*atan(x) + pi*.5d0) |
156 |
psixh = -5.d0*huol*stable + (1.d0-stable)* |
157 |
& (2.d0*log(5.e-1*(1.d0+xsq))) |
158 |
|
159 |
c Update the transfer coefficients |
160 |
|
161 |
rd = rdn/(1.d0 + rdn*(aln-psimh)/xkar) |
162 |
rh = rhn/(1.d0 + rhn*(aln-psixh)/xkar) |
163 |
re = rh |
164 |
c Update ustar, tstar, qstar using updated, shifted coefficients. |
165 |
ustar = rd*usm |
166 |
qstar = re*delq |
167 |
tstar = rh*deltap |
168 |
enddo |
169 |
c |
170 |
tau = rhoa*ustar**2 |
171 |
tau = tau*us/usm |
172 |
csha = rhoa*cpair*us*rh*rd |
173 |
clha = rhoa*lath*us*re*rd |
174 |
c |
175 |
fsha = csha*deltap |
176 |
flha = clha*delq |
177 |
evp = flha/lath/rhosw |
178 |
c |
179 |
c up long wave radiation |
180 |
cQQ mixratio=Qa/(1-Qa) |
181 |
cQQ ea=p0*mixratio/(0.62197+mixratio) |
182 |
cQQ flwupa=-0.985*stefan*tsf**4 |
183 |
cQQ & *(0.39-0.05*sqrt(ea)) |
184 |
cQQ & *(1-0.6*nc**2) |
185 |
if (iceornot.eq.0) then |
186 |
flwupa=-ocean_emissivity*stefan*tsf**4 |
187 |
else |
188 |
if (iceornot.eq.2) then |
189 |
flwupa=-snow_emissivity*stefan*tsf**4 |
190 |
else |
191 |
flwupa=-ice_emissivity*stefan*tsf**4 |
192 |
endif |
193 |
endif |
194 |
#ifdef ALLOW_THERM_SEAICE |
195 |
cswdcou -- remove --- |
196 |
cQQ if (iceornot.gt.0) then |
197 |
cswdcou --- end remove --- |
198 |
c derivatives with respect to Tsf |
199 |
cQQ dflhdT = -clha*Tf0kel*ssq*22.47/(tsf**2) |
200 |
c dflhdT = -clha*Lath*ssq/(Rvap*tsf**2) |
201 |
dflhdT = -clha*ssq*Lath*2.166847E-3/(Tsf**2) |
202 |
dfshdT = -csha |
203 |
if (iceornot.eq.0) then |
204 |
dflwupdT=-4.*ocean_emissivity*stefan*tsf**3 |
205 |
else |
206 |
if (iceornot.eq.2) then |
207 |
dflwupdT=-4.*snow_emissivity*stefan*tsf**3 |
208 |
else |
209 |
dflwupdT=-4.*ice_emissivity*stefan*tsf**3 |
210 |
endif |
211 |
endif |
212 |
cQQ dflwupdT=-4.*0.985*stefan*tsf**3 |
213 |
cQQ & *(0.39-0.05*sqrt(ea)) |
214 |
cQQ & *(1-0.6*nc**2) |
215 |
c total derivative with respect to surface temperature |
216 |
df0dT=dflwupdT+dfshdT+dflhdT |
217 |
cBB |
218 |
cBB print*,'derivatives:',df0dT,dflwupdT,dfshdT,dflhdT |
219 |
cswdcou -- remove --- |
220 |
cQQ endif |
221 |
cswdcou |
222 |
#endif /*ALLOW_THERM_SEAICE*/ |
223 |
c |
224 |
c wind stress at center points |
225 |
if (.NOT.windread) then |
226 |
ust = rhoa*exf_BulkCdn(usm)*us*uw |
227 |
vst = rhoa*exf_BulkCdn(usm)*us*vw |
228 |
endif |
229 |
#endif /*ALLOW_BULK_FORCE*/ |
230 |
|
231 |
return |
232 |
end |