1 |
c swd -- bulkf formula pkg |
2 |
|
3 |
#include "CPP_OPTIONS.h" |
4 |
|
5 |
subroutine bulkf_forcing( bi,bj, |
6 |
I mycurrenttime, |
7 |
I mycurrentiter, |
8 |
I mythid |
9 |
I ) |
10 |
|
11 |
c ================================================================== |
12 |
c SUBROUTINE BULKF_FORCING |
13 |
c ================================================================== |
14 |
c |
15 |
c o Get the surface fluxes used to force ocean model |
16 |
c Output: |
17 |
c ------ |
18 |
c ustress, vstress - wind stress |
19 |
c qnet - net heat flux |
20 |
c empmr - freshwater flux |
21 |
c --------- |
22 |
c |
23 |
c Input: |
24 |
c ------ |
25 |
c uwind, vwind - mean wind speed (m/s) at height hu (m) |
26 |
c Tair - mean air temperature (K) at height ht (m) |
27 |
c Qair - mean air humidity (kg/kg) at height hq (m) |
28 |
c theta(k=1) - sea surface temperature (K) |
29 |
c rain - precipitation |
30 |
c runoff - river(ice) runoff |
31 |
c |
32 |
c ================================================================== |
33 |
c SUBROUTINE bulkf_forcing |
34 |
c ================================================================== |
35 |
|
36 |
implicit none |
37 |
|
38 |
c == global variables == |
39 |
|
40 |
#include "EEPARAMS.h" |
41 |
#include "SIZE.h" |
42 |
#include "PARAMS.h" |
43 |
#include "DYNVARS.h" |
44 |
#include "GRID.h" |
45 |
#include "FFIELDS.h" |
46 |
#ifdef ALLOW_BULK_FORCE |
47 |
#include "BULKF.h" |
48 |
#include "BULKF_DIAG.h" |
49 |
#endif |
50 |
#ifdef ALLOW_THERM_SEAICE |
51 |
#include "ICE.h" |
52 |
#endif |
53 |
c == routine arguments == |
54 |
|
55 |
integer mythid |
56 |
integer mycurrentiter |
57 |
_RL mycurrenttime |
58 |
integer bi,bj |
59 |
|
60 |
C == Local variables == |
61 |
integer i,j,k |
62 |
|
63 |
#ifdef ALLOW_BULK_FORCE |
64 |
|
65 |
_RL evap(1-olx:snx+olx,1-oly:sny+oly) |
66 |
_RL ustress(1-olx:snx+olx,1-oly:sny+oly) |
67 |
_RL vstress(1-olx:snx+olx,1-oly:sny+oly) |
68 |
_RL savssq(1-olx:snx+olx,1-oly:sny+oly) |
69 |
_RL fsh(1-olx:snx+olx,1-oly:sny+oly) |
70 |
_RL flh(1-olx:snx+olx,1-oly:sny+oly) |
71 |
_RL flwup(1-olx:snx+olx,1-oly:sny+oly) |
72 |
_RL fswnet(1-olx:snx+olx,1-oly:sny+oly) |
73 |
|
74 |
_RL df0dT,hfl |
75 |
|
76 |
c variables to include seaice effect |
77 |
_RL tmp |
78 |
_RL albedo |
79 |
integer iceornot |
80 |
|
81 |
c == external functions == |
82 |
|
83 |
#endif /* ALLOW_BULK_FORCE */ |
84 |
|
85 |
c determine forcing field records |
86 |
|
87 |
#ifdef ALLOW_BULK_FORCE |
88 |
|
89 |
|
90 |
k = 1 |
91 |
do j = 1-oly,sny+oly |
92 |
do i = 1-olx,snx+olx |
93 |
if (hFacC(i,j,1,bi,bj).ne.0.0) then |
94 |
#ifdef ALLOW_THERM_SEAICE |
95 |
if (ICEMASK(i,j,bi,bj).gt.0.d0) then |
96 |
tmp=Tsrf(i,j,bi,bj) |
97 |
if (snowheight(i,j,bi,bj).gt.3.d-1) then |
98 |
iceornot=2 |
99 |
else |
100 |
iceornot=1 |
101 |
endif |
102 |
else |
103 |
tmp=theta(i,j,1,bi,bj) |
104 |
iceornot=0 |
105 |
endif |
106 |
#else |
107 |
tmp=theta(i,j,1,bi,bj) |
108 |
iceornot=0 |
109 |
#endif |
110 |
call BULKF_FORMULA_LANL( |
111 |
& uwind(i,j,bi,bj), vwind(i,j,bi,bj), |
112 |
& wspeed(i,j,bi,bj), Tair(i,j,bi,bj), Qair(i,j,bi,bj), |
113 |
& cloud(i,j,bi,bj), tmp, |
114 |
& flwup(i,j), flh(i,j), |
115 |
& fsh(i,j), df0dT, |
116 |
& ustress(i,j), vstress(i,j), |
117 |
& evap(i,j), savssq(i,j), |
118 |
& iceornot, readwindstress) |
119 |
cQQ use down long wave data |
120 |
flwup(i,j)=flwup(i,j)-flw(i,j,bi,bj) |
121 |
cQQ using down solar, need to have water albedo -- .1 |
122 |
#ifdef ALLOW_THERM_SEAICE |
123 |
if (ICEMASK(i,j,bi,bj).gt.0.d0) then |
124 |
call sfc_albedo(ICEHEIGHT(i,j,bi,bj), |
125 |
& SNOWHEIGHT(i,j,bi,bj), |
126 |
& Tsrf(i,j,bi,bj), |
127 |
& sage(i,j,bi,bj), albedo) |
128 |
else |
129 |
albedo=.1 |
130 |
endif |
131 |
#else |
132 |
albedo=.1 |
133 |
#endif |
134 |
fswnet(i,j)=solar(i,j,bi,bj)*(1.d0-albedo) |
135 |
else |
136 |
ustress(i,j) = 0. _d 0 |
137 |
vstress(i,j) = 0. _d 0 |
138 |
fsh(i,j) = 0. _d 0 |
139 |
flh(i,j) = 0. _d 0 |
140 |
flwup(i,j) =0. _d 0 |
141 |
evap(i,j) =0. _d 0 |
142 |
fswnet(i,j) =0. _d 0 |
143 |
savssq(i,j) =0. _d 0 |
144 |
endif |
145 |
enddo |
146 |
enddo |
147 |
|
148 |
|
149 |
if ( .NOT.readwindstress) then |
150 |
cswd move wind stresses to u and v points |
151 |
do j = 1,sny |
152 |
do i = 1,snx |
153 |
fu(i,j,bi,bj)= |
154 |
& (ustress(i,j)+ustress(i-1,j))/2* |
155 |
& maskW(i,j,1,bi,bj) |
156 |
fv(i,j,bi,bj)= |
157 |
& (vstress(i,j)+vstress(i,j-1))/2* |
158 |
& maskS(i,j,1,bi,bj) |
159 |
enddo |
160 |
enddo |
161 |
endif |
162 |
c |
163 |
c |
164 |
c Add all contributions. |
165 |
k = 1 |
166 |
do j = 1,sny |
167 |
do i = 1,snx |
168 |
if (hFacC(i,j,1,bi,bj).ne.0.0) then |
169 |
c Net surface heat flux. |
170 |
hfl = 0. _d 0 |
171 |
hfl = hfl + fsh(i,j) |
172 |
hfl = hfl + flh(i,j) |
173 |
hfl = hfl + flwup(i,j) |
174 |
c QQ minus solar for ncep data |
175 |
c QQ plus solar for dasilva |
176 |
c hfl = hfl - solar(i,j,bi,bj) |
177 |
cQQ for calculated sw net |
178 |
hfl = hfl - fswnet(i,j) |
179 |
c Heat flux: |
180 |
Qnet(i,j,bi,bj) = -hfl*maskc(i,j,k,bi,bj) |
181 |
cswdcou -- add --- |
182 |
#ifdef COUPLE_MODEL |
183 |
dFdT(i,j,bi,bj) = df0dT |
184 |
#endif |
185 |
cswdcou -- end add --- |
186 |
c Salt flux from Precipitation and Evaporation. |
187 |
EmPmR(i,j,bi,bj) = (-evap(i,j)+rain(i,j,bi,bj) |
188 |
& - runoff(i,j,bi,bj))*maskc(i,j,k,bi,bj) |
189 |
|
190 |
cccccccccccCHEATccccccccccccccccccccccccc |
191 |
c Qnet(i,j,bi,bj) = Qnetch(i,j,bi,bj) |
192 |
c EmPmR(i,j,bi,bj) = EmPch(i,j,bi,bj) |
193 |
cccccccccccccccccccccccccccccccccccccccccc |
194 |
else |
195 |
Qnet(i,j,bi,bj) =0.d0 |
196 |
EmPmR(i,j,bi,bj)=0.d0 |
197 |
cswdcou -- add --- |
198 |
#ifdef COUPLE_MODEL |
199 |
dFdT(i,j,bi,bj) = 0.d0 |
200 |
#endif |
201 |
cswdcou -- end add --- |
202 |
endif |
203 |
enddo |
204 |
enddo |
205 |
|
206 |
|
207 |
cc Update the tile edges. |
208 |
c _EXCH_XY_R8(Qnet, mythid) |
209 |
c _EXCH_XY_R8(EmPmR, mythid) |
210 |
c _EXCH_XY_R8(fu , mythid) |
211 |
c _EXCH_XY_R8(fv , mythid) |
212 |
|
213 |
|
214 |
#ifdef ALLOW_TIMEAVE |
215 |
call BULKF_AVE(bi,bj,mythid, fswnet, |
216 |
& flh, fsh, flwup, evap, savssq) |
217 |
#endif /*ALLOW_TIMEAVE*/ |
218 |
|
219 |
#endif /*ALLOW_BULK_FORCE*/ |
220 |
|
221 |
return |
222 |
end |