Alamo
Dendrite.H
Go to the documentation of this file.
1 //
2 // This is a BASIC phase field method for dendritic growth as described in
3 // "A Numerical Approach to Three-Dimensional Dendritic Solidification" by Ryo Kobayashi (1994)
4 // Reference: https://doi.org/10.1080/10586458.1994.10504577
5 //
6 
7 #ifndef INTEGRATOR_DENDRITE_H // Include guards
8 #define INTEGRATOR_DENDRITE_H //
9 
10 #include <iostream>
11 #include <fstream>
12 #include <iomanip>
13 
14 #include "AMReX.H"
15 #include "AMReX_ParallelDescriptor.H"
16 #include "AMReX_ParmParse.H"
17 
18 #include "IO/ParmParse.H"
19 #include "Integrator/Integrator.H"
20 #include "BC/Constant.H"
21 #include "IC/IC.H"
22 #include "IC/Sphere.H"
23 #include "IC/Constant.H"
24 #include "IC/Expression.H"
25 #include "Numeric/Stencil.H"
26 
27 namespace Integrator
28 {
29 class Dendrite : virtual public Integrator
30 {
31 public:
32  // Empty constructor
34  {}
35 
36  // Constructor that triggers parse
38  {
39  Parse(*this, pp);
40  }
41 
42  static void Parse(Dendrite& value, IO::ParmParse& pp)
43  {
44  pp_query_required("alpha", value.alpha); // Pre-multiplier of "m" barrier height
45  pp_query_required("delta", value.delta); // Anisotropy factor
46  pp_query_required("gamma", value.gamma); // Anisotropic temperature coupling factor
47 
48  pp_query_required("eps", value.eps); // Diffuse boundary width
49  pp_query_required("tau", value.tau); // Diffusive timescale
50 
51  // Refinement criteria for temperature
52  pp_query_default("heat.refinement_threshold", value.refinement_threshold_temp, 0.01);
53  // Refinement criteria for phi
54  pp_query_default("phi.refinement_threshold", value.refinement_threshold_phi, 0.01);
55 
56  // Expressions for ICs
57  value.ic_temp = new IC::Expression(value.geom, pp, "ic.temp");
58  value.ic_phi = new IC::Expression(value.geom, pp, "ic.phi");
59 
60  // Use a constant BC object for temperature
61  value.bc_temp = new BC::Constant(1);
62  value.bc_phi = new BC::Constant(1);
63  pp_queryclass("bc.temp", *static_cast<BC::Constant*>(value.bc_temp));
64  pp_queryclass("bc.phi", *static_cast<BC::Constant*>(value.bc_phi));
65 
66  value.RegisterNewFab(value.temp_mf, value.bc_temp, value.number_of_components, value.number_of_ghost_cells, "Temp", true);
67  value.RegisterNewFab(value.temp_old_mf, value.bc_temp, value.number_of_components, value.number_of_ghost_cells, "Temp_old", false);
68  value.RegisterNewFab(value.phi_mf, value.bc_phi, value.number_of_components, value.number_of_ghost_cells, "phi", true);
69  value.RegisterNewFab(value.phi_old_mf, value.bc_phi, value.number_of_components, value.number_of_ghost_cells, "phi_old", false);
70  }
71 protected:
72 
73  // Use the ic object to initialize the temperature field
74  void Initialize(int lev)
75  {
78  ic_temp->Initialize(lev, temp_mf);
79  ic_phi->Initialize(lev, phi_mf);
80  }
81 
82 
83  // Integrate the heat equation
84  void Advance(int lev, Set::Scalar /*time*/, Set::Scalar dt)
85  {
86  // Swap the old temp fab and the new temp fab so we use
87  // the new one.
88  std::swap(*temp_mf[lev], *temp_old_mf[lev]);
89  std::swap(*phi_mf[lev], *phi_old_mf[lev]);
90 
91  // Get the cell size corresponding to this level
92  const Set::Scalar* DX = geom[lev].CellSize();
93 
94  // Iterate over all of the patches on this level
95  for (amrex::MFIter mfi(*temp_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
96  {
97  // Get the box (index dimensions) for this patch
98  const amrex::Box& bx = mfi.tilebox();
99 
100  // Get an array-accessible handle to the data on this patch.
101  amrex::Array4<const Set::Scalar> const& temp = (*temp_old_mf[lev]).array(mfi);
102  amrex::Array4<Set::Scalar> const& temp_new = (*temp_mf[lev]).array(mfi);
103  amrex::Array4<const Set::Scalar> const& phi = (*phi_old_mf[lev]).array(mfi);
104  amrex::Array4<Set::Scalar> const& phi_new = (*phi_mf[lev]).array(mfi);
105 
106  // Iterate over the grid on this patch
107  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
108  {
109 
110  Set::Vector gradphi = Numeric::Gradient(phi, i, j, k, 0, DX);
111 
112  //
113  // Calculate anisotropy function
114  // [ equation 2.12 in reference]
115  //
116  Set::Scalar v44 = AMREX_D_TERM(
117  gradphi(0) * gradphi(0) * gradphi(0) * gradphi(0),
118  +gradphi(1) * gradphi(1) * gradphi(1) * gradphi(1),
119  +gradphi(2) * gradphi(2) * gradphi(2) * gradphi(2));
120  Set::Scalar v22 = gradphi.squaredNorm(); v22 *= v22;
121  Set::Scalar sigma = 1.0 - delta * (1.0 - v44 / (v22 + 1E-12));
122 
123  //
124  // Calculate anisotropic barrier height
125  // [ unnumbered equation on page 67 in the reference ]
126  //
127  Set::Scalar m = -(alpha / Set::Constant::Pi) * std::atan(gamma * sigma * temp(i, j, k));
128 
129 
130  //
131  // Evolve order paramter
132  // [ equation 2.10 in the reference ]
133  //
134  Set::Scalar lapphi = Numeric::Laplacian(phi, i, j, k, 0, DX);
135  Set::Scalar Dphi = (eps * eps * lapphi + phi(i, j, k) * (1.0 - phi(i, j, k)) * (phi(i, j, k) - 0.5 + m)) / tau;
136  phi_new(i, j, k) = phi(i, j, k) + dt * Dphi;
137 
138  //
139  // Evolve temperature
140  // [ equation 2.11 in the reference ]
141  //
142  Set::Scalar laptemp = Numeric::Laplacian(temp, i, j, k, 0, DX);
143  temp_new(i, j, k) = temp(i, j, k) + dt * (laptemp + Dphi);
144  });
145  }
146  }
147 
148  // Tag cells for mesh refinement based on temperature gradient
149  void TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar /*time*/, int /*ngrow*/)
150  {
151  // Get cell dimensions as done above.
152  const Set::Scalar* DX = geom[lev].CellSize();
153  // Calculate the diagonal.
154  Set::Scalar dr = sqrt(AMREX_D_TERM(DX[0] * DX[0], +DX[1] * DX[1], +DX[2] * DX[2]));
155 
156  // Iterate over the patches on this level
157  for (amrex::MFIter mfi(*temp_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
158  {
159  // Get the box and handles as done above.
160  const amrex::Box& bx = mfi.tilebox();
161  amrex::Array4<char> const& tags = a_tags.array(mfi);
162  amrex::Array4<Set::Scalar> const& temp = (*temp_mf[lev]).array(mfi);
163  amrex::Array4<Set::Scalar> const& phi = (*phi_mf[lev]).array(mfi);
164 
165  // Iterate over the grid as done above.
166  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
167  {
168  // Calculate the temperature gradient.
169  Set::Vector grad_temp = Numeric::Gradient(temp, i, j, k, 0, DX);
170  Set::Vector grad_phi = Numeric::Gradient(phi, i, j, k, 0, DX);
171 
172  // Is the gradient * cell_size too big? If so, then
173  // mark this cell as needing refinement.
174  if (grad_temp.lpNorm<2>() * dr > refinement_threshold_temp)
175  tags(i, j, k) = amrex::TagBox::SET;
176  if (grad_phi.lpNorm<2>() * dr > refinement_threshold_phi)
177  tags(i, j, k) = amrex::TagBox::SET;
178  });
179  }
180  }
181 
182 protected:
187 
188 private:
189  int number_of_components = 1; // Number of components
190  int number_of_ghost_cells = 1; // Number of ghost cells
191 
192 
193  // calculation of m
195  // evolution of p
197 
198  Set::Scalar refinement_threshold_temp = NAN; // Criterion for cell refinement
199  Set::Scalar refinement_threshold_phi = NAN; // Criterion for cell refinement
200 
203 };
204 } // namespace Integrator
205 #endif
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::Dendrite::phi_old_mf
Set::Field< Set::Scalar > phi_old_mf
Definition: Dendrite.H:186
Integrator::Dendrite::Dendrite
Dendrite(IO::ParmParse &pp)
Definition: Dendrite.H:37
Sphere.H
Integrator.H
Set::Field< Set::Scalar >
Definition: Set.H:236
Integrator::Dendrite::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int)
Definition: Dendrite.H:149
Integrator::Dendrite::ic_temp
IC::IC * ic_temp
Definition: Dendrite.H:201
Expression.H
Integrator::Dendrite::number_of_ghost_cells
int number_of_ghost_cells
Definition: Dendrite.H:190
Set::Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition: Base.H:20
Constant.H
ParmParse.H
Integrator::Dendrite::gamma
Set::Scalar gamma
Definition: Dendrite.H:194
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
Integrator::Dendrite::Dendrite
Dendrite()
Definition: Dendrite.H:33
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::Dendrite::Parse
static void Parse(Dendrite &value, IO::ParmParse &pp)
Definition: Dendrite.H:42
Integrator::Dendrite::eps
Set::Scalar eps
Definition: Dendrite.H:196
Integrator::Dendrite::refinement_threshold_phi
Set::Scalar refinement_threshold_phi
Definition: Dendrite.H:199
Integrator::Dendrite::Advance
void Advance(int lev, Set::Scalar, Set::Scalar dt)
Definition: Dendrite.H:84
Integrator::Dendrite::Initialize
void Initialize(int lev)
Definition: Dendrite.H:74
Integrator::Dendrite
Definition: Dendrite.H:29
Integrator::Dendrite::number_of_components
int number_of_components
Definition: Dendrite.H:189
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::Dendrite::alpha
Set::Scalar alpha
Definition: Dendrite.H:194
pp_query_default
#define pp_query_default(...)
Definition: ParmParse.H:100
Integrator
Collection of numerical integrator objects.
Definition: AllenCahn.H:39
BC::Constant
Definition: Constant.H:89
Integrator::Dendrite::bc_temp
BC::BC< Set::Scalar > * bc_temp
Definition: Dendrite.H:202
Integrator::Dendrite::tau
Set::Scalar tau
Definition: Dendrite.H:196
IC::Expression
Definition: Expression.H:44
Integrator::Dendrite::bc_phi
BC::BC< Set::Scalar > * bc_phi
Definition: Dendrite.H:202
Integrator::Dendrite::delta
Set::Scalar delta
Definition: Dendrite.H:194
Set::Constant::Pi
static const Set::Scalar Pi
Definition: Set.H:286
pp_query_required
#define pp_query_required(...)
Definition: ParmParse.H:99
IO::ParmParse
Definition: ParmParse.H:110
Integrator::Dendrite::refinement_threshold_temp
Set::Scalar refinement_threshold_temp
Definition: Dendrite.H:198
Integrator::Dendrite::phi_mf
Set::Field< Set::Scalar > phi_mf
Definition: Dendrite.H:185
Integrator::Dendrite::ic_phi
IC::IC * ic_phi
Definition: Dendrite.H:201
Integrator::Dendrite::temp_old_mf
Set::Field< Set::Scalar > temp_old_mf
Definition: Dendrite.H:184
Integrator::Integrator::dt
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition: Integrator.H:303
IC.H
Integrator::Dendrite::temp_mf
Set::Field< Set::Scalar > temp_mf
Definition: Dendrite.H:183