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