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  std::string type;
60  pp_query_validate("psi.ic.type", type, {"ellipse", "constant"}); // Read IC type for the eta field
61  if (type == "ellipse") value.ic_psi = new IC::Ellipse(value.geom, pp, "psi.ic.ellipse");
62  else if (type == "constant") value.ic_psi = new IC::Constant(value.geom, pp, "psi.ic.constant");
63  else Util::Abort(INFO, "Invalid value for psi.ic.type: ", type);
64 
65  value.bc_psi = new BC::Constant(1, pp, "psi.bc");
66  value.RegisterNewFab(value.psi_mf, value.bc_psi, 1, 2, "psi", true);
67  value.RegisterNewFab(value.psi_old_mf, value.bc_psi, 1, 2, "psiold", false);
68  value.psi_on = true;
69  }
70  pp_query_default("eta_ref_threshold", value.m_eta_ref_threshold, 0.01); // Refinement threshold based on eta
71  pp_query_default("alpha", value.alpha, 1.0); // :math:`\alpha` parameter
72  pp_query_default("beta", value.beta, 1.0); // :math:`\beta` parameter
73  pp_query_default("gamma", value.gamma, 1.0); // :math:`\gamma` parameter
74  pp_queryclass("L", value.L); // Mobility
75  if (pp.contains("volume0frac"))
76  {
77  if (pp.contains("volume0")) Util::Abort(INFO, "Cannot specify volume0frac and volume0");
78  Set::Scalar volumefrac;
79  pp_query_default("volume0frac", volumefrac, 0.5); // Prescribed volume fraction
80  value.volume0 = volumefrac *
81  AMREX_D_TERM((value.geom[0].ProbHi()[0] - value.geom[0].ProbLo()[0]),
82  *(value.geom[0].ProbHi()[1] - value.geom[0].ProbLo()[1]),
83  *(value.geom[0].ProbHi()[2] - value.geom[0].ProbLo()[2]));
84  }
85  else
86  pp_query_default("volume0", value.volume0, 0.5); // Prescribed total vlume
87 
88  pp_queryclass("lambda", value.lambda); // Lagrange multiplier (can be interplated)
89 
90 
91  value.RegisterIntegratedVariable(&value.volume, "volume");
92  value.RegisterIntegratedVariable(&value.w_chem_potential, "chem_potential");
93  value.RegisterIntegratedVariable(&value.w_bndry, "bndry");
94  value.RegisterIntegratedVariable(&value.w_elastic, "elastic");
95  }
96 
97  void Initialize(int lev) override
98  {
100  if (psi_on) ic_psi->Initialize(lev, psi_mf);
101  if (psi_on) ic_psi->Initialize(lev, psi_old_mf);
102  }
103 
104  virtual void UpdateModel(int a_step, Set::Scalar /*a_time*/) override
105  {
107 
108  if (a_step > 0) return;
109 
110  for (int lev = 0; lev <= finest_level; ++lev)
111  {
112  model_mf[lev]->setVal(model);
113  Util::RealFillBoundary(*model_mf[lev], geom[lev]);
114  Util::RealFillBoundary(*psi_mf[lev], geom[lev]);
115  Util::RealFillBoundary(*psi_old_mf[lev], geom[lev]);
116  }
117 
118  }
119 
120 
121  void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
122  {
123  BL_PROFILE("TopOp::Advance");
125  std::swap(psi_old_mf[lev], psi_mf[lev]);
126  const Set::Scalar* DX = geom[lev].CellSize();
127  amrex::Box domain = geom[lev].Domain();
128 
129  Set::Scalar Lnow = L(time);
130  Set::Scalar lambdaT = lambda(time);
131 
132  for (amrex::MFIter mfi(*psi_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
133  {
134  amrex::Box bx = mfi.tilebox();
135  bx.grow(1);
136  bx = bx & domain;
137  amrex::Array4<const Set::Matrix> const& sig = (*stress_mf[lev]).array(mfi);
138  amrex::Array4<const Set::Matrix> const& eps = (*strain_mf[lev]).array(mfi);
139  amrex::Array4<const Set::Scalar> const& psi = (*psi_old_mf[lev]).array(mfi);
140  amrex::Array4<Set::Scalar> const& psinew = (*psi_mf[lev]).array(mfi);
141 
142  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
143  {
144  Set::Scalar driving_force = 0.0;
145 
146  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);
147  driving_force += -beta * Numeric::Laplacian(psi, i, j, k, 0, DX);
148 
149  Set::Matrix sig_avg = Numeric::Interpolate::NodeToCellAverage(sig, i, j, k, 0);
150  Set::Matrix eps_avg = Numeric::Interpolate::NodeToCellAverage(eps, i, j, k, 0);
151 
152  driving_force += -gamma * 0.5 * (sig_avg.transpose() * eps_avg).trace() * psi(i, j, k);
153 
154  driving_force += lambdaT * (volume - volume0);
155 
156  psinew(i, j, k) = psi(i, j, k) - Lnow * dt * driving_force;
157  if (psinew(i, j, k) < 0.0) psinew(i, j, k) = 0.0;
158  if (psinew(i, j, k) > 1.0) psinew(i, j, k) = 1.0;
159  });
160  }
161  }
162 
163  void Integrate(int amrlev, Set::Scalar time, int step,
164  const amrex::MFIter& mfi, const amrex::Box& box) override
165  {
166  BL_PROFILE("TopOp::Integrate");
167  Base::Mechanics<MODEL>::Integrate(amrlev, time, step, mfi, box);
168 
169  const amrex::Real* DX = geom[amrlev].CellSize();
170  Set::Scalar dv = AMREX_D_TERM(DX[0], *DX[1], *DX[2]);
171 
172  amrex::Array4<amrex::Real> const& psi = (*psi_mf[amrlev]).array(mfi);
173  amrex::Array4<const Set::Matrix> const& sig = (*stress_mf[amrlev]).array(mfi);
174  amrex::Array4<const Set::Matrix> const& eps = (*strain_mf[amrlev]).array(mfi);
175  amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
176  {
177  volume += psi(i, j, k, 0) * dv;
178  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;
179  w_bndry += beta * 0.5 * Numeric::Gradient(psi, i, j, k, 0, DX).squaredNorm() * dv;
180 
181  Set::Matrix sig_avg = Numeric::Interpolate::NodeToCellAverage(sig, i, j, k, 0);
182  Set::Matrix eps_avg = Numeric::Interpolate::NodeToCellAverage(eps, i, j, k, 0);
183  w_elastic += gamma * 0.5 * (sig_avg.transpose() * eps_avg).trace() * psi(i, j, k) * dv;
184  });
185  }
186 
187 
188  void TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar a_time, int a_ngrow) override
189  {
191  Base::Mechanics<MODEL>::TagCellsForRefinement(lev, a_tags, a_time, a_ngrow);
192 
193  Set::Vector DX(geom[lev].CellSize());
194  Set::Scalar DXnorm = DX.lpNorm<2>();
195  for (amrex::MFIter mfi(*model_mf[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
196  {
197  amrex::Box bx = mfi.tilebox();
198  amrex::Array4<char> const& tags = a_tags.array(mfi);
199  if (psi_on)
200  {
201  amrex::Array4<Set::Scalar> const& psi = psi_mf[lev]->array(mfi);
202  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
203  {
204  auto sten = Numeric::GetStencil(i, j, k, bx);
205  {
206  Set::Vector gradpsi = Numeric::Gradient(psi, i, j, k, 0, DX.data(), sten);
207  if (gradpsi.lpNorm<2>() * DXnorm > m_eta_ref_threshold)
208  tags(i, j, k) = amrex::TagBox::SET;
209  }
210  });
211  }
212  }
213  }
214 
215 
216 protected:
217  MODEL model;
218  IC::IC* ic_psi = nullptr;
222 
223 
232 
233 
237  //Set::Scalar L = 1.0;
240  //Set::Scalar lambda = 1.0;
242 
247 };
248 
249 
250 
251 
252 
253 
254 
255 
256 
257 
258 } // namespace Integrator
259 #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:26
Util::RealFillBoundary
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T >> &a_mf, const amrex::Geometry &)
Definition: Util.H:307
BC::BC< Set::Scalar >
IC::IC
Pure abstract IC object from which all other IC objects inherit.
Definition: IC.H:12
Integrator::Mechanics
Definition: Mechanics.H:101
IC::Ellipse
Definition: Ellipse.H:21
Integrator::Base::Mechanics
Definition: Mechanics.H:23
Integrator::Integrator::box
std::vector< amrex::Box > box
Definition: Integrator.H:367
Integrator::TopOp::w_chem_potential
Set::Scalar w_chem_potential
Definition: TopOp.H:244
Integrator::TopOp::lambda
Numeric::Interpolator::Linear< Set::Scalar > lambda
Definition: TopOp.H:241
Integrator::TopOp::model
MODEL model
Definition: TopOp.H:217
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:104
Integrator::TopOp::L
Numeric::Interpolator::Linear< Set::Scalar > L
Definition: TopOp.H:238
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:97
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:163
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:655
Constant.H
Integrator::Base::Mechanics::Initialize
void Initialize(int lev) override
Use the #ic object to initialize::Temp.
Definition: Mechanics.H:140
Integrator::TopOp::beta
Set::Scalar beta
Definition: TopOp.H:235
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:1294
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
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:234
Integrator::TopOp::m_eta_ref_threshold
Set::Scalar m_eta_ref_threshold
Definition: TopOp.H:220
Integrator::TopOp::w_bndry
Set::Scalar w_bndry
Definition: TopOp.H:245
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::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:121
IO::ParmParse::contains
bool contains(std::string name)
Definition: ParmParse.H:151
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:218
Numeric::Interpolator::Linear< Set::Scalar >
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:165
Integrator::TopOp::gamma
Set::Scalar gamma
Definition: TopOp.H:236
Integrator
Collection of numerical integrator objects.
Definition: AllenCahn.H:39
Integrator::Base::Mechanics::Integrate
void Integrate(int amrlev, Set::Scalar, int, const amrex::MFIter &mfi, const amrex::Box &a_box) override
Definition: Mechanics.H:377
IC::Constant
Definition: Constant.H:18
Integrator::TopOp::psi_old_mf
Set::Field< Set::Scalar > psi_old_mf
Definition: TopOp.H:221
Newton.H
BMP.H
BC::Constant
Definition: Constant.H:89
Integrator::TopOp::w_elastic
Set::Scalar w_elastic
Definition: TopOp.H:246
Voronoi.H
Integrator::Base::Mechanics::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int) override
Definition: Mechanics.H:429
pp_query_validate
#define pp_query_validate(...)
Definition: ParmParse.H:101
BC.H
IO::ParmParse
Definition: ParmParse.H:110
Integrator::TopOp::TagCellsForRefinement
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar a_time, int a_ngrow) override
Definition: TopOp.H:188
Integrator::Base::Mechanics::Advance
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition: Mechanics.H:238
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:219
Integrator::TopOp::volume0
Set::Scalar volume0
Definition: TopOp.H:239
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:303
IC.H
Integrator::Integrator::RegisterIntegratedVariable
void RegisterIntegratedVariable(Set::Scalar *integrated_variable, std::string name, bool extensive=true)
Definition: Integrator.cpp:320
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:243
Ellipse.H