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 <random>
8
9#include "AMReX_ParmParse.H"
10#include <AMReX_MLMG.H>
11
13#include "BC/BC.H"
14#include "BC/Constant.H"
15#include "IC/Constant.H"
17#include "IC/Voronoi.H"
18#include "IC/Sphere.H"
19#include "IC/Expression.H"
20#include "IC/Ellipse.H"
21#include "IC/Random.H"
22#include "Model/Interface/GB/GB.H"
28#include "Base/Mechanics.H"
29
30namespace Integrator
31{
32
36
37///
38/// Solve the Allen-Cahn evolution equation for microstructure with parameters \f$\eta_1\ldots\eta_n\f$,
39/// where n corresponds to the number of grains.
40///
41template<class model_type>
42class PhaseFieldMicrostructure : public Base::Mechanics<model_type>
43{
44public:
45 const std::string name = "phasefieldmicrostructure." + std::string(model_type::name);
46
51 {
52 delete boundary;
53 delete ic;
54 delete mybc;
55 }
57 {
58 BL_PROFILE("PhaseFieldMicrostructure::Parse");
59
60 // Number of grain fields (may be more if using different IC)
61 pp_query_default("pf.number_of_grains", value.number_of_grains,2);
62 pp_query_required("pf.M", value.pf.M); // Mobility
63 pp_query_required("pf.gamma", value.pf.gamma); // Phase field :math:`\gamma`
64 pp_query_required("pf.sigma0", value.pf.sigma0); // Initial GB energy if not using anisotropy
65 pp_query_required("pf.l_gb", value.pf.l_gb); // Mobility
66 pp_query_default("pf.elastic_df",value.pf.elastic_df,false); // Determine whether to use elastic driving force
67 pp_query_default("pf.elastic_mult",value.pf.elastic_mult,1.0); // Multiplier of elastic energy
68
69 pp_query_default("pf.threshold.value",value.pf.threshold.value,0.0); // Value used for thresholding kinetic relation
70 pp_query_default("pf.threshold.chempot",value.pf.threshold.chempot,false); // Whether to include chemical potential in threshold
71 pp_query_default("pf.threshold.boundary",value.pf.threshold.boundary,false); // Whether to include boundary energy in threshold
72 pp_query_default("pf.threshold.corner",value.pf.threshold.corner,false); // Whether to include corner regularization in threshold
73 pp_query_default("pf.threshold.lagrange",value.pf.threshold.lagrange,false); // Whether to include lagrange multiplier in threshold
74 pp_query_default("pf.threshold.mechanics",value.pf.threshold.mechanics,false);// Whether to include mechanical driving force in threshold
75 {
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;
80
81 std::string type = "continuous";
82 pp_query_validate("pf.threshold.type",type,{"continuous","chop"}); // Type of thresholding to use
83 if (type == "continuous") value.pf.threshold.type = ThresholdType::Continuous;
84 else if (type == "chop") value.pf.threshold.type = ThresholdType::Chop;
85 }
86
87 value.pf.L = (4./3.)*value.pf.M / value.pf.l_gb;
88
89 pp_query_required("amr.max_level", value.max_level); // Maximum AMR level
90 pp_query_default("amr.ref_threshold", value.ref_threshold, 0.1); // Phase field refinement threshold
91
92
93 std::string type_str;
94 // Reading this is redundant but necessary because of the way the code was originally structured
95 // (need to fix eventually)
96 pp.query_validate("mechanics.type", type_str, {"disable","static","dynamic"});
97 value.m_type = Base::Mechanics<model_type>::Disable; // Turn mechanics off by default
98
99 if (type_str != "disable") // do this only if mechanics is activated
100 {
101 // Elasticity
102 pp_query_default("mechanics.tstart",value.mechanics.tstart, 0.0);
103
104 // Read in models - either one model for all grains, or
105 // individual models, one for each grain.
106 pp.queryclass_enumerate("mechanics.model",value.mechanics.model, value.number_of_grains);
107
108 // Mixing order
109 pp_query_validate("mechanics.mix_order",value.mechanics.model_mix_order,{1,2});
110 // Force Neumann BCs on the model
111 pp_query_default("mechanics.model_neuman_boundary",value.mechanics.model_neumann_boundary,false);
112
115 value.number_of_ghost_cells = std::max(value.number_of_ghost_cells, 3);
116 }
117
118
119 // Lagrange multiplier method for enforcing volumes
120 pp_query_default("lagrange.on", value.lagrange.on,false);
121 if (value.lagrange.on)
122 {
123 pp_query_required("lagrange.lambda", value.lagrange.lambda); // Lagrange multiplier value
124 pp_query_required("lagrange.vol0", value.lagrange.vol0); // Prescribed volume
125 pp_query_default("lagrange.tstart", value.lagrange.tstart,0.0); // Lagrange multipler start time
126 value.SetThermoInt(1);
127 }
128
129 pp_query_default("sdf.on",value.sdf.on,false); // synthetic driving force (SDF)
130 if(value.sdf.on)
131 {
132 std::vector<std::string> vals;
133 pp_queryarr("sdf.val",vals); // value of SDF for each grain
134 int nvals = static_cast<int>(vals.size());
135 if (nvals == 1)
136 for (int i = 0; i < value.number_of_grains; i++)
137 value.sdf.val.push_back(Numeric::Interpolator::Linear<Set::Scalar>(vals[0]));
138 else if (nvals == value.number_of_grains)
139 for (int i = 0; i < value.number_of_grains; i++)
140 value.sdf.val.push_back(Numeric::Interpolator::Linear<Set::Scalar>(vals[i]));
141 else
142 Util::Abort(INFO,"sdf.val received ", vals.size(), " but requires 1 or ", value.number_of_grains);
143
144 pp_query_default("sdf.tstart",value.sdf.tstart,0.0); // time to begin applying SDF
145 }
146
147 // Anisotropic grain boundary energy parameters
148
149 pp_query_default("anisotropy.on", value.anisotropy.on,false); // Turn on
150 if (value.anisotropy.on)
151 {
152 // Regularization para m
153 pp_query_required("anisotropy.beta", value.anisotropy.beta);
154 // Time to turn on anisotropy
155 pp_query_required("anisotropy.tstart", value.anisotropy.tstart);
156 value.anisotropy.timestep = value.timestep;
157 // Modify timestep when turned on
158 pp_query_required("anisotropy.timestep", value.anisotropy.timestep);
159 value.anisotropy.plot_int = value.plot_int;
160 // Modify plot_int when turned on
161 pp_query_default("anisotropy.plot_int", value.anisotropy.plot_int, -1);
162 value.anisotropy.plot_dt = value.plot_dt;
163 // Modify plot_dt when turned on
164 pp_query_default("anisotropy.plot_dt", value.anisotropy.plot_dt, -1.0);
165 // Modify thermo int when turned on
166 pp_query_default("anisotropy.thermo_int", value.anisotropy.thermo_int, -1);
167 // Modify thermo plot int when turned on
168 pp_query_default("anisotropy.thermo_plot_int", value.anisotropy.thermo_plot_int, -1);
169 // Frequency of elastic calculation
170 pp_query_default("anisotropy.elastic_int",value.anisotropy.elastic_int, -1);
171 if (value.anisotropy.on)
172 value.number_of_ghost_cells = std::max(value.number_of_ghost_cells,2);
173
174 // Determine the kind of regularization to use
175 std::map<std::string, RegularizationType> regularization_type;
176 regularization_type["wilmore"] = RegularizationType::Wilmore;
177 regularization_type["k12"] = RegularizationType::K12;
178
179 std::string regularization_type_input;
180 // Type of regularization to use
181 pp_query_validate("anisotropy.regularization", regularization_type_input,{"k12","wilmore"});
182 value.regularization = regularization_type[regularization_type_input];
183
184 pp_forbid("anisotropy.gb_type"," --> anisotropy.type");
185 // Type of GB to use
187 }
188 // Thermal fluctuations
189 pp.query("fluctuation.on",value.fluctuation.on);
190 if (value.fluctuation.on)
191 {
192 pp.query("fluctuation.amp",value.fluctuation.amp); // fluctuation amplitude
193 pp.query("fluctuation.sd",value.fluctuation.sd); // fluctuation stadard deviation
194 pp.query("fluctuation.tstart", value.fluctuation.tstart); // time to start applying fluctuation
195 value.fluctuation.norm_dist = std::normal_distribution<double>(0.0,value.fluctuation.sd);
196 }
197
198 // Disconnection generation
199 pp.query("disconnection.on",value.disconnection.on);
200 if (value.disconnection.on)
201 {
202 // Read in nucleation parameters from disconnection class
203 pp.queryclass<Model::Defect::Disconnection>("disconnection",value.disconnection.model);
204 }
205
206 // Shear coupling matrices
207 pp_query_default("shearcouple.on",value.shearcouple.on,false);
208 if (value.shearcouple.on)
209 {
210 Util::AssertException(INFO,TEST(value.m_time_evolving == true), " mechanics.time_evolving must be true when using shearcouple");
211 Util::AssertException(INFO,TEST(value.mechanics.model_mix_order == 2), " mechanics.model.mix_order must be 2 when using shearcouple");
212
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++)
216 {
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);
219 if ( i==j && pp.contains(name.data()))
220 Util::Abort(INFO,"Cannot specify self FGB ", name);
221 if (pp.contains(name.data()) && pp.contains(namerev.data()))
222 Util::Abort(INFO,"Cannot specify both ",name," and ",namerev);
223 if (pp.contains(name.data()))
224 {
225 pp.queryarr(name.data(),value.shearcouple.Fgb[i*value.number_of_grains + j]);
226 value.shearcouple.Fgb[j*value.number_of_grains + i] = - value.shearcouple.Fgb[i*value.number_of_grains + j];
227 }
228 }
229 }
230
231 // Boundary condition for eta
232 pp.select<BC::Constant>("bc.eta",value.mybc,value.number_of_grains);
233
234 // Initial condition for the order parameter eta
236
237 // Anisotropic mobility
238 pp_query_default("anisotropic_kinetics.on",value.anisotropic_kinetics.on, 0);
239 if (value.anisotropic_kinetics.on)
240 {
241 // simulation time when anisotropic kinetics is activated
242 pp_query_default("anisotropic_kinetics.tstart",value.anisotropic_kinetics.tstart, 0.0);
243 std::string mobility_filename, threshold_filename;
244 // file containing anisotropic mobility data
245 pp_query_file("anisotropic_kinetics.mobility",mobility_filename);
246 value.anisotropic_kinetics.mobility = Numeric::Interpolator::Linear<Set::Scalar>::Read(mobility_filename);
247 // file containing anisotropic mobility data
248 pp_query_file("anisotropic_kinetics.threshold",threshold_filename);
249 value.anisotropic_kinetics.threshold = Numeric::Interpolator::Linear<Set::Scalar>::Read(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);
252 }
253
254
255
256
257 value.RegisterNewFab(value.eta_mf, value.mybc, value.number_of_grains, value.number_of_ghost_cells, "Eta",true);
258 value.RegisterNewFab(value.eta_old_mf, value.mybc, value.number_of_grains, value.number_of_ghost_cells, "EtaOld",false);
259 //value.RegisterNewFab(value.driving_force_mf, value.mybc, value.number_of_grains, value.number_of_ghost_cells, "DrivingForce",false);
260 // if (value.pf.threshold.on)
261 // value.RegisterNewFab(value.driving_force_threshold_mf, value.mybc, value.number_of_grains, value.number_of_ghost_cells, "DrivingForceThreshold",false);
262 // if (value.disconnection.on)
263 // value.RegisterNewFab(value.disc_mf, new BC::Nothing(), 1, value.number_of_ghost_cells, "disc",true); // see box
264
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");
270
271 }
272
273
274protected:
275
276 /// \fn Advance
277 /// \brief Evolve phase field in time
278 void Advance (int lev, Real time, Real dt) override;
279 void Initialize (int lev) override;
280
281 void TagCellsForRefinement (int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow) override;
282
283 virtual void TimeStepBegin(amrex::Real time, int iter) override;
284 virtual void TimeStepComplete(amrex::Real time, int iter) override;
285 void Integrate(int amrlev, Set::Scalar time, int step,
286 const amrex::MFIter &mfi, const amrex::Box &box) override;
287
288 virtual void UpdateEigenstrain(int lev);
289 virtual void UpdateEigenstrain()
290 {
291 for (int lev = 0; lev <= this->max_level; lev++)
293 }
294
295 virtual void UpdateModel(int /*a_step*/, Set::Scalar /*a_time*/) override;
296
297private:
298
302
303 // Cell fab
304 Set::Field<Set::Scalar> eta_mf; // Multicomponent field variable storing \t$\eta_i\t$ for the __current__ timestep
309 //Set::Field<Set::Scalar> disc_mf; //see box
310 // Node fab
311 //Set::Field<Set::Scalar> elasticdf_mf;
313
315
316 //amrex::Real M, mu, gamma, sigma0, l_gb, beta;
321 struct {
328 bool elastic_df = false;
330 struct {
331 bool on = false;
332 bool chempot = false;
333 bool boundary = false;
334 bool corner = false;
335 bool lagrange = false;
336 bool mechanics = false;
337 bool sdf = false;
341 } pf;
342
343 struct {
344 int on = 0;
351
352 struct {
353 int on = 0;
357 int plot_int = -1;
362 int elastic_int = -1;
364
365 struct {
366 bool on = 0;
367 Set::Scalar tstart = NAN;
371
372 struct {
373 int on = 0;
374 std::vector<Numeric::Interpolator::Linear<Set::Scalar>> val;
375 Set::Scalar tstart = 0.0;
377
378 struct {
379 int on = 0;
382 Set::Scalar tstart = 0.0;
383 std::normal_distribution<double> norm_dist;
384 std::default_random_engine rand_num_gen;
386
387 struct {
388 int on = 0;
391
392 struct {
393 int on = 0;
394 std::vector<Set::Matrix> Fgb;
396
397 std::string gb_type, filename;
398
400
402
408
409 struct
410 {
411 Set::Scalar tstart = 0.0;
412 std::vector<model_type> model;
416
417 using Base::Mechanics<model_type>::model_mf;
418 using Base::Mechanics<model_type>::stress_mf;
419
420};
421}
422#endif
#define pp_query_validate(...)
Definition ParmParse.H:101
#define pp_queryarr(...)
Definition ParmParse.H:103
#define pp_query_required(...)
Definition ParmParse.H:99
#define pp_query_file(...)
Definition ParmParse.H:102
#define pp_query_default(...)
Definition ParmParse.H:100
#define pp_forbid(...)
Definition ParmParse.H:108
#define pp_queryclass(...)
Definition ParmParse.H:107
#define TEST(x)
Definition Util.H:21
#define INFO
Definition Util.H:20
Definition BC.H:42
Pure abstract IC object from which all other IC objects inherit.
Definition IC.H:23
Set each point to a random value.
Definition Random.H:20
bool contains(std::string name)
Definition ParmParse.H:154
void queryclass(std::string name, T *value, std::string file="", std::string func="", int line=-1)
Definition ParmParse.H:604
void select(std::string name, PTRTYPE *&ic_eta, Args &&... args)
Definition ParmParse.H:692
int queryclass_enumerate(std::string a_name, std::vector< T > &value, int number=1, std::string file="", std::string func="", int line=__LINE__)
Definition ParmParse.H:416
int queryarr(std::string name, std::vector< T > &value, std::string, std::string, int)
Definition ParmParse.H:336
int query_validate(std::string name, int &value, std::vector< int > possibleintvals, std::string file="", std::string func="", int line=-1)
Definition ParmParse.H:195
Set::Field< Set::Matrix > stress_mf
Definition Mechanics.H:474
Set::Field< model_type > model_mf
Definition Mechanics.H:463
std::vector< amrex::Box > box
Definition Integrator.H:446
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition Integrator.H:381
Solve the Allen-Cahn evolution equation for microstructure with parameters , where n corresponds to t...
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
struct Integrator::PhaseFieldMicrostructure::@16::@25 threshold
static void Parse(PhaseFieldMicrostructure &value, IO::ParmParse &pp)
struct Integrator::PhaseFieldMicrostructure::@21 fluctuation
void TagCellsForRefinement(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
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
struct Integrator::PhaseFieldMicrostructure::@23 shearcouple
struct Integrator::PhaseFieldMicrostructure::@22 disconnection
Numeric::Interpolator::Linear< Set::Scalar > threshold
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
std::vector< Numeric::Interpolator::Linear< Set::Scalar > > val
virtual void TimeStepComplete(amrex::Real time, int iter) override
Reads the data from a file and computes energies and its derivates.
Definition Read.H:24
A 2D interface model class.
Definition SH.H:36
static Linear< T > Read(std::string filename, int derivative=0)
Collection of numerical integrator objects.
Definition AllenCahn.H:41
amrex::Real Scalar
Definition Base.H:19
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:20
AMREX_FORCE_INLINE void AssertException(std::string file, std::string func, int line, std::string smt, bool pass, Args const &... args)
Definition Util.H:233
void Abort(const char *msg)
Definition Util.cpp:170