Alamo
Laminate.H
Go to the documentation of this file.
1 //
2 // Create a single laminate with specified orientation, thickness, and offset.
3 //
4 
5 #ifndef IC_LAMINATE_H_
6 #define IC_LAMINATE_H_
7 
8 #include "IC/IC.H"
9 #include "Util/Util.H"
10 #include "IO/ParmParse.H"
11 
12 /// \class Laminate
13 /// \brief Initialize Laminates in a matrix
14 namespace IC
15 {
16 class Laminate: public IC
17 {
18 public:
20 
21  Laminate(amrex::Vector<amrex::Geometry>& _geom): IC(_geom) {}
22  Laminate(amrex::Vector<amrex::Geometry>& _geom, IO::ParmParse& pp, std::string name): IC(_geom)
23  {
24  pp_queryclass(name, *this);
25  }
26 
27  void Add(const int& lev, Set::Field<Set::Scalar>& a_field, Set::Scalar)
28  {
29  for (amrex::MFIter mfi(*a_field[lev], amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
30  {
31  amrex::Box bx = mfi.tilebox();
32  bx.grow(a_field[lev]->nGrow());
33  amrex::IndexType type = a_field[lev]->ixType();
34 
35  amrex::Array4<Set::Scalar> const& field = a_field[lev]->array(mfi);
36  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
37 
38  Set::Vector x = Set::Position(i, j, k, geom[lev], type);
39 
40  // The base matrix is always field(i,j,k,0)
41  // Inclusions start with field(i,j,k,1)
42 
43  if (!singlefab)
44  {
45  Set::Scalar value = 0.0;
46  for (int m = 0; m < number_of_inclusions; m++)
47  {
48  Set::Scalar t = std::abs((x - center[m]).transpose() * orientation[m]) - (thickness[m] / 2.0);
49  field(i, j, k, m + 1) = 0.5 - 0.5 * std::erf(t / eps[m]);
50  if (field(i, j, k, m + 1) < 0.) field(i, j, k, m + 1) = 0.;
51  if (field(i, j, k, m + 1) > 1.) field(i, j, k, m + 1) = 1.;
52  value += field(i, j, k, m + 1);
53  }
54  field(i, j, k, 0) = 1.0 - value;
55  if (field(i, j, k, 0) < 0.) field(i, j, k, 0) = 0.;
56  if (field(i, j, k, 0) > 1.) field(i, j, k, 0) = 1.;
57  }
58  else
59  {
60  Set::Scalar t = std::abs((x - center[0]).transpose() * orientation[0]) - (thickness[0] / 2.0);
61  field(i, j, k, 0) = 0.5 - 0.5 * std::erf(t / eps[0]);
62  if (invert) field(i, j, k, 0) = 1.0 - field(i, j, k, 0);
63  }
64  });
65  }
66  a_field[lev]->FillBoundary();
67  };
68 
69 private:
71  amrex::Vector<Set::Vector> center;
72  amrex::Vector<Set::Vector> orientation;
73  amrex::Vector<Set::Vector> normal;
74  amrex::Vector<Set::Scalar> eps;
75  amrex::Vector<Set::Scalar> thickness;
77  bool singlefab = false;
78  bool invert = false;
79 
80 public:
81  static void Parse(Laminate& value, IO::ParmParse& pp)
82  {
83  // How many laminates (MUST be greater than or equal to 1).
84  pp_query_default("number_of_inclusions", value.number_of_inclusions, 1);
85 
86  if (value.number_of_inclusions < 1)
87  Util::Abort(INFO, "Number of inclusions must be at least 1. Aborting.");
88 
89  amrex::Vector<Set::Scalar> a_center;
90  // (x,y,[z]) values for the center point of the laminate
91  if (pp.contains("center")) pp_queryarr("center", a_center);
92 
93  if (a_center.size() != value.number_of_inclusions * AMREX_SPACEDIM) value.center.push_back(Set::Vector::Zero());
94  else
95  {
96  for (int i = 0; i < a_center.size(); i += AMREX_SPACEDIM)
97  value.center.push_back(Set::Vector(AMREX_D_DECL(a_center[i], a_center[i + 1], a_center[i + 2])));
98  }
99 
100  amrex::Vector<Set::Scalar> a_thickness;
101  // thickness of the laminate
102  if (pp.contains("thickness")) pp_queryarr("thickness", a_thickness);
103 
104  if (a_thickness.size() != value.number_of_inclusions && a_thickness.size() != 1)
105  Util::Abort(INFO, "Thickness of each inclusion must be specified");
106 
107  if (a_thickness.size() == 1)
108  {
109  if (a_thickness[0] <= 0.0) Util::Abort(INFO, "Invalid value of inclusion thickness");
110  for (int i = 0; i < value.number_of_inclusions; i++) value.thickness.push_back(a_thickness[0]);
111  }
112 
113  else
114  {
115  for (int i = 0; i < value.number_of_inclusions; i++)
116  {
117  if (a_thickness[i] <= 0.0) Util::Abort(INFO, "Invalid value of inclusion ", i + 1, " thickness");
118  value.thickness.push_back(a_thickness[i]);
119  }
120  }
121 
122  amrex::Vector<Set::Scalar> a_orientation;
123  // Vector normal to the interface of the laminate
124  pp_queryarr("orientation", a_orientation);
125 
126  if (a_orientation.size() != value.number_of_inclusions * AMREX_SPACEDIM && a_orientation.size() != AMREX_SPACEDIM)
127  Util::Abort(INFO, "Orientation of each inclusion must be specified");
128 
129  if (a_orientation.size() == AMREX_SPACEDIM)
130  for (int i = 0; i < value.number_of_inclusions; i++)
131  value.orientation.push_back(Set::Vector(AMREX_D_DECL(a_orientation[0], a_orientation[1], a_orientation[2])));
132 
133  else
134  for (int i = 0; i < a_orientation.size(); i += AMREX_SPACEDIM)
135  value.orientation.push_back(Set::Vector(AMREX_D_DECL(a_orientation[i], a_orientation[i + 1], a_orientation[i + 2])));
136 
137  for (int i = 0; i < value.orientation.size(); i++)
138  if (value.orientation[i].lpNorm<2>() <= 0.) value.orientation[i] = Set::Vector::Random();
139 
140  for (int i = 0; i < value.orientation.size(); i++)
141  value.orientation[i] = value.orientation[i] / value.orientation[i].lpNorm<2>();
142 
143  value.normal.resize(value.number_of_inclusions);
144  for (int i = 0; i < value.orientation.size(); i++)
145  {
146  value.normal[i] = Set::Vector::Zero();
147  if (value.orientation[i](0) != 0.)
148  {
149  AMREX_D_TERM(value.normal[i](0) = 1.;, value.normal[i](1) = 1.;, value.normal[i](2) = 1.;);
150  value.normal[i](0) = -(AMREX_D_TERM(0., +value.orientation[i](1), +value.orientation[i](2))) / value.orientation[i](0);
151  value.normal[i] = value.normal[i] / value.normal[i].lpNorm<2>();
152  }
153  else if (value.orientation[i](1) != 0.)
154  {
155  AMREX_D_TERM(value.normal[i](0) = 1.;, value.normal[i](1) = 1.;, value.normal[i](2) = 1.;);
156  value.normal[i](1) = -(AMREX_D_TERM(value.orientation[i](0), +0.0, +value.orientation[i](2))) / value.orientation[i](1);
157  value.normal[i] = value.normal[i] / value.normal[i].lpNorm<2>();
158  }
159  }
160 
161  amrex::Vector<Set::Scalar> a_eps;
162  // Diffuse thickness
163  pp_queryarr("eps", a_eps);
164 
165  if (a_eps.size() < 1)
166  for (int i = 0; i < value.number_of_inclusions; i++)
167  value.eps.push_back(1.e-5);
168  if (a_eps.size() == 1)
169  {
170  if (a_eps[0] < 0.0)
171  {
172  Util::Warning(INFO, "Invalid value of laminate.eps. Resetting to 1e-5");
173  a_eps[0] = 1.e-5;
174  }
175  for (int i = 0; i < value.number_of_inclusions; i++)
176  value.eps.push_back(a_eps[0]);
177  }
178  else
179  {
180  for (int i = 0; i < value.number_of_inclusions; i++)
181  value.eps.push_back(a_eps[i]);
182  }
183 
184  std::string mollifier;
185  pp_query("mollifier", mollifier); // Type of mollifer to use (options: dirac, [gaussian])
186  if (mollifier == "Dirac" || mollifier == "dirac")
187  value.moll = Mollifier::Dirac;
188  else
189  value.moll = Mollifier::Gaussian;
190 
191  // Switch to mode where only one component is used.
192  pp_query("singlefab", value.singlefab);
193 
194  // Take the complement of the laminate
195  pp_query("invert", value.invert);
196  }
197 };
198 }
199 #endif
IC::IC::geom
amrex::Vector< amrex::Geometry > & geom
Definition: IC.H:45
Set::Position
AMREX_FORCE_INLINE Vector Position(const int &i, const int &j, const int &k, const amrex::Geometry &geom, const amrex::IndexType &ixType)
Definition: Base.H:121
IC::Laminate::Mollifier
Mollifier
Definition: Laminate.H:19
BC::AMREX_D_DECL
@ AMREX_D_DECL
Definition: BC.H:34
Util.H
IC::Laminate::Parse
static void Parse(Laminate &value, IO::ParmParse &pp)
Definition: Laminate.H:81
IC::Laminate::Laminate
Laminate(amrex::Vector< amrex::Geometry > &_geom, IO::ParmParse &pp, std::string name)
Definition: Laminate.H:22
t
std::time_t t
Definition: FileNameParse.cpp:12
Set::Field< Set::Scalar >
Definition: Set.H:236
IC::Laminate::invert
bool invert
Definition: Laminate.H:78
IC::Laminate::eps
amrex::Vector< Set::Scalar > eps
Definition: Laminate.H:74
Set::Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition: Base.H:20
ParmParse.H
pp_query
#define pp_query(...)
Definition: ParmParse.H:104
IC::Laminate::normal
amrex::Vector< Set::Vector > normal
Definition: Laminate.H:73
Util::Gaussian
Set::Scalar Gaussian(amrex::Real mean, amrex::Real std_deviation)
Definition: Set.cpp:13
Util::Random
Set::Scalar Random()
Definition: Set.cpp:9
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
IC::Laminate::Dirac
@ Dirac
Definition: Laminate.H:19
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:105
IC::Laminate::moll
Mollifier moll
Definition: Laminate.H:76
IC::Laminate::singlefab
bool singlefab
Definition: Laminate.H:77
IC::Laminate::Add
void Add(const int &lev, Set::Field< Set::Scalar > &a_field, Set::Scalar)
Definition: Laminate.H:27
IC::Laminate::Laminate
Laminate(amrex::Vector< amrex::Geometry > &_geom)
Definition: Laminate.H:21
IO::ParmParse::contains
bool contains(std::string name)
Definition: ParmParse.H:151
pp_query_default
#define pp_query_default(...)
Definition: ParmParse.H:100
IC::Laminate::center
amrex::Vector< Set::Vector > center
Definition: Laminate.H:71
IC::Laminate::number_of_inclusions
int number_of_inclusions
Definition: Laminate.H:70
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:165
IC::Laminate::Gaussian
@ Gaussian
Definition: Laminate.H:19
IC::Laminate::orientation
amrex::Vector< Set::Vector > orientation
Definition: Laminate.H:72
IC
Definition: Affine.H:18
Util::Warning
void Warning(std::string file, std::string func, int line, Args const &... args)
Definition: Util.H:171
IO::ParmParse
Definition: ParmParse.H:110
IC::Laminate::thickness
amrex::Vector< Set::Scalar > thickness
Definition: Laminate.H:75
IC::Laminate
Definition: Laminate.H:16
INFO
#define INFO
Definition: Util.H:20
IC.H
pp_queryarr
#define pp_queryarr(...)
Definition: ParmParse.H:103