Alamo
Ellipsoid.H
Go to the documentation of this file.
1// This IC initializes a single-component fab using the formula
2//
3// .. math::
4//
5// \phi = \begin{cases}
6// \alpha_{in} & (\mathbf{x}-\mathbf{x}_0)^T\mathbf{A}(\mathbf{x}-\mathbf{x}_0) \le r \\
7// \alpha_{out} & \text{else}
8// \end{cases}
9//
10// The boundary can be mollified using an error function with parameter :math:`\varepsilon`
11//
12// .. WARNING::
13//
14// This IC is redundant with :ref:`IC::Ellipse` and :ref:`IC::Expression` and will probably be consolidated
15//
16
17#ifndef IC_ELLIPSOID_H_
18#define IC_ELLIPSOID_H_
19
20#include "IC/IC.H"
21#include "Util/Util.H"
22#include "IO/ParmParse.H"
23
24/// \class Ellipsoid
25/// \brief Initialize an ellipsoidal inclusion
26namespace IC
27{
28class Ellipsoid : public IC<Set::Scalar>
29{
30public:
31 static constexpr const char* name = "ellipsoid";
32
34
35 Ellipsoid (amrex::Vector<amrex::Geometry> &_geom) : IC(_geom) {}
36 Ellipsoid (amrex::Vector<amrex::Geometry>& _geom, IO::ParmParse& pp, std::string name) : Ellipsoid(_geom)
37 {
38 pp_queryclass(name, *this);
39 }
40
41
42 void Add(const int &lev,Set::Field<Set::Scalar> &a_field, Set::Scalar)
43 {
44 for (amrex::MFIter mfi(*a_field[lev],amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
45 {
46 amrex::Box bx = mfi.tilebox();
47 bx.grow(a_field[lev]->nGrow());
48 amrex::IndexType type = a_field[lev]->ixType();
49
50 amrex::Array4<Set::Scalar> const& field = a_field[lev]->array(mfi);
51 amrex::ParallelFor (bx,[=] AMREX_GPU_DEVICE(int i, int j, int k) {
52
53 Set::Vector x = Set::Position(i,j,k,geom[lev],type);
54
55 Set::Scalar min_value = field(i,j,k);
56 for (int i=0 ; i < center.size(); i++)
57 {
58 Set::Scalar value = 0.;
59 Set::Scalar norm = (A[i]*(x-center[i])).lpNorm<2>();
60 value = 0.5 + 0.5*std::erf(((x-center[i]).transpose() * A[i] * (x-center[i]) - 1.0) / eps[i] / norm);
61 value = in_value + value*(out_value - in_value);
62 min_value = value < min_value ? value : min_value;
63 }
64 field(i,j,k) = min_value;
65 if (field(i,j,k) < std::min(in_value,out_value)) field(i,j,k) = std::min(in_value,out_value);
66 if (field(i,j,k) > std::max(in_value,out_value)) field(i,j,k) = std::max(in_value,out_value);
67 });
68 }
69 a_field[lev]->FillBoundary();
70 };
71
72private:
74 amrex::Vector<Set::Vector> center;
75 amrex::Vector<Set::Scalar> eps;
76 amrex::Vector<Set::Matrix> A;
78
79public:
80 static void Parse(Ellipsoid & value, IO::ParmParse & pp)
81 {
82 amrex::Vector<Set::Scalar> center;
83 if(pp.contains("center")) pp_queryarr("center",center); // Center of the ellipse :math:`\mathbf{x}_0`
84
85 if(center.size() < AMREX_SPACEDIM) value.center.push_back(Set::Vector::Zero());
86 else
87 {
88 for (int i = 0; i<center.size(); i+=AMREX_SPACEDIM)
89 value.center.push_back(Set::Vector(AMREX_D_DECL(center[i],center[i+1],center[i+2])));
90 }
91
92 amrex::Vector<Set::Scalar> A;
93 if (pp.contains("A"))
94 {
95 pp_queryarr("A", A); // Matrix defining elipse radii and orientation
96 if(A.size() < AMREX_SPACEDIM*AMREX_SPACEDIM) value.A.push_back(Set::Matrix::Identity());
97 else
98 {
99 int i=0, j=0;
100 Set::Matrix temp = Set::Matrix::Zero();
101 for (i=0; i<A.size(); i+=AMREX_SPACEDIM)
102 {
103 for (j=0; j<AMREX_SPACEDIM; j++)
104 temp(i,j) = A[i+j];
105
106 if((i+j) % (AMREX_SPACEDIM*AMREX_SPACEDIM) == 0)
107 {
108 value.A.push_back(temp);
109 temp = Set::Matrix::Zero();
110 }
111 }
112 }
113 }
114 if (pp.contains("radius"))
115 {
116 amrex::Vector<Set::Scalar> a_radius;
117 pp_queryarr("radius",a_radius); // "Vector of radii (use instead of A)"
118
119 if (a_radius.size() < AMREX_SPACEDIM) value.A.push_back(Set::Matrix::Identity());
120 else
121 {
122 for (int i = 0; i<a_radius.size(); i+=AMREX_SPACEDIM)
123 {
124 Set::Matrix temp = Set::Matrix::Zero();
125 AMREX_D_TERM( temp(0,0) = 1./(a_radius[i]*a_radius[i]) ;, temp(1,1) = 1./(a_radius[i+1]*a_radius[i+1]) ;, temp(2,2) = 1./(a_radius[i+2]*a_radius[i+2]) ; );
126 value.A.push_back(temp);
127 }
128 }
129 }
130
131 amrex::Vector<Set::Scalar> eps;
132 pp_queryarr("eps", eps); // Mollifying value for erf
133
134 if(eps.size() < 1) value.eps.push_back(1.e-5);
135 for (int i =0; i < value.A.size(); i++)
136 {
137 if (eps[i] <= 0.0) eps[i] = 1.e-5;
138 value.eps.push_back(eps[i]);
139 }
140 //if(value.eps <= 0.) value.eps = 1.e-5;
141
142 pp_query_default("in_value", value.in_value, 0.0); // Value of field inside ellipse
143 pp_query_default("out_value", value.out_value, 1.0); // Value of field outside ellipse
144
145 std::string mollifier;
146 pp_query("mollifier",mollifier); // Type of mollifier to use (options: dirac, [gaussian])
147 if(mollifier == "Dirac" || mollifier == "dirac")
148 value.moll = Mollifier::Dirac;
149 else
151 }
152};
153}
154#endif
#define pp_queryarr(...)
Definition ParmParse.H:105
#define pp_query(...)
Definition ParmParse.H:108
#define pp_query_default(...)
Definition ParmParse.H:102
#define pp_queryclass(...)
Definition ParmParse.H:109
Ellipsoid(amrex::Vector< amrex::Geometry > &_geom)
Definition Ellipsoid.H:35
amrex::Vector< Set::Scalar > eps
Definition Ellipsoid.H:75
Mollifier moll
Definition Ellipsoid.H:77
static constexpr const char * name
Definition Ellipsoid.H:31
void Add(const int &lev, Set::Field< Set::Scalar > &a_field, Set::Scalar)
Definition Ellipsoid.H:42
Set::Scalar out_value
Definition Ellipsoid.H:73
Set::Scalar in_value
Definition Ellipsoid.H:73
Ellipsoid(amrex::Vector< amrex::Geometry > &_geom, IO::ParmParse &pp, std::string name)
Definition Ellipsoid.H:36
static void Parse(Ellipsoid &value, IO::ParmParse &pp)
Definition Ellipsoid.H:80
amrex::Vector< Set::Vector > center
Definition Ellipsoid.H:74
amrex::Vector< Set::Matrix > A
Definition Ellipsoid.H:76
amrex::Vector< amrex::Geometry > & geom
Definition IC.H:61
bool contains(std::string name)
Definition ParmParse.H:173
Initialize a spherical inclusion.
Definition BMP.H:20
amrex::Real Scalar
Definition Base.H:18
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:19
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:120
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, AMREX_SPACEDIM > Matrix
Definition Base.H:22