28 std::unique_ptr<MultiFab>& mf_fcor =
vec_fcor[lev];
29 auto geomdata = Geom(lev).data();
32#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
34 for (MFIter mfi(*
cons_new[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
36 auto fcor_arr = (mf_fcor)->array(mfi);
39 Real Esize = geomdata.ProbHi()[1] - geomdata.ProbLo()[1];
40 Real prob_lo = geomdata.ProbLo()[1];
41 Real dx = geomdata.CellSize()[1];
43 ParallelFor(Box(fcor_arr), [=] AMREX_GPU_DEVICE (
int i,
int j,
int )
45 Real
y = prob_lo + (j + 0.5_rt) * dx;
46 fcor_arr(i,j,0) = coriolis_f0 + coriolis_beta * (
y - 0.5_rt * Esize);
50 vec_fcor[lev]->FillBoundary(geom[lev].periodicity());
59 std::unique_ptr<MultiFab>& mf_zeta =
vec_zeta[lev];
60 std::unique_ptr<MultiFab>& mf_Zt_avg1 =
vec_Zt_avg1[lev];
61 for ( MFIter mfi(*
cons_new[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )
64 Array4<Real>
const& Zt_avg1 = (mf_Zt_avg1)->array(mfi);
65 Array4<const Real>
const& zeta = mf_zeta->const_array(mfi);
67 Box bx3 = mfi.tilebox() ; bx3.grow(IntVect(
NGROW+1,
NGROW+1,0));
69 ParallelFor(makeSlab(bx3,2,0), [=] AMREX_GPU_DEVICE (
int i,
int j,
int )
71 Zt_avg1(i,j,0) = zeta(i,j,0,nstp);
83 auto N = Geom(lev).Domain().size()[2]-1;
90 std::unique_ptr<MultiFab>& mf_ubar =
vec_ubar[lev];
91 std::unique_ptr<MultiFab>& mf_vbar =
vec_vbar[lev];
92 std::unique_ptr<MultiFab>& mf_Hz =
vec_Hz[lev];
95 for ( MFIter mfi(*
cons_new[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )
97 Array4<Real>
const& ubar = (mf_ubar)->array(mfi);
98 Array4<Real>
const& vbar = (mf_vbar)->array(mfi);
100 Array4<const Real>
const& Hz = mf_Hz->const_array(mfi);
101 Array4<const Real>
const& u = U_old->const_array(mfi);
102 Array4<const Real>
const& v = V_old->const_array(mfi);
104 Box bx2 = mfi.tilebox() ; bx2.grow(IntVect(
NGROW ,
NGROW ,0));
105 Box ubx2 = mfi.nodaltilebox(0); ubx2.grow(IntVect(
NGROW ,
NGROW ,0));
106 Box vbx2 = mfi.nodaltilebox(1); vbx2.grow(IntVect(
NGROW ,
NGROW ,0));
108 ParallelFor(makeSlab(ubx2,2,0), [=] AMREX_GPU_DEVICE (
int i,
int j,
int )
113 for (
int k=0; k<=N; k++) {
114 Real avg_hz = 0.5_rt*(Hz(i,j,k)+Hz(i-1,j,k));
116 CF += avg_hz*u(i,j,k,nstp);
118 ubar(i,j,0,0) = CF / sum_of_hz;
121 ParallelFor(makeSlab(vbx2,2,0), [=] AMREX_GPU_DEVICE (
int i,
int j,
int )
126 for(
int k=0; k<=N; k++) {
127 Real avg_hz = 0.5_rt*(Hz(i,j,k)+Hz(i,j-1,k));
129 CF += avg_hz*v(i,j,k,nstp);
131 vbar(i,j,0,0) = CF / sum_of_hz;
154 auto N = Geom(lev).Domain().size()[2]-1;
157#pragma omp parallel if (Gpu::notInLaunchRegion())
159 for (MFIter mfi(*
vec_Akk[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi) {
160 Array4<Real>
const& Akk =
vec_Akk[lev]->array(mfi);
161 Array4<Real>
const& Akp =
vec_Akp[lev]->array(mfi);
162 Array4<Real>
const& Akt =
vec_Akt[lev]->array(mfi);
163 Array4<Real>
const& Akv =
vec_Akv[lev]->array(mfi);
165 ParallelFor(makeSlab(Box(Akk),2,0), [=] AMREX_GPU_DEVICE (
int i,
int j,
int )
167 Akk(i,j, 0) = 0.0_rt;
168 Akk(i,j, N+1) = 0.0_rt;
170 Akp(i,j, 0) = 0.0_rt;
171 Akp(i,j, N+1) = 0.0_rt;
173 Akv(i,j, 0) = 0.0_rt;
174 Akv(i,j, N+1) = 0.0_rt;
176 Akt(i,j, 0) = 0.0_rt;
177 Akt(i,j, N+1) = 0.0_rt;
196#ifdef REMORA_USE_NETCDF
198 amrex::Print() <<
"Calling init_clim_nudg_coeff_from_netcdf \n " << std::endl;
200 amrex::Print() <<
"Climatology weights loaded from netcdf file \n " << std::endl;
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
void init_gls_vmix(int lev, SolverChoice solver_choice)
Initialize GLS variables.
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive scalar
std::unique_ptr< ProblemBase > prob
Pointer to container of analytical functions for problem definition.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_tke
Turbulent kinetic energy.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_gls
Turbulent generic length scale.
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
void set_grid_scale(int lev)
Set pm and pn arrays on level lev. Only works if using analytic initialization.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Lscale
Vertical mixing turbulent length scale.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hz
Width of cells in the vertical (z-) direction (3D, Hz in ROMS)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Vertical diffusion coefficient (3D)
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
void init_beta_plane_coriolis(int lev)
Calculate Coriolis parameters from beta plane parametrization.
void FillPatch(int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bccomp, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const bool fill_set=true, const int n_not_fill=0, const int icomp_calc=0, const amrex::Real dt=amrex::Real(0.0), const amrex::MultiFab &mf_calc=amrex::MultiFab())
Fill a new MultiFab by copying in phi from valid region and filling ghost cells.
void set_zeta_average(int lev)
Set Zt_avg1 to zeta.
amrex::Vector< amrex::Real > t_new
new time at each level
static SolverChoice solverChoice
Container for algorithmic choices.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akk
Turbulent kinetic energy vertical diffusion coefficient.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
free surface height (2D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
barotropic y velocity (2D)
void set_2darrays(int lev)
Set 2D momentum arrays from 3D momentum.
void init_analytic(int lev)
Initialize initial problem data from analytic functions.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
barotropic x velocity (2D)
void init_clim_nudg_coeff_from_netcdf(int lev)
Climatology nudging coefficient initialization from NetCDF file.
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > vec_nudg_coeff
Climatology nudging coefficients.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akv
Vertical viscosity coefficient (3D)
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Average of the free surface, zeta (2D)
void init_clim_nudg_coeff(int lev)
Wrapper to initialize climatology nudging coefficient.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akp
Turbulent length scale vertical diffusion coefficient.
amrex::Real coriolis_beta
amrex::Vector< amrex::Real > nudg_coeff