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