Alamo
TabulatedInterface.H
Go to the documentation of this file.
1//
2// Initialize a 2D interface in a 2 component field defind by (x,y)
3// which are an array of points along the interface.
4//
5//
6
7#ifndef IC_TABULATEDINTERFACE_H_
8#define IC_TABULATEDINTERFACE_H_
9
10#include "IO/ParmParse.H"
11#include "IC/IC.H"
12
13#include "IO/ParmParse.H"
15
16namespace IC
17{
18/// Initialize a perturbed interface using a linear interpolation
19class TabulatedInterface : public IC<Set::Scalar>
20{
21public:
23 TabulatedInterface (amrex::Vector<amrex::Geometry> &_geom) :
24 IC(_geom) {}
25 TabulatedInterface (amrex::Vector<amrex::Geometry> &_geom, IO::ParmParse &pp) :
26 IC(_geom) {pp_queryclass(*this);}
27 TabulatedInterface (amrex::Vector<amrex::Geometry> &_geom, IO::ParmParse &pp, std::string name) :
28 IC(_geom) {pp_queryclass(name,*this);}
29
30 TabulatedInterface (amrex::Vector<amrex::Geometry> &_geom,
31 std::vector<Set::Scalar> a_xs, std::vector<Set::Scalar> a_ys) :
32 IC(_geom)
33 {
34 Define(a_xs,a_ys);
35 }
36
37 TabulatedInterface (amrex::Vector<amrex::Geometry> &_geom,
38 std::vector<Set::Scalar> a_xs, std::vector<Set::Scalar> a_ys,
39 Set::Scalar a_alpha1, Set::Scalar a_alpha2, Type a_type = Type::Values) :
40 IC(_geom)
41 {
42 Define(a_xs,a_ys,a_type,a_alpha1,a_alpha2);
43 }
44
45 void Define(std::vector<Set::Scalar> a_xs, std::vector<Set::Scalar> a_ys,
46 Type a_type = Type::Partition, Set::Scalar a_alpha1 = 1.0, Set::Scalar a_alpha2 = 1.0)
47 {
48 type = a_type;
49 alpha1 = a_alpha1;
50 alpha2 = a_alpha2;
51 xs = a_xs;
52 ys = a_ys;
53 }
54
55 virtual void Add(const int &lev, Set::Field<Set::Scalar> &a_field,Set::Scalar) override
56 {
58
59 for (amrex::MFIter mfi(*a_field[lev],amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
60 {
61 amrex::Box bx = mfi.tilebox();
62 bx.grow(a_field[lev]->nGrow());
63 amrex::Array4<Set::Scalar> const& field = a_field[lev]->array(mfi);
64 amrex::ParallelFor (bx,[=] AMREX_GPU_DEVICE(int i, int j, int k) {
65 amrex::Real x = geom[lev].ProbLo()[0] + ((amrex::Real)(i) + 0.5) * geom[lev].CellSize()[0];
66 amrex::Real y = geom[lev].ProbLo()[1] + ((amrex::Real)(j) + 0.5) * geom[lev].CellSize()[1];
67 if (y > f(x))
68 {
69 if (type == Type::Partition)
70 {
71 field(i,j,k,0) += 0.;
72 field(i,j,k,1) += alpha2;
73 }
74 else
75 field(i,j,k,0) += alpha2;
76 }
77 else
78 {
79 if (type == Type::Partition)
80 {
81 field(i,j,k,0) += alpha1;
82 field(i,j,k,1) += 0.;
83 }
84 else
85 field(i,j,k,0) += alpha1;
86 }
87 });
88 }
89 };
90
91private:
93 std::vector<Set::Scalar> xs;
94 std::vector<Set::Scalar> ys;
96
97public:
98 static void Parse(TabulatedInterface & value, IO::ParmParse & pp)
99 {
100 pp_queryarr("xs",value.xs); // x location of points
101 pp_queryarr("ys",value.ys); // y location of points
102 value.Define(value.xs,value.ys);
103 }
104};
105}
106#endif
#define pp_queryarr(...)
Definition ParmParse.H:103
#define pp_queryclass(...)
Definition ParmParse.H:107
amrex::Vector< amrex::Geometry > & geom
Definition IC.H:61
Initialize a perturbed interface using a linear interpolation.
virtual void Add(const int &lev, Set::Field< Set::Scalar > &a_field, Set::Scalar) override
TabulatedInterface(amrex::Vector< amrex::Geometry > &_geom, std::vector< Set::Scalar > a_xs, std::vector< Set::Scalar > a_ys)
static void Parse(TabulatedInterface &value, IO::ParmParse &pp)
TabulatedInterface(amrex::Vector< amrex::Geometry > &_geom, IO::ParmParse &pp, std::string name)
std::vector< Set::Scalar > xs
TabulatedInterface(amrex::Vector< amrex::Geometry > &_geom)
void Define(std::vector< Set::Scalar > a_xs, std::vector< Set::Scalar > a_ys, Type a_type=Type::Partition, Set::Scalar a_alpha1=1.0, Set::Scalar a_alpha2=1.0)
std::vector< Set::Scalar > ys
TabulatedInterface(amrex::Vector< amrex::Geometry > &_geom, std::vector< Set::Scalar > a_xs, std::vector< Set::Scalar > a_ys, Set::Scalar a_alpha1, Set::Scalar a_alpha2, Type a_type=Type::Values)
TabulatedInterface(amrex::Vector< amrex::Geometry > &_geom, IO::ParmParse &pp)
Initialize a spherical inclusion.
Definition BMP.H:19
amrex::Real Scalar
Definition Base.H:19
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:20