VertexCFD  0.0-dev
VertexCFD_EquationSet_SolidInductionLessMHD_impl.hpp
1 #ifndef VERTEXCFD_EQUATIONSET_SOLIDINDUCTIONLESSMHD_IMPL_HPP
2 #define VERTEXCFD_EQUATIONSET_SOLIDINDUCTIONLESSMHD_IMPL_HPP
3 
4 #include <Panzer_BasisIRLayout.hpp>
5 #include <Panzer_EvaluatorStyle.hpp>
6 #include <Panzer_IntegrationRule.hpp>
7 #include <Panzer_Integrator_BasisTimesScalar.hpp>
8 #include <Panzer_Integrator_GradBasisDotVector.hpp>
9 
10 #include <stdexcept>
11 #include <vector>
12 
13 namespace VertexCFD
14 {
15 namespace EquationSet
16 {
17 //---------------------------------------------------------------------------//
18 template<class EvalType>
19 SolidInductionLessMHD<EvalType>::SolidInductionLessMHD(
20  const Teuchos::RCP<Teuchos::ParameterList>& params,
21  const int& default_integration_order,
22  const panzer::CellData& cell_data,
23  const Teuchos::RCP<panzer::GlobalData>& global_data,
24  const bool build_transient_support)
25  : panzer::EquationSet_DefaultImpl<EvalType>(params,
26  default_integration_order,
27  cell_data,
28  global_data,
29  build_transient_support)
30  , _dof_name("electric_potential")
31  , _equ_name("electric_potential_equation")
32 {
33  // Set default parameter values and validate the inputs.
34  Teuchos::ParameterList valid_parameters;
35  this->setDefaultValidParameters(valid_parameters);
36  valid_parameters.set(
37  "Model ID", "", "Closure model id associated with this equation set");
38  valid_parameters.set("Basis Order", 1, "Order of the basis");
39  valid_parameters.set("Integration Order", 2, "Order of the integration");
40  valid_parameters.set(
41  "Build Source Residual", false, "Build source residual");
42  params->validateParametersAndSetDefaults(valid_parameters);
43 
44  // Extract parameters.
45  const int basis_order = params->get<int>("Basis Order");
46  const int integration_order
47  = params->get<int>("Integration Order", basis_order + 1);
48  const auto model_id = params->get<std::string>("Model ID");
49  _build_source_term = params->get<bool>("Build Source Residual");
50 
51  // Set up degrees of freedom for solid induction-less MHD equation
52  const std::string residual_name = "RESIDUAL_" + _equ_name;
53  const std::string scatter_name = "SCATTER_" + _equ_name;
54  const std::string basis_type = "HGrad";
55  this->addDOF(_dof_name,
56  basis_type,
57  basis_order,
58  integration_order,
59  residual_name,
60  scatter_name);
61  this->addDOFGrad(_dof_name);
62 
63  // Add closure models and set-up DOFs
64  this->addClosureModel(model_id);
65  this->setupDOFs();
66 }
67 
68 //---------------------------------------------------------------------------//
69 template<class EvalType>
70 void SolidInductionLessMHD<EvalType>::buildAndRegisterEquationSetEvaluators(
71  PHX::FieldManager<panzer::Traits>& fm,
72  const panzer::FieldLibrary&,
73  const Teuchos::ParameterList&) const
74 {
75  // Get integration rule and basis function for field
76  const auto ir = this->getIntRuleForDOF(_dof_name);
77  const auto basis = this->getBasisIRLayoutForDOF(_dof_name);
78 
79  // Add grad-basis time residuals
80  auto add_grad_basis_time_residual = [&, this](
81  const std::string& equ_name,
82  const std::string& residual_name,
83  const double& multiplier,
84  std::vector<std::string>& op_names) {
85  const std::string closure_model_residual = residual_name + "_"
86  + equ_name;
87  const std::string full_residual = "RESIDUAL_" + closure_model_residual;
88  auto op = Teuchos::rcp(
89  new panzer::Integrator_GradBasisDotVector<EvalType, panzer::Traits>(
90  panzer::EvaluatorStyle::EVALUATES,
91  full_residual,
92  closure_model_residual,
93  *basis,
94  *ir,
95  multiplier));
96  this->template registerEvaluator<EvalType>(fm, op);
97  op_names.push_back(full_residual);
98  };
99 
100  // Collect each term for residual to form solid induction-less MHD equation
101  std::vector<std::string> residual_operator_names;
102 
103  // Second-order flux: heat flux
104  add_grad_basis_time_residual(
105  _equ_name, "ELECTRIC_POTENTIAL_FLUX", 1.0, residual_operator_names);
106 
107  // Build total residual for solid induction-less MHD equation
108  this->buildAndRegisterResidualSummationEvaluator(fm,
109  _dof_name,
110  residual_operator_names,
111  "RESIDUAL_electric_"
112  "potential_equation");
113 }
114 
115 //---------------------------------------------------------------------------/
116 template<class EvalType>
117 std::string SolidInductionLessMHD<EvalType>::fieldName(const int dof) const
118 {
119  if (dof > 0)
120  {
121  throw std::logic_error(
122  "SolidInductionLessMHD equation contributes a single DOF");
123  }
124 
125  return _dof_name;
126 }
127 
128 //---------------------------------------------------------------------------//
129 
130 } // end namespace EquationSet
131 } // end namespace VertexCFD
132 
133 #endif // end VERTEXCFD_EQUATIONSET_SOLIDINDUCTIONLESSMHD_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23