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"
32
33#include "Operator/Operator.H"
34
35
36namespace Integrator
37{
38template<class MODEL>
39class TopOp: virtual public Base::Mechanics<MODEL>
40{
41public:
42
43 TopOp(): Base::Mechanics<MODEL>() {}
44 TopOp(IO::ParmParse& pp): Base::Mechanics<MODEL>()
45 {
46 Parse(*this, pp);
47 }
48
50 {
51 delete ic_psi;
52 delete bc_psi;
53 }
54
55 // The mechanics integrator manages the solution of an elastic
56 // solve using the MLMG solver.
57 static void Parse(TopOp& value, IO::ParmParse& pp)
58 {
60
61 pp.queryclass<MODEL>("model", value.model);
62 // Read in IC for psi
63 if (pp.contains("psi.ic.type"))
64 {
65 // Initial condition for psi field
66 pp.select<IC::Ellipse,IC::Constant>("psi.ic",value.ic_psi,value.geom);
67
68 value.bc_psi = new BC::Constant(1, pp, "psi.bc");
69 value.RegisterNewFab(value.psi_mf, value.bc_psi, 1, 2, "psi", true);
70 value.RegisterNewFab(value.psi_old_mf, value.bc_psi, 1, 2, "psiold", false);
71 value.psi_on = true;
72 }
73 pp_query_default("eta_ref_threshold", value.m_eta_ref_threshold, 0.01); // Refinement threshold based on eta
74 pp_query_default("alpha", value.alpha, 1.0); // :math:`\alpha` parameter
75 pp_query_default("beta", value.beta, 1.0); // :math:`\beta` parameter
76 pp_query_default("gamma", value.gamma, 1.0); // :math:`\gamma` parameter
77 pp_queryclass("L", value.L); // Mobility
78 if (pp.contains("volume0frac"))
79 {
80 if (pp.contains("volume0")) Util::Abort(INFO, "Cannot specify volume0frac and volume0");
81 Set::Scalar volumefrac;
82 pp_query_default("volume0frac", volumefrac, 0.5); // Prescribed volume fraction
83 value.volume0 = volumefrac *
84 AMREX_D_TERM((value.geom[0].ProbHi()[0] - value.geom[0].ProbLo()[0]),
85 *(value.geom[0].ProbHi()[1] - value.geom[0].ProbLo()[1]),
86 *(value.geom[0].ProbHi()[2] - value.geom[0].ProbLo()[2]));
87 }
88 else
89 pp_query_default("volume0", value.volume0, 0.5); // Prescribed total vlume
90
91 pp_queryclass("lambda", value.lambda); // Lagrange multiplier (can be interplated)
92
93
94 value.RegisterIntegratedVariable(&value.volume, "volume");
95 value.RegisterIntegratedVariable(&value.w_chem_potential, "chem_potential");
96 value.RegisterIntegratedVariable(&value.w_bndry, "bndry");
97 value.RegisterIntegratedVariable(&value.w_elastic, "elastic");
98 }
99
100 void Initialize(int lev) override
101 {
103 if (psi_on) ic_psi->Initialize(lev, psi_mf);
105 }
106
107 virtual void UpdateModel(int a_step, Set::Scalar /*a_time*/) override
108 {
110
111 if (a_step > 0) return;
112
113 for (int lev = 0; lev <= finest_level; ++lev)
114 {
115 model_mf[lev]->setVal(model);
116 Util::RealFillBoundary(*model_mf[lev], geom[lev]);
117 Util::RealFillBoundary(*psi_mf[lev], geom[lev]);
118 Util::RealFillBoundary(*psi_old_mf[lev], geom[lev]);
119 }
120
121 }
122
123
124 void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
125 {
126 BL_PROFILE("TopOp::Advance");
128 std::swap(psi_old_mf[lev], psi_mf[lev]);
129 const Set::Scalar* DX = geom[lev].CellSize();
130 amrex::Box domain = geom[lev].Domain();
131
132 Set::Scalar Lnow = L(time);
133 Set::Scalar lambdaT = lambda(time);
134
135 for (amrex::MFIter mfi(*psi_mf[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
136 {
137 amrex::Box bx = mfi.tilebox();
138 bx.grow(1);
139 bx = bx & domain;
140 amrex::Array4<const Set::Matrix> const& sig = (*stress_mf[lev]).array(mfi);
141 amrex::Array4<const Set::Matrix> const& eps = (*strain_mf[lev]).array(mfi);
142 amrex::Array4<const Set::Scalar> const& psi = (*psi_old_mf[lev]).array(mfi);
143 amrex::Array4<Set::Scalar> const& psinew = (*psi_mf[lev]).array(mfi);
144
145 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
146 {
147 Set::Scalar driving_force = 0.0;
148
149 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);
150 driving_force += -beta * Numeric::Laplacian(psi, i, j, k, 0, DX);
151
152 Set::Matrix sig_avg = Numeric::Interpolate::NodeToCellAverage(sig, i, j, k, 0);
153 Set::Matrix eps_avg = Numeric::Interpolate::NodeToCellAverage(eps, i, j, k, 0);
154
155 driving_force += -gamma * 0.5 * (sig_avg.transpose() * eps_avg).trace() * psi(i, j, k);
156
157 driving_force += lambdaT * (volume - volume0);
158
159 psinew(i, j, k) = psi(i, j, k) - Lnow * dt * driving_force;
160 if (psinew(i, j, k) < 0.0) psinew(i, j, k) = 0.0;
161 if (psinew(i, j, k) > 1.0) psinew(i, j, k) = 1.0;
162 });
163 }
164 }
165
166 void Integrate(int amrlev, Set::Scalar time, int step,
167 const amrex::MFIter& mfi, const amrex::Box& box) override
168 {
169 BL_PROFILE("TopOp::Integrate");
170 Base::Mechanics<MODEL>::Integrate(amrlev, time, step, mfi, box);
171
172 const amrex::Real* DX = geom[amrlev].CellSize();
173 Set::Scalar dv = AMREX_D_TERM(DX[0], *DX[1], *DX[2]);
174
175 amrex::Array4<amrex::Real> const& psi = (*psi_mf[amrlev]).array(mfi);
176 amrex::Array4<const Set::Matrix> const& sig = (*stress_mf[amrlev]).array(mfi);
177 amrex::Array4<const Set::Matrix> const& eps = (*strain_mf[amrlev]).array(mfi);
178 amrex::ParallelFor(box, [=] AMREX_GPU_DEVICE(int i, int j, int k)
179 {
180 volume += psi(i, j, k, 0) * dv;
181 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;
182 w_bndry += beta * 0.5 * Numeric::Gradient(psi, i, j, k, 0, DX).squaredNorm() * dv;
183
184 Set::Matrix sig_avg = Numeric::Interpolate::NodeToCellAverage(sig, i, j, k, 0);
185 Set::Matrix eps_avg = Numeric::Interpolate::NodeToCellAverage(eps, i, j, k, 0);
186 w_elastic += gamma * 0.5 * (sig_avg.transpose() * eps_avg).trace() * psi(i, j, k) * dv;
187 });
188 }
189
190
191 void TagCellsForRefinement(int lev, amrex::TagBoxArray& a_tags, Set::Scalar a_time, int a_ngrow) override
192 {
194 Base::Mechanics<MODEL>::TagCellsForRefinement(lev, a_tags, a_time, a_ngrow);
195
196 Set::Vector DX(geom[lev].CellSize());
197 Set::Scalar DXnorm = DX.lpNorm<2>();
198 for (amrex::MFIter mfi(*model_mf[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi)
199 {
200 amrex::Box bx = mfi.tilebox();
201 amrex::Array4<char> const& tags = a_tags.array(mfi);
202 if (psi_on)
203 {
204 amrex::Array4<Set::Scalar> const& psi = psi_mf[lev]->array(mfi);
205 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k)
206 {
207 auto sten = Numeric::GetStencil(i, j, k, bx);
208 {
209 Set::Vector gradpsi = Numeric::Gradient(psi, i, j, k, 0, DX.data(), sten);
210 if (gradpsi.lpNorm<2>() * DXnorm > m_eta_ref_threshold)
211 tags(i, j, k) = amrex::TagBox::SET;
212 }
213 });
214 }
215 }
216 }
217
218
219protected:
220 MODEL model;
225
226
227 using Base::Mechanics<MODEL>::m_type;
228 using Base::Mechanics<MODEL>::finest_level;
229 using Base::Mechanics<MODEL>::geom;
230 using Base::Mechanics<MODEL>::model_mf;
231 using Base::Mechanics<MODEL>::psi_mf;
232 using Base::Mechanics<MODEL>::psi_on;
233 using Base::Mechanics<MODEL>::stress_mf;
234 using Base::Mechanics<MODEL>::strain_mf;
235
236
240 //Set::Scalar L = 1.0;
243 //Set::Scalar lambda = 1.0;
245
250};
251
252
253
254
255
256
257
258
259
260
261} // namespace Integrator
262#endif
#define pp_query_default(...)
Definition ParmParse.H:100
#define pp_queryclass(...)
Definition ParmParse.H:107
#define INFO
Definition Util.H:20
Definition BC.H:42
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
bool contains(std::string name)
Definition ParmParse.H:154
void queryclass(std::string name, T *value, std::string file="", std::string func="", int line=-1)
Definition ParmParse.H:604
void select(std::string name, PTRTYPE *&ic_eta, Args &&... args)
Definition ParmParse.H:692
Set::Field< Set::Matrix > strain_mf
Definition Mechanics.H:475
void Integrate(int amrlev, Set::Scalar, int, const amrex::MFIter &mfi, const amrex::Box &a_box) override
Definition Mechanics.H:380
static void Parse(Mechanics &value, IO::ParmParse &pp)
Definition Mechanics.H:42
void Initialize(int lev) override
Use the #ic object to initialize::Temp.
Definition Mechanics.H:143
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition Mechanics.H:241
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar, int) override
Definition Mechanics.H:439
Set::Field< Set::Matrix > stress_mf
Definition Mechanics.H:474
Set::Field< Set::Scalar > psi_mf
Definition Mechanics.H:465
Set::Field< MODEL > model_mf
Definition Mechanics.H:463
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.
std::vector< amrex::Box > box
Definition Integrator.H:446
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition Integrator.H:381
void RegisterIntegratedVariable(Set::Scalar *integrated_variable, std::string name, bool extensive=true)
Register a variable to be integrated over the spatial domain using the Integrate function.
Set::Scalar w_bndry
Definition TopOp.H:248
void Integrate(int amrlev, Set::Scalar time, int step, const amrex::MFIter &mfi, const amrex::Box &box) override
Definition TopOp.H:166
BC::BC< Set::Scalar > * bc_psi
Definition TopOp.H:222
Set::Scalar volume
Definition TopOp.H:246
Set::Scalar m_eta_ref_threshold
Definition TopOp.H:223
Set::Scalar volume0
Definition TopOp.H:242
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Definition TopOp.H:124
Set::Scalar alpha
Definition TopOp.H:237
Set::Scalar w_chem_potential
Definition TopOp.H:247
virtual void UpdateModel(int a_step, Set::Scalar) override
Definition TopOp.H:107
Set::Scalar w_elastic
Definition TopOp.H:249
Set::Scalar gamma
Definition TopOp.H:239
Set::Field< Set::Scalar > psi_old_mf
Definition TopOp.H:224
void Initialize(int lev) override
Definition TopOp.H:100
Numeric::Interpolator::Linear< Set::Scalar > lambda
Definition TopOp.H:244
Numeric::Interpolator::Linear< Set::Scalar > L
Definition TopOp.H:241
Set::Scalar beta
Definition TopOp.H:238
static void Parse(TopOp &value, IO::ParmParse &pp)
Definition TopOp.H:57
TopOp(IO::ParmParse &pp)
Definition TopOp.H:44
IC::IC< Set::Scalar > * ic_psi
Definition TopOp.H:221
void TagCellsForRefinement(int lev, amrex::TagBoxArray &a_tags, Set::Scalar a_time, int a_ngrow) override
Definition TopOp.H:191
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
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:36
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
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition Base.H:23
AMREX_FORCE_INLINE void RealFillBoundary(amrex::FabArray< amrex::BaseFab< T > > &a_mf, const amrex::Geometry &)
Definition Util.H:322
void Abort(const char *msg)
Definition Util.cpp:170
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:1331