VertexCFD  0.0-dev
VertexCFD_InitialCondition_MethodManufacturedSolution_impl.hpp
1 #ifndef VERTEXCFD_INITIALCONDITION_METHODMANUFACTUREDSOLUTION_IMPL_HPP
2 #define VERTEXCFD_INITIALCONDITION_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 
12 #include <cmath>
13 #include <string>
14 
15 namespace VertexCFD
16 {
17 namespace InitialCondition
18 {
19 //---------------------------------------------------------------------------//
20 template<class EvalType, class Traits, int NumSpaceDim>
21 MethodManufacturedSolution<EvalType, Traits, NumSpaceDim>::MethodManufacturedSolution(
22  const panzer::PureBasis& basis)
23  : _lagrange_pressure("lagrange_pressure", basis.functional)
24  , _temperature("temperature", basis.functional)
25  , _basis_name(basis.name())
26 {
27  this->addEvaluatedField(_lagrange_pressure);
28  this->addUnsharedField(_lagrange_pressure.fieldTag().clone());
29  this->addEvaluatedField(_temperature);
30  this->addUnsharedField(_temperature.fieldTag().clone());
31 
32  Utils::addEvaluatedVectorField(
33  *this, basis.functional, _velocity, "velocity_", true);
34 
35  this->setName("MethodManufacturedSolution Initial Condition"
36  + std::to_string(num_space_dim) + "D");
37 
38  _phi_coeff[0] = 0.0125;
39  _phi_coeff[1] = 1.0;
40  _phi_coeff[2] = 0.25;
41  _phi_coeff[3] = 0.5;
42  _phi_coeff[4] = 0.125;
43  _phi_coeff[5] = 0.0;
44 
45  _vel_coeff[0][0] = 0.0125;
46  _vel_coeff[0][1] = 0.08;
47  _vel_coeff[0][2] = 0.125;
48  _vel_coeff[0][3] = 0.0;
49  _vel_coeff[0][4] = 0.125;
50  _vel_coeff[0][5] = 0.25;
51 
52  _vel_coeff[1][0] = 0.0375;
53  _vel_coeff[1][1] = 1.125;
54  _vel_coeff[1][2] = 0.25;
55  _vel_coeff[1][3] = 0.0;
56  _vel_coeff[1][4] = 0.375;
57  _vel_coeff[1][5] = 0.5;
58 
59  _T_coeff[0] = 0.0625;
60  _T_coeff[1] = 1.0;
61  _T_coeff[2] = 0.375;
62  _T_coeff[3] = 0.25;
63  _T_coeff[4] = 0.25;
64  _T_coeff[5] = 0.5;
65 
66  if (num_space_dim == 3)
67  {
68  _phi_coeff[6] = 0.375;
69  _phi_coeff[7] = 1.0;
70 
71  _vel_coeff[0][6] = 0.25;
72  _vel_coeff[0][7] = 0.0;
73 
74  _vel_coeff[1][6] = 0.25;
75  _vel_coeff[1][7] = 0.5;
76 
77  _vel_coeff[2][0] = 0.025;
78  _vel_coeff[2][1] = 0.0;
79  _vel_coeff[2][2] = 0.125;
80  _vel_coeff[2][3] = 1.0;
81  _vel_coeff[2][4] = 0.25;
82  _vel_coeff[2][5] = 0.25;
83  _vel_coeff[2][6] = 0.25;
84  _vel_coeff[2][7] = 0.25;
85 
86  _T_coeff[6] = 0.125;
87  _T_coeff[7] = 0.5;
88  }
89 }
90 
91 //---------------------------------------------------------------------------//
92 template<class EvalType, class Traits, int NumSpaceDim>
93 void MethodManufacturedSolution<EvalType, Traits, NumSpaceDim>::postRegistrationSetup(
94  typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
95 {
96  _basis_index = panzer::getPureBasisIndex(
97  _basis_name, (*sd.worksets_)[0], this->wda);
98 }
99 
100 //---------------------------------------------------------------------------//
101 template<class EvalType, class Traits, int NumSpaceDim>
102 void MethodManufacturedSolution<EvalType, Traits, NumSpaceDim>::evaluateFields(
103  typename Traits::EvalData workset)
104 {
105  _basis_coords = this->wda(workset).bases[_basis_index]->basis_coordinates;
106  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
107  workset.num_cells);
108  Kokkos::parallel_for(this->getName(), policy, *this);
109 }
110 
111 //---------------------------------------------------------------------------//
112 template<class EvalType, class Traits, int NumSpaceDim>
113 KOKKOS_INLINE_FUNCTION void
114 MethodManufacturedSolution<EvalType, Traits, NumSpaceDim>::operator()(
115  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
116 {
117  const int cell = team.league_rank();
118  const int num_basis = _lagrange_pressure.extent(1);
119  using Constants::pi;
120  using std::sin;
121 
122  Kokkos::parallel_for(
123  Kokkos::TeamThreadRange(team, 0, num_basis), [&](const int basis) {
124  // function = B + A * sin(2*pi*f_x*(x-phi_x)) *
125  // sin(2*pi*f_y*(x-phi_y))
126  // * sin(2*pi*f_z*(z-phi_z)) for 3D
127  auto set_function
128  = [=](const Kokkos::Array<double, num_coeff> coeff,
129  const Kokkos::Array<double, num_space_dim> x) {
130  double return_val = coeff[0];
131  for (int i = 0; i < num_space_dim; ++i)
132  return_val *= sin(2.0 * pi * coeff[2 * (i + 1)]
133  * (x[i] - coeff[2 * (i + 1) + 1]));
134  return_val += coeff[1];
135  return return_val;
136  };
137 
138  // Note this block is all of `auto` type rather than `scalar_type`
139  // because the values coming from `_basis_coords` are never going
140  // to be FAD objects. By using `auto` instead of `double` we don't
141  // impose a precision on the mesh positions.
142 
143  // FIXME: Kokkos::Array<auto, num_space_dim> is okay??
144  Kokkos::Array<double, num_space_dim> x;
145  for (int dim = 0; dim < num_space_dim; ++dim)
146  x[dim] = _basis_coords(cell, basis, dim);
147 
148  _lagrange_pressure(cell, basis) = set_function(_phi_coeff, x);
149  _temperature(cell, basis) = set_function(_T_coeff, x);
150  for (int i = 0; i < num_space_dim; ++i)
151  _velocity[i](cell, basis) = set_function(_vel_coeff[i], x);
152  });
153 }
154 
155 //---------------------------------------------------------------------------//
156 
157 } // end namespace InitialCondition
158 } // end namespace VertexCFD
159 
160 #endif // end VERTEXCFD_INITIALCONDITION_METHODMANUFACTUREDSOLUTION_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23