1 #ifndef VERTEXCFD_CLOSURE_METHODMANUFACTUREDSOLUTION_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_METHODMANUFACTUREDSOLUTION_IMPL_HPP
4 #include "utils/VertexCFD_Utils_Constants.hpp"
5 #include "utils/VertexCFD_Utils_VectorField.hpp"
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>
18 namespace ClosureModel
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)
29 this->addEvaluatedField(_lagrange_pressure);
30 Utils::addEvaluatedVectorField(
31 *
this, ir.dl_scalar, _velocity,
"Exact_velocity_");
32 this->addEvaluatedField(_temperature);
34 this->setName(
"Method of Manufactured Solution "
35 + std::to_string(num_space_dim) +
"D");
37 _phi_coeff[0] = 0.0125;
41 _phi_coeff[4] = 0.125;
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;
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;
65 if (num_space_dim == 3)
67 _phi_coeff[6] = 0.375;
70 _vel_coeff[0][6] = 0.25;
71 _vel_coeff[0][7] = 0.0;
73 _vel_coeff[1][6] = 0.25;
74 _vel_coeff[1][7] = 0.5;
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;
91 template<
class EvalType,
class Traits,
int NumSpaceDim>
92 void MethodManufacturedSolution<EvalType, Traits, NumSpaceDim>::postRegistrationSetup(
93 typename Traits::SetupData sd, PHX::FieldManager<Traits>&)
95 _ir_index = panzer::getIntegrationRuleIndex(_ir_degree, (*sd.worksets_)[0]);
99 template<
class EvalType,
class Traits,
int NumSpaceDim>
100 void MethodManufacturedSolution<EvalType, Traits, NumSpaceDim>::evaluateFields(
101 typename Traits::EvalData workset)
103 _ip_coords = workset.int_rules[_ir_index]->ip_coordinates;
104 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
106 Kokkos::parallel_for(this->getName(), policy, *
this);
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
115 const int cell = team.league_rank();
116 const int num_point = _lagrange_pressure.extent(1);
120 Kokkos::parallel_for(
121 Kokkos::TeamThreadRange(team, 0, num_point), [&](
const int point) {
126 = [=](
const Kokkos::Array<double, num_coeff> coeff,
127 const Kokkos::Array<double, num_space_dim> x) {
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
136 * sin(2.0 * pi * coeff[6]
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);
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);
158 #endif // VERTEXCFD_CLOSURE_METHODMANUFACTUREDSOLUTION_IMPL_HPP