1#ifndef INTEGRATOR_TOPOP_H
2#define INTEGRATOR_TOPOP_H
9#include "AMReX_ParallelDescriptor.H"
10#include "AMReX_ParmParse.H"
27#include "Numeric/Stencil.H"
33#include "Operator/Operator.H"
84 AMREX_D_TERM((value.geom[0].ProbHi()[0] - value.geom[0].ProbLo()[0]),
85 *(value.geom[0].ProbHi()[1] - value.geom[0].ProbLo()[1]),
86 *(value.geom[0].ProbHi()[2] - value.geom[0].ProbLo()[2]));
111 if (a_step > 0)
return;
113 for (
int lev = 0; lev <= finest_level; ++lev)
126 BL_PROFILE(
"TopOp::Advance");
130 amrex::Box domain = geom[lev].Domain();
135 for (amrex::MFIter mfi(*
psi_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
137 amrex::Box bx = mfi.tilebox();
140 amrex::Array4<const Set::Matrix>
const& sig = (*
stress_mf[lev]).array(mfi);
141 amrex::Array4<const Set::Matrix>
const& eps = (*
strain_mf[lev]).array(mfi);
142 amrex::Array4<const Set::Scalar>
const& psi = (*
psi_old_mf[lev]).array(mfi);
143 amrex::Array4<Set::Scalar>
const& psinew = (*
psi_mf[lev]).array(mfi);
145 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
149 driving_force +=
alpha * 2.0 * psi(i, j, k) * (2.0 * psi(i, j, k) * psi(i, j, k) - 3.0 * psi(i, j, k) + 1.0);
155 driving_force += -
gamma * 0.5 * (sig_avg.transpose() * eps_avg).trace() * psi(i, j, k);
159 psinew(i, j, k) = psi(i, j, k) - Lnow *
dt * driving_force;
160 if (psinew(i, j, k) < 0.0) psinew(i, j, k) = 0.0;
161 if (psinew(i, j, k) > 1.0) psinew(i, j, k) = 1.0;
167 const amrex::MFIter& mfi,
const amrex::Box&
box)
override
169 BL_PROFILE(
"TopOp::Integrate");
172 const amrex::Real* DX = geom[amrlev].CellSize();
173 Set::Scalar dv = AMREX_D_TERM(DX[0], *DX[1], *DX[2]);
175 amrex::Array4<amrex::Real>
const& psi = (*
psi_mf[amrlev]).array(mfi);
176 amrex::Array4<const Set::Matrix>
const& sig = (*
stress_mf[amrlev]).array(mfi);
177 amrex::Array4<const Set::Matrix>
const& eps = (*
strain_mf[amrlev]).array(mfi);
178 amrex::ParallelFor(
box, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
180 volume += psi(i, j, k, 0) * dv;
181 w_chem_potential +=
alpha * psi(i, j, k, 0) * psi(i, j, k, 0) * (1. - psi(i, j, k, 0) * psi(i, j, k, 0)) * dv;
186 w_elastic +=
gamma * 0.5 * (sig_avg.transpose() * eps_avg).trace() * psi(i, j, k) * dv;
198 for (amrex::MFIter mfi(*
model_mf[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
200 amrex::Box bx = mfi.tilebox();
201 amrex::Array4<char>
const& tags = a_tags.array(mfi);
204 amrex::Array4<Set::Scalar>
const& psi =
psi_mf[lev]->array(mfi);
205 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(
int i,
int j,
int k)
211 tags(i, j, k) = amrex::TagBox::SET;
#define pp_query_default(...)
#define pp_queryclass(...)
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)
Set::Field< Set::Matrix > strain_mf
void Integrate(int amrlev, Set::Scalar, int, const amrex::MFIter &mfi, const amrex::Box &a_box) override
static void Parse(Mechanics &value, IO::ParmParse &pp)
void Initialize(int lev) override
Use the #ic object to initialize::Temp.
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int) override
Set::Field< Set::Matrix > stress_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.
std::vector< amrex::Box > box
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
void RegisterIntegratedVariable(Set::Scalar *integrated_variable, std::string name, bool extensive=true)
Register a variable to be integrated over the spatial domain using the Integrate function.
void Integrate(int amrlev, Set::Scalar time, int step, const amrex::MFIter &mfi, const amrex::Box &box) override
BC::BC< Set::Scalar > * bc_psi
Set::Scalar m_eta_ref_threshold
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Set::Scalar w_chem_potential
virtual void UpdateModel(int a_step, Set::Scalar) override
Set::Field< Set::Scalar > psi_old_mf
void Initialize(int lev) override
Numeric::Interpolator::Linear< Set::Scalar > lambda
Numeric::Interpolator::Linear< Set::Scalar > L
static void Parse(TopOp &value, IO::ParmParse &pp)
IC::IC< Set::Scalar > * ic_psi
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar a_time, int a_ngrow) override
Collection of numerical integrator objects.
AMREX_FORCE_INLINE Set::Scalar Laplacian(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)
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 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
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T > > &a_mf, const amrex::Geometry &)
void Abort(const char *msg)
static AMREX_FORCE_INLINE T NodeToCellAverage(const amrex::Array4< const T > &f, const int &i, const int &j, const int &k, const int &m)