Alamo
AllenCahn.H
Go to the documentation of this file.
1 //
2 // This is a simple implementation of an Allen-Cahn equation governed by
3 //
4 // .. math::
5 //
6 // \frac{\partial\alpha}{\partial t} =
7 // - L \Big(\frac{\lambda}{\epsilon}(2.0\alpha
8 // - 6\alpha^2 + 4\alpha^3) + \epsilon\,\kappa\,\Delta \alpha \Big)
9 //
10 // where :math:`\alpha(x,t)` is the order parameter.
11 //
12 
13 #ifndef INTEGRATOR_ALLENCAHN_H // Include guards
14 #define INTEGRATOR_ALLENCAHN_H //
15 
16 // Standard library includes
17 #include <iostream>
18 #include <fstream>
19 #include <iomanip>
20 
21 // AMReX Includes
22 #include "AMReX.H"
23 #include "AMReX_ParallelDescriptor.H"
24 #include "AMReX_ParmParse.H"
25 
26 // Alamo Includes
27 #include "IO/ParmParse.H"
28 #include "Integrator/Integrator.H"
29 #include "BC/Constant.H"
30 #include "IC/IC.H"
31 #include "IC/Sphere.H"
32 #include "IC/Constant.H"
33 #include "IC/Expression.H"
34 #include "IC/BMP.H"
35 #include "IC/Random.H"
36 #include "Numeric/Stencil.H"
37 #include "IC/Random.H"
38 
39 namespace Integrator
40 {
41 class AllenCahn : virtual public Integrator
42 {
43 public:
44  // Empty constructor
46  {}
47 
48  // Constructor that triggers parse
50  {
51  Parse(*this, pp);
52  }
53 
54  // The Parse function initializes the AllenCahn object using
55  // a parser, pp.
56  // Note that this is a static function, which means it does not have
57  // direct access to member variables. Instead, it initializes the variables
58  // inside the argument, "value", and so all references to member items are
59  // prefixed by "value."
60  static void Parse(AllenCahn& value, IO::ParmParse& pp)
61  {
62  // Criterion for mesh refinement [0.01]
63  pp_query_default("refinement_threshold", value.refinement_threshold, 0.01);
64 
65  // Value for :math:`L` (mobility)
66  pp_query_default("ch.L",value.ch.L, 1.0);
67  // Value for :math:`\epsilon` (diffuse boundary width)
68  pp_query_default("ch.eps",value.ch.eps, 0.1);
69  // Value for :math:`\kappa` (Interface energy parameter)
70  pp_query_default("ch.grad",value.ch.grad, 1.0);
71  // Value for :math:`\lambda` (Chemical potential coefficient)
72  pp_query_default("ch.chempot",value.ch.chempot, 1.0);
73 
74  std::string type = "sphere";
75  // Initial condition type ([sphere], constant, expression, bmp, random)
76  pp_query_validate("alpha.ic.type", type, {"sphere","constant","expression","bmp","random"});
77  if (type == "sphere") value.ic = new IC::Sphere(value.geom, pp, "alpha.ic.sphere");
78  else if (type == "constant") value.ic = new IC::Constant(value.geom, pp, "alpha.ic.constant");
79  else if (type == "expression") value.ic = new IC::Expression(value.geom, pp, "alpha.ic.expression");
80  else if (type == "bmp") value.ic = new IC::BMP(value.geom, pp, "alpha.ic.bmp");
81  else if (type == "random") value.ic = new IC::Random(value.geom, pp, "alpha.ic.random");
82  else Util::Abort(INFO, "Invalid ic.type ", type);
83 
84  // Use a constant BC object for temperature
85  value.bc = new BC::Constant(1);
86  // :ref:`BC::Constant` parameters
87  pp_queryclass("alpha.bc", *static_cast<BC::Constant*>(value.bc));
88 
89  // Register the temperature and old temperature fields.
90  // alpha_mf and alpha_old_mf are defined near the bottom of this Header file.
91  value.RegisterNewFab(value.alpha_mf, value.bc, value.number_of_components, value.number_of_ghost_cells, "alpha", true);
92  value.RegisterNewFab(value.alpha_old_mf, value.bc, value.number_of_components, value.number_of_ghost_cells, "alpha_old", false);
93  }
94 
95 protected:
96 
97  // Use the ic object to initialize the order parameter
98  void Initialize(int lev)
99  {
100  ic->Initialize(lev, alpha_mf);
101  ic->Initialize(lev, alpha_old_mf);
102  }
103 
104  // Integrate the Allen Cahn equation
105  void Advance(int lev, Set::Scalar /*time*/, Set::Scalar dt)
106  {
107  std::swap(*alpha_mf[lev], *alpha_old_mf[lev]);
108 
109  const Set::Scalar *DX = this->geom[lev].CellSize();
110  // Evolve alpha
111  for (amrex::MFIter mfi(*alpha_mf[lev], true); mfi.isValid(); ++mfi)
112  {
113  const amrex::Box& bx = mfi.tilebox();
114  amrex::Array4<Set::Scalar> const& alpha_new = (*alpha_mf[lev]).array(mfi);
115  amrex::Array4<const Set::Scalar> const& alpha = (*alpha_old_mf[lev]).array(mfi);
116 
117  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
118  {
119  Set::Scalar driving_force = 0.0;
120  // Chemical potential
121  Set::Scalar alpha_2 = alpha(i, j, k) * alpha(i, j, k);
122  Set::Scalar alpha_3 = alpha_2 * alpha(i, j, k);
123  driving_force += (ch.chempot / ch.eps) * (2.0 * alpha(i, j, k) - 6.0 * alpha_2 + 4.0 * alpha_3);
124  // Gradient
125  Set::Scalar alpha_lap = Numeric::Laplacian(alpha, i, j, k, 0, DX);
126  driving_force -= ch.eps * ch.grad * alpha_lap;
127  // Update
128  alpha_new(i, j, k) =
129  alpha(i, j, k) - ch.L * dt * driving_force;
130  });
131  }
132  }
133 
134  // Tag cells for mesh refinement based on temperature gradient
135  void TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar /*time*/, int /*ngrow*/)
136  {
137  // Get cell dimensions as done above.
138  const Set::Scalar* DX = geom[lev].CellSize();
139  // Calculate the diagonal.
140  Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
141 
142  // Iterate over the patches on this level
143  for (amrex::MFIter mfi(*alpha_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
144  {
145  // Get the box and handles as done above.
146  const amrex::Box& bx = mfi.tilebox();
147  amrex::Array4<char> const& tags = a_tags.array(mfi);
148  amrex::Array4<Set::Scalar> const& alpha = (*alpha_mf[lev]).array(mfi);
149 
150  // Iterate over the grid as done above.
151  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
152  {
153  // Calculate the temperature gradient.
154  Set::Vector grad = Numeric::Gradient(alpha, i, j, k, 0, DX);
155 
156  // Is the gradient * cell_size too big? If so, then
157  // mark this cell as needing refinement.
158  if (grad.lpNorm<2>() * dr > refinement_threshold)
159  tags(i, j, k) = amrex::TagBox::SET;
160  });
161  }
162  }
163 
164 protected:
165  Set::Field<Set::Scalar> alpha_mf; // Temperature field variable (current timestep)
166  Set::Field<Set::Scalar> alpha_old_mf; // Temperature field variable (previous timestep)
167 
168  struct {
169  Set::Scalar L = NAN;
173  } ch;
174 
175 private:
176  int number_of_components = 1; // Number of components
177  int number_of_ghost_cells = 2; // Number of ghost cells
178 
179  Set::Scalar refinement_threshold = NAN; // Criterion for cell refinement
180 
181  IC::IC* ic; // Object used to initialize temperature field
182  BC::BC<Set::Scalar>* bc; // Object used to update temp field boundary ghost cells
183 };
184 } // namespace Integrator
185 #endif
Integrator::AllenCahn::ic
IC::IC * ic
Definition: AllenCahn.H:181
IC::BMP
Definition: BMP.H:20
Integrator::AllenCahn
Definition: AllenCahn.H:41
IC::IC::Initialize
void Initialize(const int &a_lev, Set::Field< Set::Scalar > &a_field, Set::Scalar a_time=0.0)
Definition: IC.H:26
BC::BC< Set::Scalar >
IC::IC
Pure abstract IC object from which all other IC objects inherit.
Definition: IC.H:12
Integrator::AllenCahn::number_of_ghost_cells
int number_of_ghost_cells
Definition: AllenCahn.H:177
Integrator::AllenCahn::AllenCahn
AllenCahn()
Definition: AllenCahn.H:45
Sphere.H
Integrator::AllenCahn::number_of_components
int number_of_components
Definition: AllenCahn.H:176
Integrator.H
Integrator::AllenCahn::grad
Set::Scalar grad
Definition: AllenCahn.H:171
Integrator::AllenCahn::bc
BC::BC< Set::Scalar > * bc
Definition: AllenCahn.H:182
Set::Field< Set::Scalar >
Definition: Set.H:236
Integrator::AllenCahn::Parse
static void Parse(AllenCahn &value, IO::ParmParse &pp)
Definition: AllenCahn.H:60
Integrator::AllenCahn::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int)
Definition: AllenCahn.H:135
Expression.H
Integrator::AllenCahn::eps
Set::Scalar eps
Definition: AllenCahn.H:170
Set::Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition: Base.H:20
Constant.H
ParmParse.H
Numeric::Gradient
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)
Definition: Stencil.H:655
Constant.H
Util::Random
Set::Scalar Random()
Definition: Set.cpp:9
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
Integrator::Integrator::RegisterNewFab
void RegisterNewFab(Set::Field< Set::Scalar > &new_fab, BC::BC< Set::Scalar > *new_bc, int ncomp, int nghost, std::string name, bool writeout)
Definition: Integrator.cpp:292
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:105
Integrator::AllenCahn::AllenCahn
AllenCahn(IO::ParmParse &pp)
Definition: AllenCahn.H:49
Numeric::Laplacian
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])
Definition: Stencil.H:530
Integrator::AllenCahn::refinement_threshold
Set::Scalar refinement_threshold
Definition: AllenCahn.H:179
Integrator::AllenCahn::L
Set::Scalar L
Definition: AllenCahn.H:169
pp_query_default
#define pp_query_default(...)
Definition: ParmParse.H:100
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:165
Integrator::AllenCahn::alpha_mf
Set::Field< Set::Scalar > alpha_mf
Definition: AllenCahn.H:165
Integrator
Collection of numerical integrator objects.
Definition: AllenCahn.H:39
IC::Constant
Definition: Constant.H:18
Random.H
BMP.H
BC::Constant
Definition: Constant.H:89
IC::Expression
Definition: Expression.H:44
Integrator::AllenCahn::ch
struct Integrator::AllenCahn::@0 ch
Integrator::AllenCahn::Advance
void Advance(int lev, Set::Scalar, Set::Scalar dt)
Definition: AllenCahn.H:105
pp_query_validate
#define pp_query_validate(...)
Definition: ParmParse.H:101
IO::ParmParse
Definition: ParmParse.H:110
Integrator::AllenCahn::chempot
Set::Scalar chempot
Definition: AllenCahn.H:172
Integrator::AllenCahn::Initialize
void Initialize(int lev)
Definition: AllenCahn.H:98
INFO
#define INFO
Definition: Util.H:20
Integrator::AllenCahn::alpha_old_mf
Set::Field< Set::Scalar > alpha_old_mf
Definition: AllenCahn.H:166
Integrator::Integrator::dt
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition: Integrator.H:303
IC.H
IC::Sphere
Definition: Sphere.H:13