21 MultiFab* mf_evap, MultiFab* mf_sustr, MultiFab* mf_svstr,
22 MultiFab* mf_stflux, MultiFab* mf_lrflx, MultiFab* mf_lhflx,
26 const int IterMax = 3;
27 const BoxArray& ba = mf_cons->boxArray();
28 const DistributionMapping& dm = mf_cons->DistributionMap();
29 MultiFab mf_Taux(ba, dm, 1, IntVect(
NGROW,
NGROW,0));
30 MultiFab mf_Tauy(ba, dm, 1, IntVect(
NGROW,
NGROW,0));
32 for ( MFIter mfi(*mf_cons, TilingIfNotGPU()); mfi.isValid(); ++mfi) {
33 Array4<Real const>
const& uwind = mf_uwind->const_array(mfi);
34 Array4<Real const>
const& vwind = mf_vwind->const_array(mfi);
35 Array4<Real const>
const& cons = mf_cons->const_array(mfi);
36 Array4<Real>
const& sustr = mf_sustr->array(mfi);
37 Array4<Real>
const& svstr = mf_svstr->array(mfi);
38 Array4<Real>
const& stflux = mf_stflux->array(mfi);
39 Array4<Real>
const& lrflx = mf_lrflx->array(mfi);
40 Array4<Real>
const& lhflx = mf_lhflx->array(mfi);
41 Array4<Real>
const& shflx = mf_shflx->array(mfi);
42 Array4<Real>
const& evap = mf_evap->array(mfi);
43 Array4<Real>
const& Taux = mf_Taux.array(mfi);
44 Array4<Real>
const& Tauy = mf_Tauy.array(mfi);
46 Array4<const Real>
const& mskr =
vec_mskr[lev]->const_array(mfi);
47 Array4<const Real>
const& msku =
vec_msku[lev]->const_array(mfi);
48 Array4<const Real>
const& mskv =
vec_mskv[lev]->const_array(mfi);
49 Array4<const Real>
const& rain =
vec_rain[lev]->const_array(mfi);
65 Box bx = mfi.tilebox();
66 Box ubx = mfi.grownnodaltilebox(0,IntVect(
NGROW-1,
NGROW-1,0));
67 Box vbx = mfi.grownnodaltilebox(1,IntVect(
NGROW-1,
NGROW-1,0));
68 Box gbx1 = bx; gbx1.grow(IntVect(
NGROW,
NGROW,0));
70 ParallelFor(makeSlab(gbx1,2,0), [=] AMREX_GPU_DEVICE (
int i,
int j,
int ) {
72 Real wind_mag = std::sqrt(uwind(i,j,0)*uwind(i,j,0) + vwind(i,j,0) * vwind(i,j,0)) + eps;
73 Real TseaK = cons(i,j,N,
Temp_comp) + 273.16_rt;
79 Real LHeat = lhflx(i,j,0) * Hscale;
80 Real SHeat = shflx(i,j,0) * Hscale;
98 Real cff=(0.7859_rt+0.03477_rt*TairC)/(1.0_rt+0.00412_rt*TairC);
99 Real e_sat=std::pow(10.0_rt,cff);
101 Real cff2=TairK*TairK*TairK;
102 Real cff1=cff2*TairK;
104 (cff1*(0.39_rt-0.05_rt*std::sqrt(vap_p))*
105 (1.0_rt-0.6823_rt*cloud*cloud)+
106 cff2*4.0_rt*(TseaK-TairK));
142 Real cff_saturation_air=(1.0007_rt+3.46E-6_rt*PairM)*6.1121_rt*
143 std::exp(17.502_rt*TairC/(240.97_rt+TairC));
147 Real Qair = 0.62197_rt*(cff_saturation_air/(PairM-0.378_rt*cff_saturation_air+eps));
152 Real cff_Q = cff_saturation_air*RH;
153 Q=0.62197_rt*(cff_Q/(PairM-0.378_rt*cff_Q+eps));
160 Real cff_saturation_water=(1.0007_rt+3.46E-6_rt*PairM)*6.1121_rt*
165 Real cff_vp=cff_saturation_water*0.98_rt;
169 Real Qsea=0.62197_rt*(cff_vp/(PairM-0.378_rt*cff_vp));
178 Real rhoAir=PairM*100.0_rt/(
blk_Rgas*TairK*(1.0_rt+0.61_rt*Q));
182 Real VisAir=1.326E-5_rt*(1.0_rt+TairC*(6.542E-3_rt+TairC*
183 (8.301E-6_rt-4.84E-9_rt*TairC)));
187 Real Hlv = (2.501_rt-0.00237_rt*cons(i,j,N,
Temp_comp))*1.0e6_rt;
193 Real delW=std::sqrt(wind_mag*wind_mag+Wgus*Wgus);
199 Real u10=delW*std::log(10.0_rt/ZoW)/std::log(blk_ZW/ZoW);
200 Real Wstar=0.035_rt * u10;
201 Real Zo10=0.011_rt*Wstar*Wstar/
g+0.11_rt*VisAir/Wstar;
202 Real Cd10 =(
vonKar/std::log(10.0_rt/Zo10));
204 Real Ch10 =0.00115_rt;
205 Real Ct10 = Ch10/std::sqrt(Cd10);
206 Real ZoT10=10.0_rt/std::exp(
vonKar/Ct10);
207 Real Cd=(
vonKar/std::log(blk_ZW/Zo10));
211 Real Ct=
vonKar/std::log(blk_ZT/ZoT10);
215 Real Ri = -
g*blk_ZW*((delT-delTc)+0.61_rt*TairK*delQ)/
216 (TairK*delW*delW+eps);
219 Zetu=CC*Ri/(1.0_rt+Ri/Ribcu);
221 Zetu=CC*Ri/(1.0_rt+3.0_rt*Ri/CC);
223 Real L10 = blk_ZW/Zetu;
226 Wstar=delW*
vonKar/(std::log(blk_ZW/Zo10)-
228 Real Tstar=-(delT-delTc)*
vonKar/(std::log(blk_ZT/ZoT10)-
230 Real Qstar=-(delQ-delQc)*
vonKar/(std::log(blk_ZQ/ZoT10)-
237 if (delW > 18.0_rt) {
239 }
else if ((10.0_rt < delW) and (delW <= 18.0_rt)) {
240 charn=0.011_rt+0.125_rt*(0.018_rt-0.011_rt)*(delW-10.0_rt);
246 for (
int it=0; it<IterMax; it++) {
247 ZoW=charn*Wstar*Wstar/
g+0.11_rt*VisAir/(Wstar+eps);
248 Real Rr=ZoW*Wstar/VisAir;
250 Real ZoQ=std::min(1.15e-4_rt,5.5e-5_rt/std::pow(Rr,0.6_rt));
252 Real ZoL=
vonKar*
g*blk_ZW*(Tstar*(1.0_rt+0.61_rt*Q)+
253 0.61_rt*TairK*Qstar)/
254 (TairK*Wstar*Wstar*(1.0_rt+0.61_rt*Q)+eps);
255 Real L=blk_ZW/(ZoL+eps);
263 Wstar=std::max(eps,delW*
vonKar/(std::log(blk_ZW/ZoW)-Wpsi));
264 Tstar=-(delT-delTc)*
vonKar/(std::log(blk_ZT/ZoT)-Tpsi);
265 Qstar=-(delQ-delQc)*
vonKar/(std::log(blk_ZQ/ZoQ)-Qpsi);
268 Real Bf=-
g/TairK*Wstar*(Tstar+0.61_rt*TairK*Qstar);
274 delW=std::sqrt(wind_mag*wind_mag+Wgus*Wgus);
278 Real Wspeed=std::sqrt(wind_mag*wind_mag+Wgus*Wgus);
279 Cd=Wstar*Wstar/(Wspeed*Wspeed+eps);
282 Real Hs=-
blk_Cpa*rhoAir*Wstar*Tstar;
285 Real diffw=2.11E-5_rt*std::pow(TairK/273.16_rt,1.94_rt);
286 Real diffh=0.02411_rt*(1.0_rt+TairC*
287 (3.309E-3_rt-1.44E-6_rt*TairC))/
289 cff=Qair*Hlv/(
blk_Rgas*TairK*TairK);
290 Real wet_bulb=1.0_rt/(1.0_rt+0.622_rt*(cff*Hlv*diffw)/
292 Real Hsr=rain(i,j,0)*wet_bulb*
blk_Cpw*
294 SHeat=(Hs+Hsr) * mskr(i,j,0);
298 Real Hl=-Hlv*rhoAir*Wstar*Qstar;
301 Real upvel=-1.61_rt*Wstar*Qstar-
302 (1.0_rt+1.61_rt*Q)*Wstar*Tstar/TairK;
303 Real Hlw=rhoAir*Hlv*upvel*Q;
304 LHeat=(Hl+Hlw) * mskr(i,j,0);
307 Taur=0.85_rt*rain(i,j,0)*wind_mag;
310 cff=rhoAir*Cd*Wspeed;
311 Real sign_u = (uwind(i,j,0) >= 0.0_rt) ? 1 : -1;
312 Real sign_v = (vwind(i,j,0) >= 0.0_rt) ? 1 : -1;
313 Taux(i,j,0)=(cff*uwind(i,j,0)+Taur*sign_u) * mskr(i,j,0);
314 Tauy(i,j,0)=(cff*vwind(i,j,0)+Taur*sign_v) * mskr(i,j,0);
345 lrflx(i,j,0) = LRad*Hscale2;
346 lhflx(i,j,0) = -LHeat*Hscale2;
347 shflx(i,j,0) = -SHeat*Hscale2;
348 stflux(i,j,0,
Temp_comp)=(srflux + lrflx(i,j,0) + lhflx(i,j,0) + shflx(i,j,0)) * mskr(i,j,0);
349 evap(i,j,0) = (LHeat / Hlv+eps) * mskr(i,j,0);
350 stflux(i,j,0,
Salt_comp) = mskr(i,j,0) * (evap(i,j,0)-rain(i,j,0)) /
rhow;
354 ParallelFor(makeSlab(ubx,2,0), [=] AMREX_GPU_DEVICE (
int i,
int j,
int ) {
355 sustr(i,j,0) = cff_rho*(Taux(i-1,j,0) + Taux(i,j,0)) * msku(i,j,0);
357 ParallelFor(makeSlab(vbx,2,0), [=] AMREX_GPU_DEVICE (
int i,
int j,
int ) {
358 svstr(i,j,0) = cff_rho*(Tauy(i,j-1,0) + Tauy(i,j,0)) * mskv(i,j,0);