1 #include "lift_drag.hpp" 5 template <
int dim,
int nstate,
typename real,
typename MeshType>
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())
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;
26 template <
int dim,
int nstate,
typename real,
typename MeshType>
30 const double ref_length = this->euler_fad_fad.ref_length;
31 const double dynamic_pressure_inf = this->euler_fad_fad.dynamic_pressure_inf;
33 return 1.0 / (ref_length * dynamic_pressure_inf);
36 template <
int dim,
int nstate,
typename real,
typename MeshType>
40 dealii::Tensor<2,dim,double> rotation_matrix;
41 if constexpr (dim == 1) {
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);
50 if constexpr (dim == 3) {
51 rotation_matrix[0][2] = 0.0;
52 rotation_matrix[1][2] = 0.0;
54 rotation_matrix[2][0] = 0.0;
55 rotation_matrix[2][1] = 0.0;
56 rotation_matrix[2][2] = 1.0;
59 return rotation_matrix;
62 template <
int dim,
int nstate,
typename real,
typename MeshType>
66 dealii::Tensor<1,dim,double> lift_direction;
67 lift_direction[0] = 0.0;
68 lift_direction[1] = 1.0;
70 if constexpr (dim == 1) {
73 if constexpr (dim == 3) {
74 lift_direction[2] = 0.0;
77 dealii::Tensor<1,dim,double> vec;
78 vec = rotation_matrix * lift_direction;
83 template <
int dim,
int nstate,
typename real,
typename MeshType>
87 dealii::Tensor<1,dim,double> drag_direction;
88 drag_direction[0] = 1.0;
89 drag_direction[1] = 0.0;
91 if constexpr (dim == 1) {
94 if constexpr (dim == 3) {
95 drag_direction[2] = 0.0;
98 dealii::Tensor<1,dim,double> vec;
99 vec = rotation_matrix * drag_direction;
104 template <
int dim,
int nstate,
typename real,
typename MeshType>
107 const bool compute_dIdX,
108 const bool compute_d2I)
112 if (functional_type == Functional_types::lift) {
113 this->pcout <<
"Lift value: " << value <<
"\n";
117 if (functional_type == Functional_types::drag) {
118 this->pcout <<
"Drag value: " << value <<
"\n";
Sacado::Fad::DFad< FadType > FadFadType
Sacado AD type that allows 2nd derivatives.
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.
double initialize_force_dimensionalization_factor()
Compute force dimensionalization factor.
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.
Files for the baseline physics.
dealii::Tensor< 2, dim, double > initialize_rotation_matrix(const double angle_of_attack)
Initialize rotation matrix based on given angle of attack.
LiftDragFunctional(std::shared_ptr< DGBase< dim, real, MeshType >> dg_input, const Functional_types functional_type)
Constructor.
Functional_types
Switch between lift and drag functional types.
real evaluate_functional(const bool compute_dIdW=false, const bool compute_dIdX=false, const bool compute_d2I=false) override
Destructor.
DGBase is independent of the number of state variables.
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.