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

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

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


Revision 1.1 - (show annotations) (download)
Thu Nov 21 19:11:42 2002 UTC (21 years, 6 months ago) by cheisey
Branch: MAIN
CVS Tags: checkpoint47b_post, checkpoint51f_post, checkpoint48i_post, checkpoint47a_post, checkpoint48d_post, checkpoint50b_post, checkpoint47i_post, checkpoint48g_post, branchpoint-genmake2, checkpoint50c_pre, checkpoint50, checkpoint51j_post, branch-exfmods-tag, checkpoint47e_post, checkpoint50f_post, checkpoint50a_post, checkpoint48e_post, checkpoint47c_post, checkpoint50f_pre, checkpoint48b_post, checkpoint47j_post, checkpoint47d_pre, checkpoint50d_pre, checkpoint47h_post, checkpoint51d_post, checkpoint51, checkpoint48c_pre, checkpoint48c_post, checkpoint50d_post, checkpoint47f_post, checkpoint51b_pre, checkpoint50e_post, checkpoint47g_post, checkpoint50h_post, checkpoint50c_post, checkpoint51a_post, checkpoint47d_post, checkpoint50e_pre, checkpoint50b_pre, checkpoint48d_pre, checkpoint50i_post, checkpoint51e_post, checkpoint51b_post, checkpoint48a_post, checkpoint51h_pre, checkpoint48f_post, checkpoint51i_pre, checkpoint50g_post, checkpoint51c_post, checkpoint51g_post, checkpoint48, checkpoint49, checkpoint51f_pre, checkpoint48h_post
Branch point for: branch-genmake2, branch-exfmods-curt
Two packages:  bulk_force (Bulk forcing)
and therm_seaice (thermodynamic_seaice) - adopted from LANL CICE.v2.0.2
Earlier integration from Stephaine Dutkiewicz
and Patrick Heimbach.

Two ifdef statements for compile time,
ALLOW_THERM_SEAICE and ALLOW_BULK_FORCE

Two switches in data.pkg to turn on at run-time:

cat data.pkg
# Packages
 &PACKAGES
 useBulkForce=.TRUE.,
 useThermSeaIce=.TRUE.,
 &

WARNING:  useSEAICE and useThermSEAICE are mutually exclusive.

The bulk package requires an additional parameter file
with two namelists, data.ice and data.blk.

c ADAPTED FROM:
c LANL CICE.v2.0.2
c-----------------------------------------------------------------------
c.. thermodynamics (vertical physics) based on M. Winton 3-layer model
c.. See Bitz, C. M. and W. H. Lipscomb, 1999:  "An energy-conserving
c..       thermodynamic sea ice model for climate study."  J. Geophys.
c..       Res., 104, 15669 - 15677.
c..     Winton, M., 1999:  "A reformulated three-layer sea ice model."
c..       Submitted to J. Atmos. Ocean. Technol.

