Alamo
TopOp.H
Go to the documentation of this file.
1 #ifndef INTEGRATOR_TOPOP_H
2 #define INTEGRATOR_TOPOP_H
3 #include <iostream>
4 #include <fstream>
5 #include <iomanip>
6 #include <numeric>
7 
8 #include "AMReX.H"
9 #include "AMReX_ParallelDescriptor.H"
10 #include "AMReX_ParmParse.H"
11 
12 #include "IO/ParmParse.H"
14 
15 
16 #include "IC/IC.H"
17 #include "BC/BC.H"
21 
22 #include "IC/Ellipse.H"
23 #include "IC/Voronoi.H"
24 #include "IC/Constant.H"
25 #include "IC/BMP.H"
26 #include "BC/Constant.H"
27 #include "Numeric/Stencil.H"
28 
29 #include "Model/Solid/Solid.H"
30 #include "Solver/Nonlocal/Linear.H"
31 #include "Solver/Nonlocal/Newton.H"
32 
33 #include "Operator/Operator.H"
34 
35 
36 namespace Integrator
37 {
38 template<class MODEL>
39 class TopOp: virtual public Base::Mechanics<MODEL>
40 {
41 public:
42 
43  TopOp(): Base::Mechanics<MODEL>() {}
44  TopOp(IO::ParmParse& pp): Base::Mechanics<MODEL>()
45  {
46  Parse(*this, pp);
47  }
48 
49  // The mechanics integrator manages the solution of an elastic
50  // solve using the MLMG solver.
51  static void Parse(TopOp& value, IO::ParmParse& pp)
52  {
54 
55  pp_queryclass("model", value.model);
56  // Read in IC for psi
57  if (pp.contains("psi.ic.type"))
58  {
59  // Initial condition for psi field
60  pp.select<IC::Ellipse,IC::Constant>("psi.ic",value.ic_psi,value.geom);
61 
62  value.bc_psi = new BC::Constant(1, pp, "psi.bc");
63  value.RegisterNewFab(value.psi_mf, value.bc_psi, 1, 2, "psi", true);
64  value.RegisterNewFab(value.psi_old_mf, value.bc_psi, 1, 2, "psiold", false);
65  value.psi_on = true;
66  }
67  pp_query_default("eta_ref_threshold", value.m_eta_ref_threshold, 0.01); // Refinement threshold based on eta
68  pp_query_default("alpha", value.alpha, 1.0); // :math:`\alpha` parameter
69  pp_query_default("beta", value.beta, 1.0); // :math:`\beta` parameter
70  pp_query_default("gamma", value.gamma, 1.0); // :math:`\gamma` parameter
71  pp_queryclass("L", value.L); // Mobility
72  if (pp.contains("volume0frac"))
73  {
74  if (pp.contains("volume0")) Util::Abort(INFO, "Cannot specify volume0frac and volume0");
75  Set::Scalar volumefrac;
76  pp_query_default("volume0frac", volumefrac, 0.5); // Prescribed volume fraction
77  value.volume0 = volumefrac *
78  AMREX_D_TERM((value.geom[0].ProbHi()[0] - value.geom[0].ProbLo()[0]),
79  *(value.geom[0].ProbHi()[1] - value.geom[0].ProbLo()[1]),
80  *(value.geom[0].ProbHi()[2] - value.geom[0].ProbLo()[2]));
81  }
82  else
83  pp_query_default("volume0", value.volume0, 0.5); // Prescribed total vlume
84 
85  pp_queryclass("lambda", value.lambda); // Lagrange multiplier (can be interplated)
86 
87 
88  value.RegisterIntegratedVariable(&value.volume, "volume");
89  value.RegisterIntegratedVariable(&value.w_chem_potential, "chem_potential");
90  value.RegisterIntegratedVariable(&value.w_bndry, "bndry");
91  value.RegisterIntegratedVariable(&value.w_elastic, "elastic");
92  }
93 
94  void Initialize(int lev) override
95  {
97  if (psi_on) ic_psi->Initialize(lev, psi_mf);
98  if (psi_on) ic_psi->Initialize(lev, psi_old_mf);
99  }
100 
101  virtual void UpdateModel(int a_step, Set::Scalar /*a_time*/) override
102  {
104 
105  if (a_step > 0) return;
106 
107  for (int lev = 0; lev <= finest_level; ++lev)
108  {
109  model_mf[lev]->setVal(model);
110  Util::RealFillBoundary(*model_mf[lev], geom[lev]);
111  Util::RealFillBoundary(*psi_mf[lev], geom[lev]);
112  Util::RealFillBoundary(*psi_old_mf[lev], geom[lev]);
113  }
114 
115  }
116 
117 
118  void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
119  {
120  BL_PROFILE("TopOp::Advance");
122  std::swap(psi_old_mf[lev], psi_mf[lev]);
123  const Set::Scalar* DX = geom[lev].CellSize();
124  amrex::Box domain = geom[lev].Domain();
125 
126  Set::Scalar Lnow = L(time);
127  Set::Scalar lambdaT = lambda(time);
128 
129  for (amrex::MFIter mfi(*psi_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
130  {
131  amrex::Box bx = mfi.tilebox();
132  bx.grow(1);
133  bx = bx & domain;
134  amrex::Array4<const Set::Matrix> const& sig = (*stress_mf[lev]).array(mfi);
135  amrex::Array4<const Set::Matrix> const& eps = (*strain_mf[lev]).array(mfi);
136  amrex::Array4<const Set::Scalar> const& psi = (*psi_old_mf[lev]).array(mfi);
137  amrex::Array4<Set::Scalar> const& psinew = (*psi_mf[lev]).array(mfi);
138 
139  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
140  {
141  Set::Scalar driving_force = 0.0;
142 
143  driving_force += alpha * 2.0 * psi(i, j, k) * (2.0 * psi(i, j, k) * psi(i, j, k) - 3.0 * psi(i, j, k) + 1.0);
144  driving_force += -beta * Numeric::Laplacian(psi, i, j, k, 0, DX);
145 
146  Set::Matrix sig_avg = Numeric::Interpolate::NodeToCellAverage(sig, i, j, k, 0);
147  Set::Matrix eps_avg = Numeric::Interpolate::NodeToCellAverage(eps, i, j, k, 0);
148 
149  driving_force += -gamma * 0.5 * (sig_avg.transpose() * eps_avg).trace() * psi(i, j, k);
150 
151  driving_force += lambdaT * (volume - volume0);
152 
153  psinew(i, j, k) = psi(i, j, k) - Lnow * dt * driving_force;
154  if (psinew(i, j, k) < 0.0) psinew(i, j, k) = 0.0;
155  if (psinew(i, j, k) > 1.0) psinew(i, j, k) = 1.0;
156  });
157  }
158  }
159 
160  void Integrate(int amrlev, Set::Scalar time, int step,
161  const amrex::MFIter& mfi, const amrex::Box& box) override
162  {
163  BL_PROFILE("TopOp::Integrate");
164  Base::Mechanics<MODEL>::Integrate(amrlev, time, step, mfi, box);
165 
166  const amrex::Real* DX = geom[amrlev].CellSize();
167  Set::Scalar dv = AMREX_D_TERM(DX[0], *DX[1], *DX[2]);
168 
169  amrex::Array4<amrex::Real> const& psi = (*psi_mf[amrlev]).array(mfi);
170  amrex::Array4<const Set::Matrix> const& sig = (*stress_mf[amrlev]).array(mfi);
171  amrex::Array4<const Set::Matrix> const& eps = (*strain_mf[amrlev]).array(mfi);
172  amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
173  {
174  volume += psi(i, j, k, 0) * dv;
175  w_chem_potential += alpha * psi(i, j, k, 0) * psi(i, j, k, 0) * (1. - psi(i, j, k, 0) * psi(i, j, k, 0)) * dv;
176  w_bndry += beta * 0.5 * Numeric::Gradient(psi, i, j, k, 0, DX).squaredNorm() * dv;
177 
178  Set::Matrix sig_avg = Numeric::Interpolate::NodeToCellAverage(sig, i, j, k, 0);
179  Set::Matrix eps_avg = Numeric::Interpolate::NodeToCellAverage(eps, i, j, k, 0);
180  w_elastic += gamma * 0.5 * (sig_avg.transpose() * eps_avg).trace() * psi(i, j, k) * dv;
181  });
182  }
183 
184 
185  void TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar a_time, int a_ngrow) override
186  {
188  Base::Mechanics<MODEL>::TagCellsForRefinement(lev, a_tags, a_time, a_ngrow);
189 
190  Set::Vector DX(geom[lev].CellSize());
191  Set::Scalar DXnorm = DX.lpNorm<2>();
192  for (amrex::MFIter mfi(*model_mf[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
193  {
194  amrex::Box bx = mfi.tilebox();
195  amrex::Array4<char> const& tags = a_tags.array(mfi);
196  if (psi_on)
197  {
198  amrex::Array4<Set::Scalar> const& psi = psi_mf[lev]->array(mfi);
199  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
200  {
201  auto sten = Numeric::GetStencil(i, j, k, bx);
202  {
203  Set::Vector gradpsi = Numeric::Gradient(psi, i, j, k, 0, DX.data(), sten);
204  if (gradpsi.lpNorm<2>() * DXnorm > m_eta_ref_threshold)
205  tags(i, j, k) = amrex::TagBox::SET;
206  }
207  });
208  }
209  }
210  }
211 
212 
213 protected:
214  MODEL model;
215  IC::IC* ic_psi = nullptr;
219 
220 
229 
230 
234  //Set::Scalar L = 1.0;
237  //Set::Scalar lambda = 1.0;
239 
244 };
245 
246 
247 
248 
249 
250 
251 
252 
253 
254 
255 } // namespace Integrator
256 #endif
TensionTest.H
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
Util::RealFillBoundary
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T >> &a_mf, const amrex::Geometry &)
Definition: Util.H:319
BC::BC< Set::Scalar >
IC::IC
Pure abstract IC object from which all other IC objects inherit.
Definition: IC.H:21
Integrator::Mechanics
Definition: Mechanics.H:101
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
IC::Ellipse
Definition: Ellipse.H:15
Integrator::Base::Mechanics
Definition: Mechanics.H:23
Integrator::Integrator::box
std::vector< amrex::Box > box
Definition: Integrator.H:463
Integrator::TopOp::w_chem_potential
Set::Scalar w_chem_potential
Definition: TopOp.H:241
Integrator::TopOp::lambda
Numeric::Interpolator::Linear< Set::Scalar > lambda
Definition: TopOp.H:238
Integrator::TopOp::model
MODEL model
Definition: TopOp.H:214
Integrator::TopOp::TopOp
TopOp(IO::ParmParse &pp)
Definition: TopOp.H:44
Expression.H
Integrator::TopOp::UpdateModel
virtual void UpdateModel(int a_step, Set::Scalar) override
Definition: TopOp.H:101
Integrator::TopOp::L
Numeric::Interpolator::Linear< Set::Scalar > L
Definition: TopOp.H:235
Integrator::Base::Mechanics::psi_mf
Set::Field< Set::Scalar > psi_mf
Definition: Mechanics.H:455
Set::Field< Set::Scalar >
Definition: Set.H:236
Integrator::TopOp::Initialize
void Initialize(int lev) override
Definition: TopOp.H:94
Integrator::Base::Mechanics::model_mf
Set::Field< MODEL > model_mf
Definition: Mechanics.H:453
Integrator::TopOp::Integrate
void Integrate(int amrlev, Set::Scalar time, int step, const amrex::MFIter &mfi, const amrex::Box &box) override
Definition: TopOp.H:160
Integrator::Base::Mechanics::m_type
Type m_type
Definition: Mechanics.H:459
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
Integrator::Base::Mechanics::Initialize
void Initialize(int lev) override
Use the #ic object to initialize::Temp.
Definition: Mechanics.H:134
Integrator::TopOp::beta
Set::Scalar beta
Definition: TopOp.H:232
Numeric::Interpolate::NodeToCellAverage
static AMREX_FORCE_INLINE T NodeToCellAverage(const amrex::Array4< const T > &f, const int &i, const int &j, const int &k, const int &m)
Definition: Stencil.H:1315
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:106
Solid.H
Set::Matrix
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition: Base.H:23
Mechanics.H
Constant.H
Integrator::TopOp::alpha
Set::Scalar alpha
Definition: TopOp.H:231
Integrator::TopOp::m_eta_ref_threshold
Set::Scalar m_eta_ref_threshold
Definition: TopOp.H:217
Integrator::TopOp::w_bndry
Set::Scalar w_bndry
Definition: TopOp.H:242
Integrator::Base::Mechanics::Parse
static void Parse(Mechanics &value, IO::ParmParse &pp)
Definition: Mechanics.H:34
Integrator::TopOp::Advance
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition: TopOp.H:118
IO::ParmParse::contains
bool contains(std::string name)
Definition: ParmParse.H:153
Integrator::TopOp::Parse
static void Parse(TopOp &value, IO::ParmParse &pp)
Definition: TopOp.H:51
pp_query_default
#define pp_query_default(...)
Definition: ParmParse.H:100
Integrator::TopOp::ic_psi
IC::IC * ic_psi
Definition: TopOp.H:215
Numeric::Interpolator::Linear< Set::Scalar >
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:166
Integrator::TopOp::gamma
Set::Scalar gamma
Definition: TopOp.H:233
Integrator
Collection of numerical integrator objects.
Definition: AllenCahn.H:40
Integrator::Base::Mechanics::Integrate
void Integrate(int amrlev, Set::Scalar, int, const amrex::MFIter &mfi, const amrex::Box &a_box) override
Definition: Mechanics.H:371
IC::Constant
Definition: Constant.H:18
Integrator::TopOp::psi_old_mf
Set::Field< Set::Scalar > psi_old_mf
Definition: TopOp.H:218
Newton.H
BMP.H
BC::Constant
Definition: Constant.H:43
Integrator::TopOp::w_elastic
Set::Scalar w_elastic
Definition: TopOp.H:243
Voronoi.H
Integrator::Base::Mechanics::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int) override
Definition: Mechanics.H:429
BC.H
IO::ParmParse
Definition: ParmParse.H:112
Integrator::TopOp::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar a_time, int a_ngrow) override
Definition: TopOp.H:185
Integrator::Base::Mechanics::Advance
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition: Mechanics.H:232
Integrator::Base::Mechanics::strain_mf
Set::Field< Set::Matrix > strain_mf
Definition: Mechanics.H:465
Integrator::Base::Mechanics::stress_mf
Set::Field< Set::Matrix > stress_mf
Definition: Mechanics.H:464
Integrator::TopOp::bc_psi
BC::BC< Set::Scalar > * bc_psi
Definition: TopOp.H:216
Integrator::TopOp::volume0
Set::Scalar volume0
Definition: TopOp.H:236
Integrator::TopOp
Definition: TopOp.H:39
INFO
#define INFO
Definition: Util.H:20
Integrator::TopOp::TopOp
TopOp()
Definition: TopOp.H:43
Integrator::Integrator::dt
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition: Integrator.H:398
IC.H
Integrator::Integrator::RegisterIntegratedVariable
void RegisterIntegratedVariable(Set::Scalar *integrated_variable, std::string name, bool extensive=true)
Definition: Integrator.cpp:342
Numeric::GetStencil
static AMREX_FORCE_INLINE std::array< StencilType, AMREX_SPACEDIM > GetStencil(const int i, const int j, const int k, const amrex::Box domain)
Definition: Stencil.H:20
Integrator::Base::Mechanics::psi_on
bool psi_on
Definition: Mechanics.H:456
Linear.H
Integrator::TopOp::volume
Set::Scalar volume
Definition: TopOp.H:240
IO::ParmParse::select
void select(std::string name, PTRTYPE *&ic_eta, Args &&... args)
Definition: ParmParse.H:541
Ellipse.H