Alamo
PhaseFieldMicrostructure.H
Go to the documentation of this file.
1 ///
2 /// \file PhaseFieldMicrostructure.H
3 ///
4 #ifndef INTEGRATOR_PHASEFIELDMICROSTRUCTURE_H
5 #define INTEGRATOR_PHASEFIELDMICROSTRUCTURE_H
6 
7 #include <iostream>
8 #include <fstream>
9 #include <iomanip>
10 
11 #include "AMReX.H"
12 #include "AMReX_ParmParse.H"
13 #include "AMReX_ParallelDescriptor.H"
14 #include <AMReX_MLMG.H>
15 
16 #include "Integrator/Integrator.H"
17 #include "Integrator/Mechanics.H"
18 
19 #include "BC/BC.H"
20 #include "BC/Constant.H"
22 #include "IC/Constant.H"
23 #include "IC/TabulatedInterface.H"
24 #include "IC/PerturbedInterface.H"
25 #include "IC/Voronoi.H"
26 #include "IC/Sphere.H"
27 #include "IC/Expression.H"
28 #include "IC/Ellipse.H"
29 #include "IC/Random.H"
30 
31 #include "Model/Interface/GB/GB.H"
32 #include "Model/Interface/GB/Sin.H"
35 #include "Model/Interface/GB/SH.H"
36 
39 #include "Operator/Elastic.H"
40 
41 namespace Integrator
42 {
43 //using model_type = Model::Solid::Affine::Cubic;
44 
46  Wilmore = 1,
47  K12 = 2};
48 
49 ///
50 /// \class PhaseFieldMicrostructure
51 /// \brief Microstructure evolution with grain boundary anisotropy
52 ///
53 /// Solve the Allen-Cahn evolution equation for microstructure with parameters \f$\eta_1\ldots\eta_n\f$,
54 /// where n corresponds to the number of grains.
55 ///
56 template<class model_type>
57 class PhaseFieldMicrostructure : public Base::Mechanics<model_type>
58 {
59 public:
62  {pp_queryclass(*this);}
64  {
65  delete boundary;
66  delete ic;
67  delete mybc;
68  }
70  {
71  BL_PROFILE("PhaseFieldMicrostructure::Parse");
72 
73  // Number of grain fields (may be more if using different IC)
74  pp_query_default("pf.number_of_grains", value.number_of_grains,2);
75  pp_query_required("pf.M", value.pf.M); // Mobility
76  pp_query_required("pf.mu", value.pf.mu); // Phase field :math:`\mu`
77  pp_query_required("pf.gamma", value.pf.gamma); // Phase field :math:`\gamma`
78  pp_query_required("pf.sigma0", value.pf.sigma0); // Initial GB energy if not using anisotropy
79  pp_query_required("pf.l_gb", value.pf.l_gb); // Mobility
80  pp_query_default("pf.elastic_df",value.pf.elastic_df,false); // Determine whether to use elastic driving force
81  pp_query_default("pf.elastic_mult",value.pf.elastic_mult,1.0); // Multiplier of elastic energy
82 
83  pp_query_default("pf.threshold.value",value.pf.threshold.value,0.0); // Value used for thresholding kinetic relation
84  pp_query_default("pf.threshold.chempot",value.pf.threshold.chempot,false); // Whether to include chemical potential in threshold
85  pp_query_default("pf.threshold.boundary",value.pf.threshold.boundary,false); // Whether to include boundary energy in threshold
86  pp_query_default("pf.threshold.corner",value.pf.threshold.corner,false); // Whether to include corner regularization in threshold
87  pp_query_default("pf.threshold.lagrange",value.pf.threshold.lagrange,false); // Whether to include lagrange multiplier in threshold
88  pp_query_default("pf.threshold.mechanics",value.pf.threshold.mechanics,false);// Whether to include mechanical driving force in threshold
89  {
90  value.pf.threshold.on =
91  value.pf.threshold.chempot || value.pf.threshold.boundary ||
92  value.pf.threshold.corner || value.pf.threshold.lagrange ||
93  value.pf.threshold.mechanics;
94 
95  std::string type = "continuous";
96  pp_query_validate("pf.threshold.type",type,{"continuous","chop"}); // Type of thresholding to use
97  if (type == "continuous") value.pf.threshold.type = ThresholdType::Continuous;
98  else if (type == "chop") value.pf.threshold.type = ThresholdType::Chop;
99  }
100 
101  value.pf.L = (4./3.)*value.pf.M / value.pf.l_gb;
102 
103  pp_query_required("amr.max_level", value.max_level); // Maximum AMR level
104  pp_query_default("amr.ref_threshold", value.ref_threshold, 0.1); // Phase field refinement threshold
105 
106  // Elasticity
107  pp_query_default("mechanics.tstart",value.mechanics.tstart, 0.0);
108  value.mechanics.model.clear();
109  value.mechanics.model.resize(value.number_of_grains,model_type::Zero());
110  for (int i=0; i < value.number_of_grains; i++)
111  {
112  // By default, read in the model specified by "mechanics.model"
113  if (pp.getEntries("mechanics.model").size())
114  pp_queryclass("mechanics.model", value.mechanics.model[i]);
115 
116  // If individual models are specified, read those in and overwrite
117  std::string name = "mechanics.model" + std::to_string(i+1);
118  if (pp.getEntries(name).size())
119  pp_queryclass(std::string(name.data()), value.mechanics.model[i]);
120  }
121 
122  value.m_type = Base::Mechanics<model_type>::Disable; // Turn mechanics off by default
125  value.number_of_ghost_cells = std::max(value.number_of_ghost_cells,3);
126 
127 
128  // Lagrange multiplier method for enforcing volumes
129  pp_query_default("lagrange.on", value.lagrange.on,false);
130  if (value.lagrange.on)
131  {
132  pp_query_required("lagrange.lambda", value.lagrange.lambda); // Lagrange multiplier value
133  pp_query_required("lagrange.vol0", value.lagrange.vol0); // Prescribed volume
134  pp_query_default("lagrange.tstart", value.lagrange.tstart,0.0); // Lagrange multipler start time
135  value.SetThermoInt(1);
136  }
137 
138  pp_query_default("sdf.on",value.sdf.on,false); // synthetic driving force (SDF)
139  if(value.sdf.on)
140  {
141  std::vector<std::string> vals;
142  pp_queryarr("sdf.val",vals); // value of SDF for each grain
143  int nvals = static_cast<int>(vals.size());
144  if (nvals == 1)
145  for (int i = 0; i < value.number_of_grains; i++)
146  value.sdf.val.push_back(Numeric::Interpolator::Linear<Set::Scalar>(vals[0]));
147  else if (nvals == value.number_of_grains)
148  for (int i = 0; i < value.number_of_grains; i++)
149  value.sdf.val.push_back(Numeric::Interpolator::Linear<Set::Scalar>(vals[i]));
150  else
151  Util::Abort(INFO,"sdf.val received ", vals.size(), " but requires 1 or ", value.number_of_grains);
152 
153  pp_query_default("sdf.tstart",value.sdf.tstart,0.0); // time to begin applying SDF
154  }
155 
156  // Anisotropic grain boundary energy parameters
157 
158  pp_query_default("anisotropy.on", value.anisotropy.on,false); // Turn on
159  if (value.anisotropy.on)
160  {
161  // Regularization para m
162  pp_query_required("anisotropy.beta", value.anisotropy.beta);
163  // Time to turn on anisotropy
164  pp_query_required("anisotropy.tstart", value.anisotropy.tstart);
165  value.anisotropy.timestep = value.timestep;
166  // Modify timestep when turned on
167  pp_query_required("anisotropy.timestep", value.anisotropy.timestep);
168  value.anisotropy.plot_int = value.plot_int;
169  // Modify plot_int when turned on
170  pp_query_default("anisotropy.plot_int", value.anisotropy.plot_int, -1);
171  value.anisotropy.plot_dt = value.plot_dt;
172  // Modify plot_dt when turned on
173  pp_query_default("anisotropy.plot_dt", value.anisotropy.plot_dt, -1.0);
174  // Modify thermo int when turned on
175  pp_query_default("anisotropy.thermo_int", value.anisotropy.thermo_int, -1);
176  // Modify thermo plot int when turned on
177  pp_query_default("anisotropy.thermo_plot_int", value.anisotropy.thermo_plot_int, -1);
178  // Frequency of elastic calculation
179  pp_query_default("anisotropy.elastic_int",value.anisotropy.elastic_int, -1);
180  if (value.anisotropy.on)
181  value.number_of_ghost_cells = std::max(value.number_of_ghost_cells,2);
182 
183  // Determine the kind of regularization to use
184  std::map<std::string, RegularizationType> regularization_type;
185  regularization_type["wilmore"] = RegularizationType::Wilmore;
186  regularization_type["k12"] = RegularizationType::K12;
187 
188  std::string regularization_type_input;
189  // Type of regularization to use
190  pp_query_validate("anisotropy.regularization", regularization_type_input,{"k12","wilmore"});
191  value.regularization = regularization_type[regularization_type_input];
192 
193  pp_forbid("anisotropy.gb_type"," --> anisotropy.type");
194  // Type of GB to use
196  }
197 
198  // Boundary condition for eta
199  pp.select<BC::Constant>("bc.eta",value.mybc,value.number_of_grains);
200 
201  // Initial condition for the order parameter eta
203 
204  // Anisotropic mobility
205  pp_query_default("anisotropic_kinetics.on",value.anisotropic_kinetics.on, 0);
206  if (value.anisotropic_kinetics.on)
207  {
208  // simulation time when anisotropic kinetics is activated
209  pp_query_default("anisotropic_kinetics.tstart",value.anisotropic_kinetics.tstart, 0.0);
210  std::string mobility_filename, threshold_filename;
211  // file containing anisotropic mobility data
212  pp_query_file("anisotropic_kinetics.mobility",mobility_filename);
213  value.anisotropic_kinetics.mobility = Numeric::Interpolator::Linear<Set::Scalar>::Read(mobility_filename);
214  // file containing anisotropic mobility data
215  pp_query_file("anisotropic_kinetics.threshold",threshold_filename);
216  value.anisotropic_kinetics.threshold = Numeric::Interpolator::Linear<Set::Scalar>::Read(threshold_filename);
217  value.RegisterNewFab(value.anisotropic_kinetics.L_mf, value.mybc, value.number_of_grains, 0, "mobility",true);
218  value.RegisterNewFab(value.anisotropic_kinetics.threshold_mf, value.mybc, value.number_of_grains, 0, "theshold",true);
219  }
220 
221 
222 
223 
224  value.RegisterNewFab(value.eta_mf, value.mybc, value.number_of_grains, value.number_of_ghost_cells, "Eta",true);
225  value.RegisterNewFab(value.driving_force_mf, value.mybc, value.number_of_grains, value.number_of_ghost_cells, "DrivingForce",false);
226  if (value.pf.threshold.on)
227  value.RegisterNewFab(value.driving_force_threshold_mf, value.mybc, value.number_of_grains, value.number_of_ghost_cells, "DrivingForceThreshold",false);
228 
229  value.RegisterIntegratedVariable(&value.volume, "volume");
230  value.RegisterIntegratedVariable(&value.area, "area");
231  value.RegisterIntegratedVariable(&value.gbenergy, "gbenergy");
232  value.RegisterIntegratedVariable(&value.realgbenergy, "realgbenergy");
233  value.RegisterIntegratedVariable(&value.regenergy, "regenergy");
234 
235  }
236 
237 
238 protected:
239 
240  /// \fn Advance
241  /// \brief Evolve phase field in time
242  void Advance (int lev, Real time, Real dt) override;
243  void Initialize (int lev) override;
244 
245  void TagCellsForRefinement (int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow) override;
246 
247  virtual void TimeStepBegin(amrex::Real time, int iter) override;
248  virtual void TimeStepComplete(amrex::Real time, int iter) override;
249  void Integrate(int amrlev, Set::Scalar time, int step,
250  const amrex::MFIter &mfi, const amrex::Box &box) override;
251 
252 
253  virtual void UpdateModel(int /*a_step*/, Set::Scalar /*a_time*/) override;
254 
255 private:
256 
260 
261  // Cell fab
262  Set::Field<Set::Scalar> eta_mf; // Multicomponent field variable storing \t$\eta_i\t$ for the __current__ timestep
265 
267 
268  //amrex::Real M, mu, gamma, sigma0, l_gb, beta;
271  Continuous = 0, Chop = 1
272  };
273  struct {
274  Set::Scalar M = NAN;
275  Set::Scalar L = NAN;
280  bool elastic_df = false;
282  struct {
283  bool on = false;
284  bool chempot = false;
285  bool boundary = false;
286  bool corner = false;
287  bool lagrange = false;
288  bool mechanics = false;
289  bool sdf = false;
291  ThresholdType type = ThresholdType::Continuous;
292  } threshold;
293  } pf;
294 
295  struct {
296  int on = 0;
303 
304  struct {
305  int on = 0;
309  int plot_int = -1;
311  int thermo_int = -1, thermo_plot_int = -1;
314  int elastic_int = -1;
315  } anisotropy;
316 
317  struct {
318  bool on = 0;
319  Set::Scalar tstart = NAN;
322  } lagrange;
323 
324  struct {
325  int on = 0;
326  std::vector<Numeric::Interpolator::Linear<Set::Scalar>> val;
327  Set::Scalar tstart = 0.0;
328  } sdf;
329 
330  std::string gb_type, filename;
331 
333 
334  IC::IC *ic = nullptr;
335 
341 
342  struct
343  {
344  Set::Scalar tstart = 0.0;
345  std::vector<model_type> model;
346  } mechanics;
347 
348 
349 };
350 }
351 #endif
Integrator::PhaseFieldMicrostructure::plot_dt
Set::Scalar plot_dt
Definition: PhaseFieldMicrostructure.H:310
Integrator::PhaseFieldMicrostructure::lambda
Set::Scalar lambda
Definition: PhaseFieldMicrostructure.H:321
SH.H
Integrator::PhaseFieldMicrostructure::threshold
Numeric::Interpolator::Linear< Set::Scalar > threshold
Definition: PhaseFieldMicrostructure.H:299
Integrator::PhaseFieldMicrostructure::value
Set::Scalar value
Definition: PhaseFieldMicrostructure.H:290
Integrator::PhaseFieldMicrostructure::mechanics
bool mechanics
Definition: PhaseFieldMicrostructure.H:288
pp_forbid
#define pp_forbid(...)
Definition: ParmParse.H:107
IC::IC
Pure abstract IC object from which all other IC objects inherit.
Definition: IC.H:21
BC::BC< Set::Scalar >
Integrator::PhaseFieldMicrostructure::filename
std::string filename
Definition: PhaseFieldMicrostructure.H:330
Integrator::K12
@ K12
Definition: PhaseFieldMicrostructure.H:47
Integrator::PhaseFieldMicrostructure::driving_force_mf
Set::Field< Set::Scalar > driving_force_mf
Definition: PhaseFieldMicrostructure.H:263
IC::Ellipse
Definition: Ellipse.H:15
Integrator::Base::Mechanics
Definition: Mechanics.H:23
Integrator::Integrator::box
std::vector< amrex::Box > box
Definition: Integrator.H:463
PerturbedInterface.H
Integrator::PhaseFieldMicrostructure::driving_force_threshold_mf
Set::Field< Set::Scalar > driving_force_threshold_mf
Definition: PhaseFieldMicrostructure.H:264
Integrator::PhaseFieldMicrostructure::thermo_int
int thermo_int
Definition: PhaseFieldMicrostructure.H:311
Mechanics.H
Model::Interface::GB::GB
Definition: GB.H:18
Integrator::PhaseFieldMicrostructure::gamma
Set::Scalar gamma
Definition: PhaseFieldMicrostructure.H:277
Integrator::PhaseFieldMicrostructure::elastic_mult
Set::Scalar elastic_mult
Definition: PhaseFieldMicrostructure.H:281
Read.H
Model::Interface::GB::Read
Definition: Read.H:23
Sphere.H
Integrator::PhaseFieldMicrostructure::mybc
BC::BC< Set::Scalar > * mybc
Definition: PhaseFieldMicrostructure.H:266
Integrator.H
Integrator::PhaseFieldMicrostructure::number_of_ghost_cells
int number_of_ghost_cells
Definition: PhaseFieldMicrostructure.H:258
Integrator::PhaseFieldMicrostructure::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Definition: PhaseFieldMicrostructure.cpp:281
Set::Field< Set::Scalar >
Definition: Set.H:236
Integrator::PhaseFieldMicrostructure::phi0
Set::Scalar phi0
Definition: PhaseFieldMicrostructure.H:313
Integrator::PhaseFieldMicrostructure::threshold
struct Integrator::PhaseFieldMicrostructure::@16::@22 threshold
Integrator::PhaseFieldMicrostructure::val
std::vector< Numeric::Interpolator::Linear< Set::Scalar > > val
Definition: PhaseFieldMicrostructure.H:326
Integrator::PhaseFieldMicrostructure::mu
Set::Scalar mu
Definition: PhaseFieldMicrostructure.H:276
Integrator::RegularizationType
RegularizationType
Definition: PhaseFieldMicrostructure.H:45
Integrator::PhaseFieldMicrostructure::ref_threshold
Set::Scalar ref_threshold
Definition: PhaseFieldMicrostructure.H:259
Expression.H
Integrator::PhaseFieldMicrostructure::eta_mf
Set::Field< Set::Scalar > eta_mf
Definition: PhaseFieldMicrostructure.H:262
Integrator::PhaseFieldMicrostructure::Chop
@ Chop
Definition: PhaseFieldMicrostructure.H:271
Integrator::PhaseFieldMicrostructure::elastic_df
bool elastic_df
Definition: PhaseFieldMicrostructure.H:280
Integrator::PhaseFieldMicrostructure::ic
IC::IC * ic
Definition: PhaseFieldMicrostructure.H:334
Integrator::PhaseFieldMicrostructure::plot_int
int plot_int
Definition: PhaseFieldMicrostructure.H:309
Constant.H
Constant.H
Integrator::PhaseFieldMicrostructure
Microstructure evolution with grain boundary anisotropy.
Definition: PhaseFieldMicrostructure.H:57
Integrator::PhaseFieldMicrostructure::UpdateModel
virtual void UpdateModel(int, Set::Scalar) override
Definition: PhaseFieldMicrostructure.cpp:312
Cubic.H
Elastic.H
Integrator::PhaseFieldMicrostructure::pf
struct Integrator::PhaseFieldMicrostructure::@16 pf
Integrator::PhaseFieldMicrostructure::M
Set::Scalar M
Definition: PhaseFieldMicrostructure.H:274
Integrator::PhaseFieldMicrostructure::L
Set::Scalar L
Definition: PhaseFieldMicrostructure.H:275
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
Integrator::PhaseFieldMicrostructure::tstart
Set::Scalar tstart
Definition: PhaseFieldMicrostructure.H:297
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:106
Integrator::PhaseFieldMicrostructure::number_of_grains
int number_of_grains
Definition: PhaseFieldMicrostructure.H:257
Integrator::PhaseFieldMicrostructure::corner
bool corner
Definition: PhaseFieldMicrostructure.H:286
Model::Interface::GB::Sin
Definition: Sin.H:18
AbsSin.H
Integrator::PhaseFieldMicrostructure::chempot
bool chempot
Definition: PhaseFieldMicrostructure.H:284
Integrator::PhaseFieldMicrostructure::elastic_int
int elastic_int
Definition: PhaseFieldMicrostructure.H:314
Integrator::PhaseFieldMicrostructure::on
bool on
Definition: PhaseFieldMicrostructure.H:283
IC::Voronoi
Definition: Voronoi.H:10
Integrator::PhaseFieldMicrostructure::lagrange
bool lagrange
Definition: PhaseFieldMicrostructure.H:287
Integrator::PhaseFieldMicrostructure::sigma1
Set::Scalar sigma1
Definition: PhaseFieldMicrostructure.H:312
Constant.H
TabulatedInterface.H
Integrator::PhaseFieldMicrostructure::Integrate
void Integrate(int amrlev, Set::Scalar time, int step, const amrex::MFIter &mfi, const amrex::Box &box) override
Definition: PhaseFieldMicrostructure.cpp:360
IC::PerturbedInterface
Definition: PerturbedInterface.H:27
Integrator::PhaseFieldMicrostructure::theta0
Set::Scalar theta0
Definition: PhaseFieldMicrostructure.H:312
Integrator::PhaseFieldMicrostructure::sdf
bool sdf
Definition: PhaseFieldMicrostructure.H:289
Integrator::PhaseFieldMicrostructure::Parse
static void Parse(PhaseFieldMicrostructure &value, IO::ParmParse &pp)
Definition: PhaseFieldMicrostructure.H:69
Sin.H
Integrator::PhaseFieldMicrostructure::l_gb
Set::Scalar l_gb
Definition: PhaseFieldMicrostructure.H:279
pp_query_file
#define pp_query_file(...)
Definition: ParmParse.H:102
Integrator::PhaseFieldMicrostructure::L_mf
Set::Field< Set::Scalar > L_mf
Definition: PhaseFieldMicrostructure.H:300
Model::Interface::GB::SH
Definition: SH.H:35
Cubic.H
Integrator::PhaseFieldMicrostructure::area
Set::Scalar area
Definition: PhaseFieldMicrostructure.H:337
pp_query_default
#define pp_query_default(...)
Definition: ParmParse.H:100
Integrator::PhaseFieldMicrostructure::type
ThresholdType type
Definition: PhaseFieldMicrostructure.H:291
Numeric::Interpolator::Linear< Set::Scalar >
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:166
Integrator
Collection of numerical integrator objects.
Definition: AllenCahn.H:40
IC::Constant
Definition: Constant.H:18
Random.H
Integrator::PhaseFieldMicrostructure::gbenergy
Set::Scalar gbenergy
Definition: PhaseFieldMicrostructure.H:338
Integrator::PhaseFieldMicrostructure::volume
Set::Scalar volume
Definition: PhaseFieldMicrostructure.H:336
Numeric::Interpolator::Linear::Read
static Linear< T > Read(std::string filename, int derivative=0)
BC::Constant
Definition: Constant.H:43
IC::Expression
Definition: Expression.H:44
IC::Random
Set each point to a random value.
Definition: Random.H:15
Integrator::Wilmore
@ Wilmore
Definition: PhaseFieldMicrostructure.H:46
Integrator::PhaseFieldMicrostructure::PhaseFieldMicrostructure
PhaseFieldMicrostructure(IO::ParmParse &pp)
Definition: PhaseFieldMicrostructure.H:61
Integrator::PhaseFieldMicrostructure::mobility
Numeric::Interpolator::Linear< Set::Scalar > mobility
Definition: PhaseFieldMicrostructure.H:298
Voronoi.H
Integrator::PhaseFieldMicrostructure::anisotropic_kinetics
struct Integrator::PhaseFieldMicrostructure::@17 anisotropic_kinetics
pp_query_validate
#define pp_query_validate(...)
Definition: ParmParse.H:101
Integrator::PhaseFieldMicrostructure::thermo_plot_int
int thermo_plot_int
Definition: PhaseFieldMicrostructure.H:311
pp_query_required
#define pp_query_required(...)
Definition: ParmParse.H:99
Integrator::PhaseFieldMicrostructure::gb_type
std::string gb_type
Definition: PhaseFieldMicrostructure.H:330
Integrator::PhaseFieldMicrostructure::Advance
void Advance(int lev, Real time, Real dt) override
Evolve phase field in time.
Definition: PhaseFieldMicrostructure.cpp:29
BC.H
IO::ParmParse
Definition: ParmParse.H:112
Integrator::PhaseFieldMicrostructure::timestep
Set::Scalar timestep
Definition: PhaseFieldMicrostructure.H:307
Integrator::PhaseFieldMicrostructure::realgbenergy
Set::Scalar realgbenergy
Definition: PhaseFieldMicrostructure.H:339
Integrator::PhaseFieldMicrostructure::regenergy
Set::Scalar regenergy
Definition: PhaseFieldMicrostructure.H:340
Integrator::PhaseFieldMicrostructure::anisotropy
struct Integrator::PhaseFieldMicrostructure::@18 anisotropy
Integrator::PhaseFieldMicrostructure::Initialize
void Initialize(int lev) override
Definition: PhaseFieldMicrostructure.cpp:273
Integrator::PhaseFieldMicrostructure::threshold_mf
Set::Field< Set::Scalar > threshold_mf
Definition: PhaseFieldMicrostructure.H:301
IO::ParmParse::queryclass
void queryclass(std::string name, T *value, std::string file="", std::string func="", int line=-1)
Definition: ParmParse.H:455
Integrator::PhaseFieldMicrostructure::model
std::vector< model_type > model
Definition: PhaseFieldMicrostructure.H:345
Integrator::PhaseFieldMicrostructure::vol0
Set::Scalar vol0
Definition: PhaseFieldMicrostructure.H:320
Integrator::PhaseFieldMicrostructure::sigma0
Set::Scalar sigma0
Definition: PhaseFieldMicrostructure.H:278
Integrator::PhaseFieldMicrostructure::beta
Set::Scalar beta
Definition: PhaseFieldMicrostructure.H:306
Integrator::PhaseFieldMicrostructure::Continuous
@ Continuous
Definition: PhaseFieldMicrostructure.H:271
Integrator::PhaseFieldMicrostructure::TimeStepBegin
virtual void TimeStepBegin(amrex::Real time, int iter) override
Definition: PhaseFieldMicrostructure.cpp:346
INFO
#define INFO
Definition: Util.H:20
Integrator::Integrator::dt
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition: Integrator.H:398
Integrator::PhaseFieldMicrostructure::TimeStepComplete
virtual void TimeStepComplete(amrex::Real time, int iter) override
Definition: PhaseFieldMicrostructure.cpp:306
Integrator::PhaseFieldMicrostructure::regularization
RegularizationType regularization
Definition: PhaseFieldMicrostructure.H:269
Model::Interface::GB::AbsSin
Definition: AbsSin.H:18
IC::Sphere
Definition: Sphere.H:13
Integrator::PhaseFieldMicrostructure::~PhaseFieldMicrostructure
virtual ~PhaseFieldMicrostructure()
Definition: PhaseFieldMicrostructure.H:63
Integrator::PhaseFieldMicrostructure::ThresholdType
ThresholdType
Definition: PhaseFieldMicrostructure.H:270
IO::ParmParse::select
void select(std::string name, PTRTYPE *&ic_eta, Args &&... args)
Definition: ParmParse.H:541
Integrator::PhaseFieldMicrostructure::boundary
bool boundary
Definition: PhaseFieldMicrostructure.H:285
Ellipse.H
Integrator::PhaseFieldMicrostructure::PhaseFieldMicrostructure
PhaseFieldMicrostructure()
Definition: PhaseFieldMicrostructure.H:60
pp_queryarr
#define pp_queryarr(...)
Definition: ParmParse.H:103