c.. authors Elizabeth C. Hunke and William Lipscomb
c..         Fluid Dynamics Group, Los Alamos National Laboratory
c-----------------------------------------------------------------------

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_exf(
6 I uw, vw, us, Ta, Qa, nc, tsf_in,
7 I flwup, 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_BULKFORMULA
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 flwup ! 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
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 integer niter_bulk, iter
86
87 #ifdef ALLOW_BULKFORMULA
88
89 c --- external functions ---
90 _RL exf_BulkCdn
91 external exf_BulkCdn
92 _RL exf_BulkqSat
93 external exf_BulkqSat
94 _RL exf_BulkRhn
95 external exf_BulkRhn
96
97 c -- compute turbulent surface fluxes
98 ht = 2.d0
99 hq = 2.d0
100 hu = 10.d0
101 zref = 10.d0
102 aln = log(ht/zref)
103 niter_bulk = 2
104 cdalton = 0.0346000
105 czol = hu*xkar*gravity
106 psim_fac=5.d0
107 umin=5.d-1
108 c
109 Tsf=Tsf_in+Tf0kel
110 c Wind speed
111 if (us.eq.0) then
112 us = sqrt(uw*uw + vw*vw)
113 endif
114 usm = max(us,umin)
115 cQQQ try to limit drag
116 cQQ usm = min(usm,5.d0)
117 c
118 t0 = Ta*(1.d0 + humid_fac*Qa)
119 ssq = saltsat*exf_BulkqSat(Tsf)/rhoa
120 cBB debugging
121 cBB print*,'ice, ssq', qcoef, ssq
122 c
123 deltap = ta - tsf + gamma_blk*ht
124 delq = Qa - ssq
125 stable = 5.d-1 + sign(5.d-1, deltap)
126 rdn=sqrt(exf_BulkCdn(usm))
127 ustar=rdn*usm
128 tstar=exf_BulkRhn(stable)*deltap
129 qstar=cdalton*delq
130 c
131 c interation with psi-functions to find transfer coefficients
132 do iter=1,niter_bulk
133 huol = czol*(tstar/t0 +
134 & qstar/(1.d0/humid_fac+Qa))/
135 & ustar**2
136 huol = max(huol,-100.d0)
137 huol = min(huol, 100.d0)
138 stable = .5d0+ sign(.5d0, huol)
139 htol = huol*ht/hu
140 hqol = huol*hq/hu
141
142 c Evaluate all stability functions assuming hq = ht.
143 xsq = max(sqrt(abs(1.d0 - 16.*huol)),1.d0)
144 x = sqrt(xsq)
145 psimh = -psim_fac*huol*stable +
146 & (1.d0 - stable)*
147 & log((1.d0 + x*(2.d0 + x))*
148 & (1.d0 + xsq)/8.d0) - 2.d0*atan(x) +
149 & pi*.5d0
150 xsq = max(sqrt(abs(1.d0 - 16.*htol)),1.d0)
151 psixh = -psim_fac*htol*stable + (1.d0 - stable)*
152 & 2.d0*log((1.d0 + xsq)/2.d0)
153
154 c Shift wind speed using old coefficient
155 ccc rd = rdn/(1.d0 + rdn/xkar*
156 ccc & (log(hu/zref) - psimh) )
157 rd = rdn/(1.d0 - rdn/xkar*psimh )
158 shn = us*rd/rdn
159 uzn = max(shn, umin)
160
161 c Update the transfer coefficients at 10 meters and neutral stability.
162
163 rdn = sqrt(exf_BulkCdn(uzn))
164
165 c Shift all coefficients to the measurement height and stability.
166 c rd = rdn/(1.d0 + rdn/xkar*(log(hu/zref) - psimh))
167 rd = rdn/(1.d0 - rdn/xkar*psimh)
168 rh = exf_BulkRhn(stable)/(1.d0 +
169 & exf_BulkRhn(stable)/
170 & xkar*(aln - psixh))
171 re = cdalton/(1.d0 + cdalton/xkar*(aln - psixh))
172
173 c Update ustar, tstar, qstar using updated, shifted coefficients.
174 ustar = rd*usm
175 qstar = re*delq
176 tstar = rh*deltap
177 enddo
178 c
179 tau = rhoa*ustar**2
180 tau = tau*us/usm
181 fsha = cpair*tau*tstar/ustar
182 flha = Lvap *tau*qstar/ustar
183 evp = tau*qstar/ustar
184 c
185
186 c up long wave radiation
187 mixratio=Qa/(1-Qa)
188 ea=p0*mixratio/(0.62197+mixratio)
189 flwup=-0.985*stefan*tsf**4
190 & *(0.39-0.05*sqrt(ea))
191 & *(1-0.6*nc**2)
192 #ifdef ALLOW_TSEAICE
193 if (iceornot.eq.1) then
194 c derivatives with respect to Tsf
195 dflhdT = -(Lvap*tau*re*ssq/ustar)/(tsf**2)
196 dfshdT = -cpair*tau*rh/ustar
197 dflwupdT=-4.*0.985*stefan*tsf**3
198 & *(0.39-0.05*sqrt(ea))
199 c total derivative with respect to surface temperature
200 df0dT=dflwupdT+dfshdT+dflhdT
201 cBB
202 cBB print*,'derivatives:',df0dT,dflwupdT,dfshdT,dflhdT
203 endif
204 #endif /*ALLOW_TSEAICE*/
205 c
206 c wind stress at center points
207 if (.NOT.windread) then
208 ust = rhoa*exf_BulkCdn(usm)*us*uw
209 vst = rhoa*exf_BulkCdn(usm)*us*vw
210 endif
211 #endif /*ALLOW_BULKFORMULA*/
212 ssq=deltap
213
214 return
215 end

  ViewVC Help
Powered by ViewVC 1.1.22