1 |
C $Header: /u/gcmpack/MITgcm/pkg/bulk_force/bulkf_forcing.F,v 1.13 2009/04/28 18:08:13 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "BULK_FORCE_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: BULKF_FORCING |
8 |
C !INTERFACE: |
9 |
SUBROUTINE BULKF_FORCING( |
10 |
I myTime, myIter, myThid ) |
11 |
|
12 |
C !DESCRIPTION: \bv |
13 |
C *==========================================================* |
14 |
C | SUBROUTINE BULKF_FORCING |
15 |
C *==========================================================* |
16 |
C \ev |
17 |
|
18 |
C o Get the surface fluxes used to force ocean model |
19 |
C Output: |
20 |
C ------ |
21 |
C ustress, vstress :: wind stress |
22 |
C qnet :: net heat flux |
23 |
C empmr :: freshwater flux |
24 |
C --------- |
25 |
C |
26 |
C Input: |
27 |
C ------ |
28 |
C uwind, vwind :: mean wind speed (m/s) at height hu (m) |
29 |
C Tair :: mean air temperature (K) at height ht (m) |
30 |
C Qair :: mean air humidity (kg/kg) at height hq (m) |
31 |
C theta(k=1) :: sea surface temperature (K) |
32 |
C rain :: precipitation |
33 |
C runoff :: river(ice) runoff |
34 |
C |
35 |
C ================================================================== |
36 |
C SUBROUTINE bulkf_forcing |
37 |
C ================================================================== |
38 |
|
39 |
C !USES: |
40 |
IMPLICIT NONE |
41 |
C == global variables == |
42 |
#include "EEPARAMS.h" |
43 |
#include "SIZE.h" |
44 |
#include "PARAMS.h" |
45 |
#include "DYNVARS.h" |
46 |
#include "GRID.h" |
47 |
#include "FFIELDS.h" |
48 |
#include "BULKF_PARAMS.h" |
49 |
#include "BULKF.h" |
50 |
#include "BULKF_INT.h" |
51 |
#include "BULKF_TAVE.h" |
52 |
|
53 |
C !INPUT/OUTPUT PARAMETERS: |
54 |
C == routine arguments == |
55 |
_RL myTime |
56 |
INTEGER myIter |
57 |
INTEGER myThid |
58 |
CEOP |
59 |
|
60 |
#ifdef ALLOW_BULK_FORCE |
61 |
C == Local variables == |
62 |
INTEGER bi,bj |
63 |
INTEGER i,j |
64 |
INTEGER ks, iceornot |
65 |
|
66 |
_RL df0dT, hfl, evp, dEvdT |
67 |
#ifdef ALLOW_FORMULA_AIM |
68 |
_RL SHF(1), EVPloc(1), SLRU(1) |
69 |
_RL dEvp(1), sFlx(0:2) |
70 |
#endif |
71 |
|
72 |
C- surface level index: |
73 |
ks = 1 |
74 |
|
75 |
C- Compute surface fluxes over ice-free ocean only: |
76 |
iceornot = 0 |
77 |
|
78 |
DO bj=myByLo(myThid),myByHi(myThid) |
79 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
80 |
|
81 |
DO j = 1-OLy,sNy+OLy |
82 |
DO i = 1-OLx,sNx+OLx |
83 |
IF ( maskC(i,j,ks,bi,bj).NE.0. _d 0 ) THEN |
84 |
|
85 |
#ifdef ALLOW_FORMULA_AIM |
86 |
IF ( useFluxFormula_AIM ) THEN |
87 |
CALL BULKF_FORMULA_AIM( |
88 |
I theta(i,j,ks,bi,bj), flwdwn(i,j,bi,bj), |
89 |
I thAir(i,j,bi,bj), Tair(i,j,bi,bj), |
90 |
I Qair(i,j,bi,bj), wspeed(i,j,bi,bj), |
91 |
O SHF, EVPloc, SLRU, |
92 |
O dEvp, sFlx, |
93 |
I iceornot, myThid ) |
94 |
|
95 |
flwup(i,j,bi,bj)= ocean_emissivity*SLRU(1) |
96 |
C- reverse sign (AIM convention -> BULKF convention): |
97 |
fsh(i,j,bi,bj) = -SHF(1) |
98 |
flh(i,j,bi,bj) = -Lvap*EVPloc(1) |
99 |
C- Convert from g/m2/s to m/s |
100 |
evap(i,j,bi,bj) = EVPloc(1) * 1. _d -3 / rhoFW |
101 |
dEvdT = dEvp(1) * 1. _d -3 |
102 |
df0dT = sFlx(2) |
103 |
|
104 |
ELSEIF ( blk_nIter.EQ.0 ) THEN |
105 |
#else /* ALLOW_FORMULA_AIM */ |
106 |
IF ( blk_nIter.EQ.0 ) THEN |
107 |
#endif /* ALLOW_FORMULA_AIM */ |
108 |
CALL BULKF_FORMULA_LANL( |
109 |
I uwind(i,j,bi,bj),vwind(i,j,bi,bj),wspeed(i,j,bi,bj), |
110 |
I Tair(i,j,bi,bj), Qair(i,j,bi,bj), |
111 |
I cloud(i,j,bi,bj),theta(i,j,ks,bi,bj), |
112 |
O flwup(i,j,bi,bj), flh(i,j,bi,bj), |
113 |
O fsh(i,j,bi,bj), df0dT, |
114 |
O ustress(i,j,bi,bj), vstress(i,j,bi,bj), |
115 |
O evp, savssq(i,j,bi,bj), dEvdT, |
116 |
I iceornot, myThid ) |
117 |
C Note that the LANL flux conventions are opposite |
118 |
C of what they are in the model. |
119 |
|
120 |
C- Convert from kg/m2/s to m/s |
121 |
evap(i,j,bi,bj) = evp/rhoFW |
122 |
|
123 |
ELSE |
124 |
CALL BULKF_FORMULA_LAY( |
125 |
I uwind(i,j,bi,bj), vwind(i,j,bi,bj), |
126 |
I wspeed(i,j,bi,bj), Tair(i,j,bi,bj), |
127 |
I Qair(i,j,bi,bj), theta(i,j,ks,bi,bj), |
128 |
O flwup(i,j,bi,bj), flh(i,j,bi,bj), |
129 |
O fsh(i,j,bi,bj), df0dT, |
130 |
O ustress(i,j,bi,bj), vstress(i,j,bi,bj), |
131 |
O evp, savssq(i,j,bi,bj), dEvdT, |
132 |
I iceornot, i,j,bi,bj,myThid ) |
133 |
|
134 |
C- Convert from kg/m2/s to m/s |
135 |
evap(i,j,bi,bj) = evp/rhoFW |
136 |
|
137 |
ENDIF |
138 |
|
139 |
C- use down long wave data |
140 |
flwupnet(i,j,bi,bj)=flwup(i,j,bi,bj)-flwdwn(i,j,bi,bj) |
141 |
C- using down solar, need to have water albedo -- .1 |
142 |
fswnet(i,j,bi,bj) = solar(i,j,bi,bj) |
143 |
& *( 1. _d 0 - ocean_albedo ) |
144 |
ElSE |
145 |
ustress(i,j,bi,bj) = 0. _d 0 |
146 |
vstress(i,j,bi,bj) = 0. _d 0 |
147 |
fsh(i,j,bi,bj) = 0. _d 0 |
148 |
flh(i,j,bi,bj) = 0. _d 0 |
149 |
flwup(i,j,bi,bj) = 0. _d 0 |
150 |
evap(i,j,bi,bj) = 0. _d 0 |
151 |
fswnet(i,j,bi,bj) = 0. _d 0 |
152 |
savssq(i,j,bi,bj) = 0. _d 0 |
153 |
ENDIF |
154 |
ENDDO |
155 |
ENDDO |
156 |
|
157 |
IF ( calcWindStress ) THEN |
158 |
C- move wind stresses to u and v points |
159 |
DO j = 1-OLy,sNy+OLy |
160 |
DO i = 1-OLx+1,sNx+OLx |
161 |
fu(i,j,bi,bj) = maskW(i,j,1,bi,bj) |
162 |
& *(ustress(i,j,bi,bj)+ustress(i-1,j,bi,bj))*0.5 _d 0 |
163 |
ENDDO |
164 |
ENDDO |
165 |
DO j = 1-OLy+1,sNy+OLy |
166 |
DO i = 1-OLx,sNx+OLx |
167 |
fv(i,j,bi,bj) = maskS(i,j,1,bi,bj) |
168 |
& *(vstress(i,j,bi,bj)+vstress(i,j-1,bi,bj))*0.5 _d 0 |
169 |
ENDDO |
170 |
ENDDO |
171 |
ENDIF |
172 |
|
173 |
C- Add all contributions. |
174 |
DO j = 1-OLy,sNy+OLy |
175 |
DO i = 1-OLx,sNx+OLx |
176 |
IF ( maskC(i,j,ks,bi,bj).NE.0. _d 0 ) THEN |
177 |
C- Net downward surface heat flux : |
178 |
hfl = 0. _d 0 |
179 |
hfl = hfl + fsh(i,j,bi,bj) |
180 |
hfl = hfl + flh(i,j,bi,bj) |
181 |
hfl = hfl - flwupnet(i,j,bi,bj) |
182 |
hfl = hfl + fswnet(i,j,bi,bj) |
183 |
C- Heat flux: |
184 |
Qnet(i,j,bi,bj) = -hfl |
185 |
Qsw (i,j,bi,bj) = -fswnet(i,j,bi,bj) |
186 |
#ifdef COUPLE_MODEL |
187 |
dFdT(i,j,bi,bj) = df0dT |
188 |
#endif |
189 |
C- Fresh-water flux from Precipitation and Evaporation. |
190 |
EmPmR(i,j,bi,bj) = (evap(i,j,bi,bj)-rain(i,j,bi,bj) |
191 |
& - runoff(i,j,bi,bj))*rhoConstFresh |
192 |
C---- cheating: now done in S/R BULKF_FLUX_ADJUST, over ice-free ocean only |
193 |
c Qnet(i,j,bi,bj) = Qnetch(i,j,bi,bj) |
194 |
c EmPmR(i,j,bi,bj) = EmPch(i,j,bi,bj) |
195 |
C---- |
196 |
ELSE |
197 |
Qnet(i,j,bi,bj) = 0. _d 0 |
198 |
Qsw (i,j,bi,bj) = 0. _d 0 |
199 |
EmPmR(i,j,bi,bj)= 0. _d 0 |
200 |
#ifdef COUPLE_MODEL |
201 |
dFdT(i,j,bi,bj) = 0. _d 0 |
202 |
#endif |
203 |
ENDIF |
204 |
ENDDO |
205 |
ENDDO |
206 |
|
207 |
IF ( temp_EvPrRn .NE. UNSET_RL ) THEN |
208 |
C-- Account for energy content of Precip + RunOff & Evap. Assumes: |
209 |
C 1) Rain has same temp as Air |
210 |
C 2) Snow has no heat capacity (consistent with seaice & thsice pkgs) |
211 |
C 3) No distinction between sea-water Cp and fresh-water Cp |
212 |
C 4) Run-Off comes at the temp of surface water (with same Cp) |
213 |
C 5) Evap is released to the Atmos @ surf-temp (=SST); should be using |
214 |
C the water-vapor heat capacity here and consistently in Bulk-Formulae; |
215 |
C Could also be put directly into Latent Heat flux. |
216 |
c IF ( SnowFile .NE. ' ' ) THEN |
217 |
C-- Melt snow (if provided) into the ocean and account for rain-temp |
218 |
c DO j = 1-OLy,sNy+OLy |
219 |
c DO i = 1-OLx,sNx+OLx |
220 |
c Qnet(i,j,bi,bj) = Qnet(i,j,bi,bj) |
221 |
c & + Lfresh*snowPrecip(i,j,bi,bj)*rhoConstFresh |
222 |
c & - HeatCapacity_Cp |
223 |
c & *( Tair(i,j,bi,bj) - Tf0kel - temp_EvPrRn ) |
224 |
c & *( rain(i,j,bi,bj)- snowPrecip(i,j,bi,bj) ) |
225 |
c & *rhoConstFresh |
226 |
c ENDDO |
227 |
c ENDDO |
228 |
c ELSE |
229 |
C-- Make snow (according to Air Temp) and melt it in the ocean |
230 |
C note: here we just use the same criteria as over seaice but would be |
231 |
C better to consider a higher altitude air temp, e.g., 850.mb |
232 |
DO j = 1-OLy,sNy+OLy |
233 |
DO i = 1-OLx,sNx+OLx |
234 |
IF ( Tair(i,j,bi,bj).LE.Tf0kel ) THEN |
235 |
Qnet(i,j,bi,bj) = Qnet(i,j,bi,bj) |
236 |
& + Lfresh*rain(i,j,bi,bj)*rhoConstFresh |
237 |
ELSE |
238 |
C-- Account for rain-temp |
239 |
Qnet(i,j,bi,bj) = Qnet(i,j,bi,bj) |
240 |
& - HeatCapacity_Cp |
241 |
& *( Tair(i,j,bi,bj) - Tf0kel - temp_EvPrRn ) |
242 |
& *rain(i,j,bi,bj)*rhoConstFresh |
243 |
ENDIF |
244 |
ENDDO |
245 |
ENDDO |
246 |
c ENDIF |
247 |
C-- Account for energy content of Evap and RunOff: |
248 |
DO j = 1-OLy,sNy+OLy |
249 |
DO i = 1-OLx,sNx+OLx |
250 |
c Qnet(i,j,bi,bj) = Qnet(i,j,bi,bj) |
251 |
c & + ( theta(i,j,ks,bi,bj) - temp_EvPrRn ) |
252 |
c & *( evap(i,j,bi,bj)*cpwv |
253 |
c & - runoff(i,j,bi,bj)*HeatCapacity_Cp |
254 |
c & )*rhoConstFresh |
255 |
Qnet(i,j,bi,bj) = Qnet(i,j,bi,bj) |
256 |
& + ( theta(i,j,ks,bi,bj) - temp_EvPrRn ) |
257 |
& *( evap(i,j,bi,bj) - runoff(i,j,bi,bj) ) |
258 |
& *HeatCapacity_Cp*rhoConstFresh |
259 |
Qnet(i,j,bi,bj) = Qnet(i,j,bi,bj)*maskC(i,j,ks,bi,bj) |
260 |
ENDDO |
261 |
ENDDO |
262 |
ENDIF |
263 |
|
264 |
IF ( blk_taveFreq.GT.0. _d 0 ) |
265 |
& CALL BULKF_AVE( bi, bj, myThid ) |
266 |
|
267 |
C-- end bi,bj loops |
268 |
ENDDO |
269 |
ENDDO |
270 |
|
271 |
C-- Update the tile edges. |
272 |
C jmc: Is it necessary ? |
273 |
c _EXCH_XY_RS(Qnet, myThid) |
274 |
c _EXCH_XY_RS(EmPmR, myThid) |
275 |
c CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid) |
276 |
|
277 |
#endif /*ALLOW_BULK_FORCE*/ |
278 |
|
279 |
RETURN |
280 |
END |