39 MultiFab& mf_uold, MultiFab& mf_vold,
40 MultiFab& mf_u, MultiFab& mf_v,
43 MultiFab& S_old, MultiFab& S_new,
44 MultiFab& mf_W, MultiFab& mf_DC,
45 const MultiFab* mf_z_r,
46 const MultiFab* mf_z_w,
48 const MultiFab* mf_pm,
49 const MultiFab* mf_pn,
50 const MultiFab* mf_sustr,
51 const MultiFab* mf_svstr,
52 const MultiFab* mf_bustr,
53 const MultiFab* mf_bvstr,
54 const MultiFab* mf_msku,
55 const MultiFab* mf_mskv,
56 const int iic,
const int ntfirst,
57 const int nnew,
int nstp,
int nrhs,
58 int N,
const Real dt_lev)
60 const BoxArray& ba = S_old.boxArray();
61 const DistributionMapping& dm = S_old.DistributionMap();
65 MultiFab mf_saltcache(ba,dm,1,IntVect(
NGROW,
NGROW,0));
66 MultiFab mf_tempcache(ba,dm,1,IntVect(
NGROW,
NGROW,0));
74 for ( MFIter mfi(S_new, TilingIfNotGPU()); mfi.isValid(); ++mfi )
76 Array4<Real>
const& DC = mf_DC.array(mfi);
77 Array4<Real>
const& Akv =
vec_Akv[lev]->array(mfi);
78 Array4<Real>
const& Akt =
vec_Akt[lev]->array(mfi);
79 Array4<Real>
const& Hz =
vec_Hz[lev]->array(mfi);
80 Array4<Real>
const& Huon =
vec_Huon[lev]->array(mfi);
81 Array4<Real>
const& Hvom =
vec_Hvom[lev]->array(mfi);
82 Array4<Real const>
const& uold = mf_uold.const_array(mfi);
83 Array4<Real const>
const& vold = mf_vold.const_array(mfi);
84 Array4<Real>
const& u = (mf_u).array(mfi);
85 Array4<Real>
const& v = (mf_v).array(mfi);
87 Array4<Real const>
const& z_r = mf_z_r->const_array(mfi);
88 Array4<Real const>
const& z_w = mf_z_w->const_array(mfi);
89 Array4<Real const>
const& h = mf_h->const_array(mfi);
90 Array4<Real const>
const& pm = mf_pm->const_array(mfi);
91 Array4<Real const>
const& pn = mf_pn->const_array(mfi);
93 Array4<Real >
const& ru = mf_ru->array(mfi);
94 Array4<Real >
const& rv = mf_rv->array(mfi);
95 Array4<Real >
const& W = (mf_W).array(mfi);
96 Array4<Real const>
const& sustr = mf_sustr->const_array(mfi);
97 Array4<Real const>
const& svstr = mf_svstr->const_array(mfi);
98 Array4<Real const>
const& bustr = mf_bustr->const_array(mfi);
99 Array4<Real const>
const& bvstr = mf_bvstr->const_array(mfi);
100 Array4<Real const>
const& msku = mf_msku->const_array(mfi);
101 Array4<Real const>
const& mskv = mf_mskv->const_array(mfi);
103 Real lambda = 1.0_rt;
105 Box bx = mfi.tilebox();
106 Box gbx = mfi.growntilebox();
107 Box gbx1 = mfi.growntilebox(IntVect(
NGROW-1,
NGROW-1,0));
108 Box gbx2 = mfi.growntilebox(IntVect(
NGROW,
NGROW,0));
134 tbxp1D.makeSlab(2,0);
136 tbxp2D.makeSlab(2,0);
138 FArrayBox fab_FC(convert(tbxp2,IntVect(0,0,1)),1,amrex::The_Async_Arena());
140 auto FC=fab_FC.array();
143 for (
int i_comp=0; i_comp <
NCONS; i_comp++) {
144 Array4<Real>
const& sstore = (
vec_sstore[lev])->array(mfi,i_comp);
145#ifdef REMORA_USE_NETCDF
146 FArrayBox* fab_river_source;
152 const Array4<const Real>& river_source = (
solverChoice.
do_rivers_cons[i_comp]) ? fab_river_source->const_array() : Array4<const Real>();
154 const Array4<const int >& river_pos = Array4<const int>();
155 const Array4<const Real>& river_source = Array4<const Real>();
158 mf_scalarcache.array(mfi,i_comp), Hz, Huon, Hvom,
159 W, DC, FC, sstore, z_w, h, pm, pn, msku, mskv, river_pos,
160 river_source, iic, ntfirst,
165 for (
int i_comp=0; i_comp <
NCONS; i_comp++) {
166 const Array4<Real const>& stflx =
vec_stflx[lev]->const_array(mfi,i_comp);
167 const Array4<Real const>& btflx =
vec_btflx[lev]->const_array(mfi,i_comp);
168 prestep_diffusion(bx,gbx,0,0,S_new.array(mfi,i_comp), S_old.array(mfi,i_comp), ru,
169 Hz, Akt, DC, FC, stflx, btflx, z_r, pm, pn, iic, iic, nnew, nstp,
170 nrhs, N, lambda, dt_lev);
180 prestep_diffusion(tbxp1, gbx, 1, 0, u, uold, ru, Hz, Akv, DC, FC,
181 sustr, bustr, z_r, pm, pn, iic, ntfirst, nnew, nstp, nrhs, N, lambda, dt_lev);
183 prestep_diffusion(tbxp1, gbx, 0, 1, v, vold, rv, Hz, Akv, DC, FC,
184 svstr, bvstr, z_r, pm, pn, iic, ntfirst, nnew, nstp, nrhs, N, lambda, dt_lev);
void prestep_diffusion(const amrex::Box &bx, const amrex::Box &gbx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &vel, const amrex::Array4< amrex::Real const > &vel_old, const amrex::Array4< amrex::Real > &rvel, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &Akv, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &sstr, const amrex::Array4< amrex::Real const > &bstr, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const int iic, const int ntfirst, const int nnew, int nstp, int nrhs, int N, const amrex::Real lambda, const amrex::Real dt_lev)
Update velocities or tracers with diffusion/viscosity as the last part of the prestep.
void prestep(int lev, amrex::MultiFab &mf_uold, amrex::MultiFab &mf_vold, amrex::MultiFab &mf_u, amrex::MultiFab &mf_v, amrex::MultiFab *mf_ru, amrex::MultiFab *mf_rv, amrex::MultiFab &S_old, amrex::MultiFab &S_new, amrex::MultiFab &mf_W, amrex::MultiFab &mf_DC, const amrex::MultiFab *mf_z_r, const amrex::MultiFab *mf_z_w, const amrex::MultiFab *mf_h, const amrex::MultiFab *mf_pm, const amrex::MultiFab *mf_pn, const amrex::MultiFab *mf_sustr, const amrex::MultiFab *mf_svstr, const amrex::MultiFab *mf_bustr, const amrex::MultiFab *mf_bvstr, const amrex::MultiFab *mf_msku, const amrex::MultiFab *mf_mskv, const int iic, const int nfirst, const int nnew, int nstp, int nrhs, int N, const amrex::Real dt_lev)
Wrapper function for prestep.
void prestep_t_advection(const amrex::Box &tbx, const amrex::Box &gbx, const amrex::Array4< amrex::Real > &tempold, const amrex::Array4< amrex::Real > &tempcache, const amrex::Array4< amrex::Real > &Hz, const amrex::Array4< amrex::Real > &Huon, const amrex::Array4< amrex::Real > &Hvom, const amrex::Array4< amrex::Real > &W, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real > &sstore, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &msku, const amrex::Array4< amrex::Real const > &mskv, const amrex::Array4< int const > &river_pos, const amrex::Array4< amrex::Real const > &river_source, int iic, int ntfirst, int nrhs, int N, const amrex::Real dt_lev)
Prestep advection calculations for the tracers.