4#ifndef INTEGRATOR_PHASEFIELDMICROSTRUCTURE_H
5#define INTEGRATOR_PHASEFIELDMICROSTRUCTURE_H
9#include "AMReX_ParmParse.H"
10#include <AMReX_MLMG.H>
22#include "Model/Interface/GB/GB.H"
41template<
class model_type>
45 const std::string
name =
"phasefieldmicrostructure." + std::string(model_type::name);
58 BL_PROFILE(
"PhaseFieldMicrostructure::Parse");
76 value.pf.threshold.on =
77 value.pf.threshold.chempot ||
value.pf.threshold.boundary ||
78 value.pf.threshold.corner ||
value.pf.threshold.lagrange ||
79 value.pf.threshold.mechanics;
81 std::string
type =
"continuous";
96 pp.
query_validate(
"mechanics.type", type_str, {
"disable",
"static",
"dynamic"});
99 if (type_str !=
"disable")
115 value.number_of_ghost_cells = std::max(
value.number_of_ghost_cells, 3);
121 if (
value.lagrange.on)
126 value.SetThermoInt(1);
132 std::vector<std::string> vals;
134 int nvals =
static_cast<int>(vals.size());
136 for (
int i = 0; i <
value.number_of_grains; i++)
138 else if (nvals ==
value.number_of_grains)
139 for (
int i = 0; i <
value.number_of_grains; i++)
142 Util::Abort(
INFO,
"sdf.val received ", vals.size(),
" but requires 1 or ",
value.number_of_grains);
150 if (
value.anisotropy.on)
171 if (
value.anisotropy.on)
172 value.number_of_ghost_cells = std::max(
value.number_of_ghost_cells,2);
175 std::map<std::string, RegularizationType> regularization_type;
179 std::string regularization_type_input;
181 pp_query_validate(
"anisotropy.regularization", regularization_type_input,{
"k12",
"wilmore"});
182 value.regularization = regularization_type[regularization_type_input];
184 pp_forbid(
"anisotropy.gb_type",
" --> anisotropy.type");
189 pp.query(
"fluctuation.on",
value.fluctuation.on);
190 if (
value.fluctuation.on)
192 pp.query(
"fluctuation.amp",
value.fluctuation.amp);
193 pp.query(
"fluctuation.sd",
value.fluctuation.sd);
194 pp.query(
"fluctuation.tstart",
value.fluctuation.tstart);
195 value.fluctuation.norm_dist = std::normal_distribution<double>(0.0,
value.fluctuation.sd);
199 pp.query(
"disconnection.on",
value.disconnection.on);
200 if (
value.disconnection.on)
208 if (
value.shearcouple.on)
213 value.shearcouple.Fgb.resize(
value.number_of_grains *
value.number_of_grains, Set::Matrix::Zero());
214 for (
int i = 0 ; i <
value.number_of_grains ; i++)
215 for (
int j = 0 ; j <
value.number_of_grains ; j++)
217 std::string
name =
"shearcouple.Fgb."+std::to_string(i)+
"."+std::to_string(j);
218 std::string namerev =
"shearcouple.Fgb."+std::to_string(j)+
"."+std::to_string(i);
226 value.shearcouple.Fgb[j*
value.number_of_grains + i] = -
value.shearcouple.Fgb[i*
value.number_of_grains + j];
239 if (
value.anisotropic_kinetics.on)
243 std::string mobility_filename, threshold_filename;
245 pp_query_file(
"anisotropic_kinetics.mobility",mobility_filename);
248 pp_query_file(
"anisotropic_kinetics.threshold",threshold_filename);
250 value.RegisterNewFab(
value.anisotropic_kinetics.L_mf,
value.mybc,
value.number_of_grains, 0,
"mobility",
true);
251 value.RegisterNewFab(
value.anisotropic_kinetics.threshold_mf,
value.mybc,
value.number_of_grains, 0,
"theshold",
true);
265 value.RegisterIntegratedVariable(&
value.volume,
"volume");
266 value.RegisterIntegratedVariable(&
value.area,
"area");
267 value.RegisterIntegratedVariable(&
value.gbenergy,
"gbenergy");
268 value.RegisterIntegratedVariable(&
value.realgbenergy,
"realgbenergy");
269 value.RegisterIntegratedVariable(&
value.regenergy,
"regenergy");
278 void Advance (
int lev, Real time, Real
dt)
override;
281 void TagCellsForRefinement (
int lev, amrex::TagBoxArray& tags, amrex::Real time,
int ngrow)
override;
283 virtual void TimeStepBegin(amrex::Real time,
int iter)
override;
286 const amrex::MFIter &mfi,
const amrex::Box &
box)
override;
291 for (
int lev = 0; lev <= this->max_level; lev++)
374 std::vector<Numeric::Interpolator::Linear<Set::Scalar>>
val;
394 std::vector<Set::Matrix>
Fgb;
#define pp_query_validate(...)
#define pp_query_required(...)
#define pp_query_file(...)
#define pp_query_default(...)
#define pp_queryclass(...)
Pure abstract IC object from which all other IC objects inherit.
Set each point to a random value.
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__)
int queryarr(std::string name, std::vector< T > &value, std::string, std::string, int)
int query_validate(std::string name, int &value, std::vector< int > possibleintvals, std::string file="", std::string func="", int line=-1)
Set::Field< Set::Matrix > stress_mf
Set::Field< model_type > model_mf
std::vector< amrex::Box > box
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Solve the Allen-Cahn evolution equation for microstructure with parameters , where n corresponds to t...
Set::Field< Set::Scalar > threshold_mf
void Initialize(int lev) override
void Advance(int lev, Real time, Real dt) override
Evolve phase field in time.
Numeric::Interpolator::Linear< Set::Scalar > mobility
virtual void UpdateModel(int, Set::Scalar) override
BC::BC< Set::Scalar > * mybc
struct Integrator::PhaseFieldMicrostructure::@16::@25 threshold
Set::Field< Set::Scalar > L_mf
RegularizationType regularization
static void Parse(PhaseFieldMicrostructure &value, IO::ParmParse &pp)
virtual void UpdateEigenstrain()
struct Integrator::PhaseFieldMicrostructure::@21 fluctuation
PhaseFieldMicrostructure()
PhaseFieldMicrostructure(IO::ParmParse &pp)
Model::Defect::Disconnection model
void TagCellsForRefinement(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Set::Field< Set::Scalar > fluct_mf
std::default_random_engine rand_num_gen
Set::Field< Set::Scalar > driving_force_mf
void Integrate(int amrlev, Set::Scalar time, int step, const amrex::MFIter &mfi, const amrex::Box &box) override
virtual void TimeStepBegin(amrex::Real time, int iter) override
IC::IC< Set::Scalar > * ic
struct Integrator::PhaseFieldMicrostructure::@23 shearcouple
Set::Field< Set::Scalar > totaldf_mf
struct Integrator::PhaseFieldMicrostructure::@22 disconnection
Set::Scalar ref_threshold
Numeric::Interpolator::Linear< Set::Scalar > threshold
std::vector< Set::Matrix > Fgb
struct Integrator::PhaseFieldMicrostructure::@18 anisotropy
std::normal_distribution< double > norm_dist
struct Integrator::PhaseFieldMicrostructure::@17 anisotropic_kinetics
Set::Field< Set::Scalar > driving_force_threshold_mf
struct Integrator::PhaseFieldMicrostructure::@16 pf
int number_of_ghost_cells
std::vector< model_type > model
std::vector< Numeric::Interpolator::Linear< Set::Scalar > > val
virtual void TimeStepComplete(amrex::Real time, int iter) override
bool model_neumann_boundary
Set::Field< Set::Scalar > eta_old_mf
virtual ~PhaseFieldMicrostructure()
Set::Field< Set::Scalar > eta_mf
Reads the data from a file and computes energies and its derivates.
A 2D interface model class.
static Linear< T > Read(std::string filename, int derivative=0)
Collection of numerical integrator objects.
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
AMREX_FORCE_INLINE void AssertException(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
void Abort(const char *msg)