VertexCFD  0.0-dev
VertexCFD_Closure_MethodManufacturedSolution_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_METHODMANUFACTUREDSOLUTION_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_METHODMANUFACTUREDSOLUTION_IMPL_HPP
3 
4 #include "utils/VertexCFD_Utils_Constants.hpp"
5 #include "utils/VertexCFD_Utils_VectorField.hpp"
6 
7 #include <Panzer_GlobalIndexer.hpp>
8 #include <Panzer_HierarchicParallelism.hpp>
9 #include <Panzer_PureBasis.hpp>
10 #include <Panzer_Workset_Utilities.hpp>
11 #include <Teuchos_Array.hpp>
12 
13 #include <cmath>
14 #include <string>
15 
16 namespace VertexCFD
17 {
18 namespace ClosureModel
19 {
20 //---------------------------------------------------------------------------//
21 template<class EvalType, class Traits, int NumSpaceDim>
22 MethodManufacturedSolution<EvalType, Traits, NumSpaceDim>::MethodManufacturedSolution(
23  const panzer::IntegrationRule& ir)
24  : _lagrange_pressure("Exact_lagrange_pressure", ir.dl_scalar)
25  , _temperature("Exact_temperature", ir.dl_scalar)
26  , _ir_degree(ir.cubature_degree)
27 
28 {
29  this->addEvaluatedField(_lagrange_pressure);
30  Utils::addEvaluatedVectorField(
31  *this, ir.dl_scalar, _velocity, "Exact_velocity_");
32  this->addEvaluatedField(_temperature);
33 
34  this->setName("Method of Manufactured Solution "
35  + std::to_string(num_space_dim) + "D");
36 
37  _phi_coeff[0] = 0.0125;
38  _phi_coeff[1] = 1.0;
39  _phi_coeff[2] = 0.25;
40  _phi_coeff[3] = 0.5;
41  _phi_coeff[4] = 0.125;
42  _phi_coeff[5] = 0.0;
43 
44  _vel_coeff[0][0] = 0.0125;
45  _vel_coeff[0][1] = 0.08;
46  _vel_coeff[0][2] = 0.125;
47  _vel_coeff[0][3] = 0.0;
48  _vel_coeff[0][4] = 0.125;
49  _vel_coeff[0][5] = 0.25;
50 
51  _vel_coeff[1][0] = 0.0375;
52  _vel_coeff[1][1] = 1.125;
53  _vel_coeff[1][2] = 0.25;
54  _vel_coeff[1][3] = 0.0;
55  _vel_coeff[1][4] = 0.375;
56  _vel_coeff[1][5] = 0.5;
57 
58  _T_coeff[0] = 0.0625;
59  _T_coeff[1] = 1.0;
60  _T_coeff[2] = 0.375;
61  _T_coeff[3] = 0.25;
62  _T_coeff[4] = 0.25;
63  _T_coeff[5] = 0.5;
64 
65  if (num_space_dim == 3)
66  {
67  _phi_coeff[6] = 0.375;
68  _phi_coeff[7] = 1.0;
69 
70  _vel_coeff[0][6] = 0.25;
71  _vel_coeff[0][7] = 0.0;
72 
73  _vel_coeff[1][6] = 0.25;
74  _vel_coeff[1][7] = 0.5;
75 
76  _vel_coeff[2][0] = 0.025;
77  _vel_coeff[2][1] = 0.0;
78  _vel_coeff[2][2] = 0.125;
79  _vel_coeff[2][3] = 1.0;
80  _vel_coeff[2][4] = 0.25;
81  _vel_coeff[2][5] = 0.25;
82  _vel_coeff[2][6] = 0.25;
83  _vel_coeff[2][7] = 0.25;
84 
85  _T_coeff[6] = 0.125;
86  _T_coeff[7] = 0.5;
87  }
88 }
89 
90 //---------------------------------------------------------------------------//
91 template<class EvalType, class Traits, int NumSpaceDim>
92 void MethodManufacturedSolution<EvalType, Traits, NumSpaceDim>::postRegistrationSetup(
93  typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
94 {
95  _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
96 }
97 
98 //---------------------------------------------------------------------------//
99 template<class EvalType, class Traits, int NumSpaceDim>
100 void MethodManufacturedSolution<EvalType, Traits, NumSpaceDim>::evaluateFields(
101  typename Traits::EvalData workset)
102 {
103  _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
104  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
105  workset.num_cells);
106  Kokkos::parallel_for(this->getName(), policy, *this);
107 }
108 
109 //---------------------------------------------------------------------------//
110 template<class EvalType, class Traits, int NumSpaceDim>
111 KOKKOS_INLINE_FUNCTION void
112 MethodManufacturedSolution<EvalType, Traits, NumSpaceDim>::operator()(
113  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
114 {
115  const int cell = team.league_rank();
116  const int num_point = _lagrange_pressure.extent(1);
117  using Constants::pi;
118  using std::sin;
119 
120  Kokkos::parallel_for(
121  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
122  // function = B + A * sin(2*pi*f_x*(x-phi_x)) *
123  // sin(2*pi*f_y*(y-phi_y))
124  // * sin(2*pi*f_z*(z-phi_z)) for 3D
125  auto set_function
126  = [=](const Kokkos::Array<double, num_coeff> coeff,
127  const Kokkos::Array<double, num_space_dim> x) {
128  const double val
129  = coeff[0]
130  * sin(2.0 * pi * coeff[2] * (x[0] - coeff[3]))
131  * sin(2.0 * pi * coeff[4] * (x[1] - coeff[5]));
132  const double return_val
133  = num_space_dim == 2
134  ? val + coeff[1]
135  : val
136  * sin(2.0 * pi * coeff[6]
137  * (x[2] - coeff[7]))
138  + coeff[1];
139  return return_val;
140  };
141 
142  Kokkos::Array<double, num_space_dim> x;
143  for (int dim = 0; dim < num_space_dim; ++dim)
144  x[dim] = _ip_coords(cell, point, dim);
145 
146  _lagrange_pressure(cell, point) = set_function(_phi_coeff, x);
147  _temperature(cell, point) = set_function(_T_coeff, x);
148  for (int i = 0; i < num_space_dim; ++i)
149  _velocity[i](cell, point) = set_function(_vel_coeff[i], x);
150  });
151 }
152 
153 //---------------------------------------------------------------------------//
154 
155 } // end namespace ClosureModel
156 } // end namespace VertexCFD
157 
158 #endif // VERTEXCFD_CLOSURE_METHODMANUFACTUREDSOLUTION_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23