[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
lift_drag.cpp
1 #include "lift_drag.hpp"
2 
3 namespace PHiLiP {
4 
5 template <int dim,int nstate,typename real,typename MeshType>
8  std::shared_ptr<DGBase<dim,real,MeshType>> dg_input,
9  const Functional_types functional_type)
10  : Functional<dim,nstate,real,MeshType>(dg_input)
11  , functional_type(functional_type)
12  , euler_fad_fad(dynamic_cast< Physics::Euler<dim,dim+2,FadFadType> &>(*(this->physics_fad_fad)))
13  , angle_of_attack(this->euler_fad_fad.angle_of_attack)
14  , rotation_matrix(initialize_rotation_matrix(this->angle_of_attack))
15  , lift_vector(initialize_lift_vector(this->rotation_matrix))
16  , drag_vector(initialize_drag_vector(this->rotation_matrix))
17  , force_dimensionalization_factor(this->initialize_force_dimensionalization_factor())
18 {
19  switch(functional_type) {
20  case Functional_types::lift : this->force_vector = lift_vector; break;
21  case Functional_types::drag : this->force_vector = drag_vector; break;
22  default: break;
23  }
24 }
25 
26 template <int dim,int nstate,typename real,typename MeshType>
29 {
30  const double ref_length = this->euler_fad_fad.ref_length;
31  const double dynamic_pressure_inf = this->euler_fad_fad.dynamic_pressure_inf;
32 
33  return 1.0 / (ref_length * dynamic_pressure_inf);
34 }
35 
36 template <int dim,int nstate,typename real,typename MeshType>
37 dealii::Tensor<2,dim,double> LiftDragFunctional<dim,nstate,real,MeshType>
38 ::initialize_rotation_matrix(const double angle_of_attack)
39 {
40  dealii::Tensor<2,dim,double> rotation_matrix;
41  if constexpr (dim == 1) {
42  assert(false);
43  }
44 
45  rotation_matrix[0][0] = cos(angle_of_attack);
46  rotation_matrix[0][1] = -sin(angle_of_attack);
47  rotation_matrix[1][0] = sin(angle_of_attack);
48  rotation_matrix[1][1] = cos(angle_of_attack);
49 
50  if constexpr (dim == 3) {
51  rotation_matrix[0][2] = 0.0;
52  rotation_matrix[1][2] = 0.0;
53 
54  rotation_matrix[2][0] = 0.0;
55  rotation_matrix[2][1] = 0.0;
56  rotation_matrix[2][2] = 1.0;
57  }
58 
59  return rotation_matrix;
60 }
61 
62 template <int dim,int nstate,typename real,typename MeshType>
63 dealii::Tensor<1,dim,double> LiftDragFunctional<dim,nstate,real,MeshType>
64 ::initialize_lift_vector (const dealii::Tensor<2,dim,double> &rotation_matrix)
65 {
66  dealii::Tensor<1,dim,double> lift_direction;
67  lift_direction[0] = 0.0;
68  lift_direction[1] = 1.0;
69 
70  if constexpr (dim == 1) {
71  assert(false);
72  }
73  if constexpr (dim == 3) {
74  lift_direction[2] = 0.0;
75  }
76 
77  dealii::Tensor<1,dim,double> vec;
78  vec = rotation_matrix * lift_direction;
79 
80  return vec;
81 }
82 
83 template <int dim,int nstate,typename real,typename MeshType>
84 dealii::Tensor<1,dim,double> LiftDragFunctional<dim,nstate,real,MeshType>
85 ::initialize_drag_vector (const dealii::Tensor<2,dim,double> &rotation_matrix)
86 {
87  dealii::Tensor<1,dim,double> drag_direction;
88  drag_direction[0] = 1.0;
89  drag_direction[1] = 0.0;
90 
91  if constexpr (dim == 1) {
92  assert(false);
93  }
94  if constexpr (dim == 3) {
95  drag_direction[2] = 0.0;
96  }
97 
98  dealii::Tensor<1,dim,double> vec;
99  vec = rotation_matrix * drag_direction;
100 
101  return vec;
102 }
103 
104 template <int dim,int nstate,typename real,typename MeshType>
106 ::evaluate_functional(const bool compute_dIdW,
107  const bool compute_dIdX,
108  const bool compute_d2I)
109 {
110  double value = Functional<dim,nstate,real,MeshType>::evaluate_functional(compute_dIdW, compute_dIdX, compute_d2I);
111 
112  if (functional_type == Functional_types::lift) {
113  this->pcout << "Lift value: " << value << "\n";
114  //std::cout << "Lift value: " << value << std::cout;
115  //std::cout << "Lift value: " << value << std::cout;
116  }
117  if (functional_type == Functional_types::drag) {
118  this->pcout << "Drag value: " << value << "\n";
119  }
120 
121  return value;
122 }
123 
124 #if PHILIP_DIM!=1
126 #endif
127 
128 } // PHiLiP namespace
129 
Sacado::Fad::DFad< FadType > FadFadType
Sacado AD type that allows 2nd derivatives.
Definition: lift_drag.hpp:23
dealii::Tensor< 1, dim, double > initialize_drag_vector(const dealii::Tensor< 2, dim, double > &rotation_matrix)
Initialize drag vector with given rotation matrix based on angle of attack.
Definition: lift_drag.cpp:85
double initialize_force_dimensionalization_factor()
Compute force dimensionalization factor.
Definition: lift_drag.cpp:28
dealii::Tensor< 1, dim, double > initialize_lift_vector(const dealii::Tensor< 2, dim, double > &rotation_matrix)
Initialize lift vector with given rotation matrix based on angle of attack.
Definition: lift_drag.cpp:64
Files for the baseline physics.
Definition: ADTypes.hpp:10
dealii::Tensor< 2, dim, double > initialize_rotation_matrix(const double angle_of_attack)
Initialize rotation matrix based on given angle of attack.
Definition: lift_drag.cpp:38
LiftDragFunctional(std::shared_ptr< DGBase< dim, real, MeshType >> dg_input, const Functional_types functional_type)
Constructor.
Definition: lift_drag.cpp:7
Functional_types
Switch between lift and drag functional types.
Definition: lift_drag.hpp:20
Functional base class.
Definition: functional.h:43
real evaluate_functional(const bool compute_dIdW=false, const bool compute_dIdX=false, const bool compute_d2I=false) override
Destructor.
Definition: lift_drag.cpp:106
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
virtual real evaluate_functional(const bool compute_dIdW=false, const bool compute_dIdX=false, const bool compute_d2I=false)
Evaluates the functional derivative with respect to the solution variable.
Definition: functional.cpp:795