Alamo
SH.H
Go to the documentation of this file.
1 #ifndef MODEL_INTERFACE_GB_SH_H
2 #define MODEL_INTERFACE_GB_SH_H
3 
4 #include <iostream>
5 #include <fstream>
6 
7 #include "AMReX.H"
8 #include "GB.H"
9 #include "Set/Set.H"
10 #include "Util/Util.H"
11 
12 #define PI 3.14159265
13 
14 namespace Model
15 {
16 namespace Interface
17 {
18 namespace GB
19 {
20 ///
21 /// A 2D interface model class. Under construction
22 ///
23 /// Here are some useful trig identities for reference:
24 ///
25 /// \f{align}{\sin(2 \arcsin x) &= 2 x \sqrt{1-x^2} &
26 /// \cos(2 \arcsin x) &= 1 - 2 x^2 \\ %
27 /// \sin(2 \arccos x) &= 2 x \sqrt{1-x^2} &
28 /// \cos(2 \arccos x) &= 2 x^2 - 1 \\ %
29 /// \sin(2 \arctan x) &= \frac{2 x}{1 + x^2} &
30 /// \cos(2 \arctan x) &= \frac{1-x^2}{1+x^2} \f}
31 /// Specializing for tangent:
32 /// \f{align}{\sin(2 \arctan y/x) &= \frac{2 x y}{x^2 + y^2} &
33 /// \cos(2 \arctan y/x) &= \frac{x^2-y^2}{x^2+y^2} \f}
34 ///
35 class SH : public GB
36 {
37 public:
38  static constexpr const char* name = "sh";
39 
40  SH(){};
42  SH(IO::ParmParse &pp,std::string name){pp_queryclass(name,*this);};
43  SH(const amrex::Real a_theta0, const amrex::Real a_phi0, const amrex::Real a_sigma0, const amrex::Real a_sigma1)
44  {
45  Define(a_theta0,a_phi0,a_sigma0,a_sigma1);
46  };
47  void Define(const amrex::Real a_theta0, const amrex::Real a_phi0, const amrex::Real a_sigma0, const amrex::Real a_sigma1)
48  {
49  theta0 = a_theta0;
50  phi0 = a_phi0;
51  sigma0 = a_sigma0;
52  sigma1 = a_sigma1;
53  };
54  void Randomize()
55  {
57  sigma0 = Util::Random();
58  sigma1 = Util::Random();
59  };
60  Set::Scalar W(const Set::Scalar a_theta, const Set::Scalar a_phi) const
61  {
62  return sigma0 + sigma1*(1.0 - (cos(2*(a_phi))*cos(2*(a_phi)) * sin(2*(a_theta))*sin(2*(a_theta))));
63  };
64  std::array<Set::Scalar,2> DW(const Set::Scalar a_theta, const Set::Scalar a_phi) const
65  {
66  return {- sigma1 * 4.0 * cos(2*a_phi) * cos(2*a_phi) * sin(2*a_theta) * cos(2*a_theta),
67  + sigma1 * 4.0 * sin(2*a_phi) * cos(2*a_phi) * sin(2*a_theta) * sin(2*a_theta)} ;
68  };
69 
70  virtual Set::Scalar W(const Set::Vector &a_n) const override
71  {
72  Set::Vector n = a_n / a_n.lpNorm<2>();
73  Set::Scalar theta = std::acos(n(0));
74  Set::Scalar phi = std::atan2(n(2),n(1));
75  return W(theta,phi);
76  };
77  virtual Set::Scalar DW(const Set::Vector &a_n, const Set::Vector &t_n) const override
78  {
79  const Set::Scalar alpha = 1E-4;
80  return (W(a_n + alpha*t_n) - W(a_n - alpha*t_n)) / 2.0 / alpha;
81  };
82  virtual Set::Scalar DDW(const Set::Vector &a_n, const Set::Vector &t_n) const override
83  {
84  const Set::Scalar alpha = 1E-4;
85  return (W(a_n + alpha*t_n) - 2.0*W(a_n) + W(a_n - alpha*t_n)) / alpha / alpha;
86  };
87 
88 
89  Set::Scalar W(const Set::Scalar ) const override {return NAN;};
90  Set::Scalar DW(const Set::Scalar ) const override {return NAN;};
91  Set::Scalar DDW(const Set::Scalar ) const override {return NAN;};
92 
93 private:
94  Set::Scalar theta0 = NAN, phi0 = NAN, sigma0 = NAN, sigma1 = NAN;
95 
96 public:
97  static void Parse(SH & value, amrex::ParmParse & pp)
98  {
99  pp_query("theta0",value.theta0); // Theta offset (degrees)
100  value.theta0 *= 0.01745329251; // convert degrees into radians
101  pp_query("phi0",value.phi0); // Phi offset (radians)
102  value.phi0 *= 0.01745329251; // convert degrees into radians
103  pp_query("sigma0",value.sigma0); // Minimum energy value
104  pp_query("sigma1",value.sigma1); // Energy multiplier
105 
106  std::string reg_str = "wilhelm";
107  pp_query("regularization",reg_str); // Type of regularization to use: {wilhelm,k23}
108  if (reg_str == "wilhelm") value.regularization = Regularization::Wilhelm;
109  else if (reg_str == "k23") value.regularization = Regularization::K23;
110  else Util::Abort(INFO,"Invalid regularization string ",reg_str);
111  }
112 
113 };
114 }
115 }
116 }
117 #endif
Model::Interface::GB::SH::SH
SH(const amrex::Real a_theta0, const amrex::Real a_phi0, const amrex::Real a_sigma0, const amrex::Real a_sigma1)
Definition: SH.H:43
Model::Interface::GB::SH::sigma1
Set::Scalar sigma1
Definition: SH.H:94
Model::Interface::GB::SH::Define
void Define(const amrex::Real a_theta0, const amrex::Real a_phi0, const amrex::Real a_sigma0, const amrex::Real a_sigma1)
Definition: SH.H:47
Model::Interface::GB::SH::DW
virtual Set::Scalar DW(const Set::Vector &a_n, const Set::Vector &t_n) const override
Definition: SH.H:77
Util.H
Model::Interface::GB::GB
Definition: GB.H:18
Model::Interface::GB::GB::regularization
Regularization regularization
Definition: GB.H:141
Model::Interface::GB::SH::Parse
static void Parse(SH &value, amrex::ParmParse &pp)
Definition: SH.H:97
Model::Interface::GB::SH::W
Set::Scalar W(const Set::Scalar) const override
Definition: SH.H:89
Model::Interface::GB::SH::SH
SH(IO::ParmParse &pp, std::string name)
Definition: SH.H:42
Model::Interface::GB::SH::W
Set::Scalar W(const Set::Scalar a_theta, const Set::Scalar a_phi) const
Definition: SH.H:60
Set::Vector
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition: Base.H:20
pp_query
#define pp_query(...)
Definition: ParmParse.H:105
Model::Interface::GB::SH::Randomize
void Randomize()
Definition: SH.H:54
Util::Random
Set::Scalar Random()
Definition: Set.cpp:9
Set::Scalar
amrex::Real Scalar
Definition: Base.H:19
pp_queryclass
#define pp_queryclass(...)
Definition: ParmParse.H:106
Model::Interface::GB::SH::DDW
Set::Scalar DDW(const Set::Scalar) const override
Definition: SH.H:91
Model::Interface::GB::SH::DW
std::array< Set::Scalar, 2 > DW(const Set::Scalar a_theta, const Set::Scalar a_phi) const
Definition: SH.H:64
Model::Interface::GB::SH::W
virtual Set::Scalar W(const Set::Vector &a_n) const override
Definition: SH.H:70
GB.H
Model::Interface::GB::SH::theta0
Set::Scalar theta0
Definition: SH.H:94
Model::Interface::GB::SH::name
static constexpr const char * name
Definition: SH.H:38
Model::Interface::GB::SH
Definition: SH.H:35
Model::Interface::GB::SH::SH
SH()
Definition: SH.H:40
Util::Abort
void Abort(const char *msg)
Definition: Util.cpp:166
Model::Interface::GB::SH::DW
Set::Scalar DW(const Set::Scalar) const override
Definition: SH.H:90
Set.H
Set::Constant::Pi
static const Set::Scalar Pi
Definition: Set.H:286
IO::ParmParse
Definition: ParmParse.H:112
Model::Interface::GB::SH::sigma0
Set::Scalar sigma0
Definition: SH.H:94
INFO
#define INFO
Definition: Util.H:20
Model::Interface::GB::SH::DDW
virtual Set::Scalar DDW(const Set::Vector &a_n, const Set::Vector &t_n) const override
Definition: SH.H:82
Model
Definition: Constant.H:12
Model::Interface::GB::SH::phi0
Set::Scalar phi0
Definition: SH.H:94
Model::Interface::GB::SH::SH
SH(IO::ParmParse &pp)
Definition: SH.H:41