/[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.4 - (show annotations) (download)
Sun Nov 23 01:36:55 2003 UTC (20 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint52l_pre, hrcube4, hrcube5, checkpoint52d_pre, checkpoint52j_pre, checkpoint52k_post, checkpoint52f_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint52e_pre, checkpoint52e_post, checkpoint52b_post, checkpoint52c_post, checkpoint52f_pre, checkpoint52d_post, checkpoint52i_post, checkpoint52h_pre, checkpoint52j_post, branch-netcdf, checkpoint52l_post
Branch point for: netcdf-sm0
Changes since 1.3: +59 -77 lines
light cleaning to make it work with new thSIce pkg.

1 C $Header: $
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,
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)
43 _RL ssq ! surface specific humidity (kg/kg)
44
45 #ifdef ALLOW_BULK_FORCE
46
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 c --- external functions ---
92 _RL exf_BulkCdn
93 external exf_BulkCdn
94 c _RL exf_BulkqSat
95 c external exf_BulkqSat
96 c _RL exf_BulkRhn
97 c external exf_BulkRhn
98
99 cQQQQQQQQQQQQQQQQQQQQq
100 c -- compute turbulent surface fluxes
101 ht = 2. _d 0
102 hq = 2. _d 0
103 hu = 10. _d 0
104 zref = 10. _d 0
105 zice = 0.0005 _d 0
106 aln = log(ht/zref)
107 niter_bulk = 5
108 cdalton = 0.0346000 _d 0
109 czol = zref*xkar*gravity
110 psim_fac=5. _d 0
111 umin=1. _d 0
112 c
113 lath=Lvap
114 if (iceornot.gt.0) lath=Lvap+Lfresh
115 Tsf=Tsf_in+Tf0kel
116 c Wind speed
117 if (us.eq.0. _d 0) then
118 us = sqrt(uw*uw + vw*vw)
119 endif
120 usm = max(us,umin)
121 cQQQ try to limit drag
122 cQQ usm = min(usm,5. _d 0)
123 c
124 t0 = Ta*(1. _d 0 + humid_fac*Qa)
125 cQQ ssq = 0.622*6.11*exp(22.47*(1.d0-Tf0kel/tsf))/p0
126 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
127 ssq = 3.797915 _d 0*exp(
128 & lath*(7.93252 _d -6 - 2.166847 _d -3/Tsf)
129 & )/p0
130 cBB debugging
131 cBB print*,'ice, ssq', ssq
132 c
133 deltap = ta - tsf + gamma_blk*ht
134 delq = Qa - ssq
135 c
136 c initialize estimate exchange coefficients
137 rdn=xkar/(log(zref/zice))
138 rhn=rdn
139 ren=rdn
140 c calculate turbulent scales
141 ustar=rdn*usm
142 tstar=rhn*deltap
143 qstar=ren*delq
144 c
145 c interation with psi-functions to find transfer coefficients
146 do iter=1,niter_bulk
147 huol = czol/ustar**2 *(tstar/t0 +
148 & qstar/(1. _d 0/humid_fac+Qa))
149 huol = sign( min(abs(huol),10. _d 0), huol)
150 stable = 5. _d -1 + sign(5. _d -1 , huol)
151 xsq = max(sqrt(abs(1. _d 0 - 16. _d 0*huol)),1. _d 0)
152 x = sqrt(xsq)
153 psimh = -5. _d 0*huol*stable + (1. _d 0-stable)*
154 & (2. _d 0*log(5. _d -1*(1. _d 0+x)) +
155 & 2. _d 0*log(5. _d -1*(1. _d 0+xsq)) -
156 & 2. _d 0*atan(x) + pi*.5 _d 0)
157 psixh = -5. _d 0*huol*stable + (1. _d 0-stable)*
158 & (2. _d 0*log(5. _d -1*(1. _d 0+xsq)))
159
160 c Update the transfer coefficients
161
162 rd = rdn/(1. _d 0 + rdn*(aln-psimh)/xkar)
163 rh = rhn/(1. _d 0 + rhn*(aln-psixh)/xkar)
164 re = rh
165 c Update ustar, tstar, qstar using updated, shifted coefficients.
166 ustar = rd*usm
167 qstar = re*delq
168 tstar = rh*deltap
169 enddo
170 c
171 tau = rhoa*ustar**2
172 tau = tau*us/usm
173 csha = rhoa*cpair*us*rh*rd
174 clha = rhoa*lath*us*re*rd
175 c
176 fsha = csha*deltap
177 flha = clha*delq
178 evp = -flha/lath/rhofw
179 c
180 c up long wave radiation
181 cQQ mixratio=Qa/(1-Qa)
182 cQQ ea=p0*mixratio/(0.62197+mixratio)
183 cQQ flwupa=-0.985*stefan*tsf**4
184 cQQ & *(0.39-0.05*sqrt(ea))
185 cQQ & *(1-0.6*nc**2)
186 if (iceornot.eq.0) then
187 flwupa=ocean_emissivity*stefan*tsf**4
188 dflwupdT=4. _d 0*ocean_emissivity*stefan*tsf**3
189 else
190 if (iceornot.eq.2) then
191 flwupa=snow_emissivity*stefan*tsf**4
192 dflwupdT=4. _d 0*snow_emissivity*stefan*tsf**3
193 else
194 flwupa=ice_emissivity*stefan*tsf**4
195 dflwupdT=4. _d 0*ice_emissivity*stefan*tsf**3
196 endif
197 endif
198 cQQ dflhdT = -clha*Tf0kel*ssq*22.47/(tsf**2)
199 c dflhdT = -clha*Lath*ssq/(Rvap*tsf**2)
200 dflhdT = -clha*ssq*Lath*2.166847 _d -3/(Tsf**2)
201 dfshdT = -csha
202 cQQ dflwupdT= 4.*0.985*stefan*tsf**3
203 cQQ & *(0.39-0.05*sqrt(ea))
204 cQQ & *(1-0.6*nc**2)
205 c total derivative with respect to surface temperature
206 df0dT=-dflwupdT+dfshdT+dflhdT
207 c
208 c wind stress at center points
209 if (.NOT.windread) then
210 ust = rhoa*exf_BulkCdn(usm)*us*uw
211 vst = rhoa*exf_BulkCdn(usm)*us*vw
212 endif
213 #endif /*ALLOW_BULK_FORCE*/
214
215 RETURN
216 END

  ViewVC Help
Powered by ViewVC 1.1.22