59#ifndef INTEGRATOR_MECHANICS_H
60#define INTEGRATOR_MECHANICS_H
67#include "AMReX_ParallelDescriptor.H"
68#include "AMReX_ParmParse.H"
89#include "Numeric/Stencil.H"
95#include "Operator/Operator.H"
104 static constexpr const char *
name =
"mechanics";
141 if (nmodels > 1 && pp.
contains(
"ic.type"))
170 if (pp.
contains(
"trac_normal.ic.type"))
190 else eta_mf[lev]->setVal(1.0);
202 for (
int lev = 0; lev <= finest_level; ++lev)
205 psi_mf[lev]->FillBoundary();
206 amrex::Box domain = this->geom[lev].Domain();
207 domain.convert(amrex::IntVect::TheNodeVector());
210 for (MFIter mfi(*
model_mf[lev],
false); mfi.isValid(); ++mfi)
212 amrex::Box bx = mfi.nodaltilebox();
215 amrex::Array4<Set::Vector>
const& rhs =
rhs_mf[lev]->array(mfi);
216 amrex::Array4<const Set::Scalar>
const& psi =
psi_mf[lev]->array(mfi);
217 amrex::Array4<const Set::Scalar>
const& trac_normal =
trac_normal_mf[lev]->array(mfi);
219 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
221 rhs(i,j,k) = trac_normal(i,j,k) * grad;
228 if (a_step > 0)
return;
230 for (
int lev = 0; lev <= finest_level; ++lev)
232 eta_mf[lev]->FillBoundary();
234 amrex::Box domain = this->geom[lev].Domain();
235 domain.convert(amrex::IntVect::TheNodeVector());
239 for (MFIter mfi(*
model_mf[lev],
false); mfi.isValid(); ++mfi)
241 amrex::Box bx = mfi.grownnodaltilebox();
244 amrex::Array4<MODEL>
const& model =
model_mf[lev]->array(mfi);
245 amrex::Array4<const Set::Scalar>
const& eta =
eta_mf[lev]->array(mfi);
247 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k) {
248 model(i, j, k) = MODEL::Zero();
249 for (
unsigned int n = 0; n <
models.size(); n++)
250 model(i, j, k) += eta(i, j, k, n) *
models[n];
258 for (MFIter mfi(*
model_mf[lev],
false); mfi.isValid(); ++mfi)
260 amrex::Box bx = mfi.grownnodaltilebox() & domain;
261 amrex::Array4<MODEL>
const& model = this->
model_mf[lev]->array(mfi);
262 const Dim3 lo = amrex::lbound(domain), hi = amrex::ubound(domain);
264 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
266 if (i==lo.x && j==lo.y)
267 model(i,j,k) = 0.5*(model(i+1,j,k)+model(i,j+1,k));
268 else if (i==lo.x && j==hi.y)
269 model(i,j,k) = 0.5*(model(i+1,j,k)+model(i,j-1,k));
270 else if (i==hi.x && j==lo.y)
271 model(i,j,k) = 0.5*(model(i-1,j,k)+model(i,j+1,k));
272 else if (i==hi.x && j==hi.y)
273 model(i,j,k) = 0.5*(model(i-1,j,k)+model(i,j-1,k));
276 model(i,j,k) = model(i+1,j,k);
278 model(i,j,k) = model(i-1,j,k);
280 model(i,j,k) = model(i,j+1,k);
282 model(i,j,k) = model(i,j-1,k);
298 for (amrex::MFIter mfi(*
eta_mf[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
300 amrex::Box bx = mfi.nodaltilebox();
301 amrex::Array4<char>
const& tags = a_tags.array(mfi);
302 amrex::Array4<Set::Scalar>
const& eta =
eta_mf[lev]->array(mfi);
303 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
306 for (
int n = 0; n <
eta_mf[lev]->nComp(); n++)
310 tags(i, j, k) = amrex::TagBox::SET;
315 amrex::Array4<Set::Scalar>
const& psi =
psi_mf[lev]->array(mfi);
316 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
322 tags(i, j, k) = amrex::TagBox::SET;
#define pp_query_default(...)
Pure abstract IC object from which all other IC objects inherit.
void Initialize(const int &a_lev, Set::Field< T > &a_field, Set::Scalar a_time=0.0)
bool contains(std::string name)
void queryclass(std::string name, T *value, std::string file="", std::string func="", int line=-1)
void select(std::string name, PTRTYPE *&ic_eta, Args &&... args)
int queryclass_enumerate(std::string a_name, std::vector< T > &value, int number=1, std::string file="", std::string func="", int line=__LINE__)
Set::Scalar m_elastic_ref_threshold
void Initialize(int lev) override
Use the #ic object to initialize::Temp.
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int) override
Set::Field< Set::Vector > rhs_mf
Set::Field< Set::Scalar > psi_mf
Set::Field< MODEL > model_mf
void RegisterNewFab(Set::Field< Set::Scalar > &new_fab, BC::BC< Set::Scalar > *new_bc, int ncomp, int nghost, std::string name, bool writeout, std::vector< std::string > suffix={})
Add a new cell-based scalar field.
void RegisterNodalFab(Set::Field< Set::Scalar > &new_fab, int ncomp, int nghost, std::string name, bool writeout, std::vector< std::string > suffix={})
Add a new node-based scalar field.
IC::IC< Set::Scalar > * ic_trac_normal
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar a_time, int a_ngrow) override
IC::IC< Set::Scalar > * ic_eta
Set::Field< Set::Scalar > eta_mf
static void Parse(Mechanics &value, IO::ParmParse &pp)
void Regrid(int lev, Set::Scalar time) override
Mechanics(IO::ParmParse &pp)
std::vector< MODEL > models
BC::BC< Set::Scalar > * bc_trac_normal
void Initialize(int lev) override
bool model_neumann_boundary
BC::BC< Set::Scalar > * bc_psi
Set::Scalar m_eta_ref_threshold
static constexpr const char * name
IC::IC< Set::Scalar > * ic_psi
Set::Field< Set::Scalar > trac_normal_mf
virtual void UpdateModel(int a_step, Set::Scalar time) override
Collection of numerical integrator objects.
static AMREX_FORCE_INLINE std::array< StencilType, AMREX_SPACEDIM > GetStencil(const int i, const int j, const int k, const amrex::Box domain)
AMREX_FORCE_INLINE Set::Vector CellGradientOnNode(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM])
AMREX_FORCE_INLINE Set::Vector Gradient(const amrex::Array4< const Set::Scalar > &f, const int &i, const int &j, const int &k, const int &m, const Set::Scalar dx[AMREX_SPACEDIM], std::array< StencilType, AMREX_SPACEDIM > stencil=DefaultType)
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T > > &a_mf, const amrex::Geometry &)
AMREX_FORCE_INLINE void AssertException(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
AMREX_FORCE_INLINE void Assert(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
void Message(std::string file, std::string func, int line, Args const &... args)