/[MITgcm]/MITgcm/pkg/bulk_force/bulkf_formula_lanl.F
ViewVC logotype

Contents of /MITgcm/pkg/bulk_force/bulkf_formula_lanl.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.6 - (show annotations) (download)
Fri Nov 4 01:28:39 2005 UTC (18 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57y_post, checkpoint58, checkpoint57z_post, checkpoint57y_pre, checkpoint57x_post
Changes since 1.5: +3 -7 lines
remove unused variables (reduces number of compiler warning)

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

  ViewVC Help
Powered by ViewVC 1.1.22