https://mooseframework.inl.gov
ConservativeAdvection.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "ConservativeAdvection.h"
11 #include "SystemBase.h"
12 
13 using namespace libMesh;
14 
17 
18 template <bool is_ad>
21 {
23  params.addClassDescription("Conservative form of $\\nabla \\cdot \\vec{v} u$ which in its weak "
24  "form is given by: $(-\\nabla \\psi_i, \\vec{v} u)$.");
25  MooseEnum upwinding_type("none full", "none");
26  params.addParam<MooseEnum>("upwinding_type",
27  upwinding_type,
28  "Type of upwinding used. None: Typically results in overshoots and "
29  "undershoots, but numerical diffusion is minimized. Full: Overshoots "
30  "and undershoots are avoided, but numerical diffusion is large");
31  params.addParam<MaterialPropertyName>("advected_quantity",
32  "An optional material property to be advected. If not "
33  "supplied, then the variable will be used.");
34  return params;
35 }
36 
37 template <>
40 {
41  InputParameters params = generalParams();
42  params.addCoupledVar("velocity", "Velocity vector");
43  params.deprecateCoupledVar("velocity", "velocity_variable", "12/31/2025");
44  params.addParam<MaterialPropertyName>("velocity_material", "Velocity vector given as a material");
45  return params;
46 }
47 
48 template <>
51 {
52  InputParameters params = generalParams();
53  params.addCoupledVar("velocity_variable", "Velocity vector given as a variable");
54  params.addParam<MaterialPropertyName>("velocity", "Velocity vector given as a material");
55  params.deprecateParam("velocity", "velocity_material", "12/31/2025");
56  return params;
57 }
58 
59 template <bool is_ad>
61  : GenericKernel<is_ad>(parameters),
62  _velocity(this->isParamValid("velocity_variable")
63  ? &this->template coupledGenericVectorValue<is_ad>("velocity_variable")
64  : (this->isParamValid("velocity_material")
65  ? &this->template getGenericMaterialProperty<RealVectorValue, is_ad>(
66  "velocity_material")
67  .get()
68  : nullptr)),
69  _adv_quant(
70  isParamValid("advected_quantity")
71  ? this->template getGenericMaterialProperty<Real, is_ad>("advected_quantity").get()
72  : _u),
73  _upwinding(
74  this->template getParam<MooseEnum>("upwinding_type").template getEnum<UpwindingType>()),
75  _u_nodal(_var.template genericDofValues<is_ad>()),
76  _upwind_node(0),
77  _dtotal_mass_out(0)
78 {
79  if (_upwinding != UpwindingType::none && this->isParamValid("advected_quantity"))
80  paramError(
81  "advected_quantity",
82  "Upwinding is not compatable with an advected quantity that is not the primary variable.");
83 
84  if (!_velocity ||
85  (this->isParamValid("velocity_variable") && this->isParamValid("velocity_material")))
86  paramError("velocity_variable",
87  "Either velocity_variable or velocity_material must be specified");
88 
89  if (this->_has_diag_save_in)
90  paramError("diag_save_in",
91  "_local_ke not computed for global AD indexing. Save-in is deprecated anyway. Use "
92  "the tagging system instead.");
93 }
94 
95 template <bool is_ad>
98 {
99  return -_grad_test[_i][_qp] * (*_velocity)[_qp];
100 }
101 
102 template <bool is_ad>
105 {
106  // This is the no-upwinded version
107  // It gets called via GenericKernel<is_ad>::computeResidual()
108  return negSpeedQp() * _adv_quant[_qp];
109 }
110 
111 template <>
112 Real
114 {
115  // This is the no-upwinded version
116  // It gets called via GenericKernel<false>::computeJacobian()
117  return negSpeedQp() * _phi[_j][_qp];
118 }
119 
120 template <>
121 Real
123 {
124  mooseError("Internal error, should never get here when using AD");
125  return 0.0;
126 }
127 
128 template <bool is_ad>
129 void
131 {
132  switch (_upwinding)
133  {
134  case UpwindingType::none:
136  break;
137  case UpwindingType::full:
138  fullUpwind(JacRes::CALCULATE_RESIDUAL);
139  break;
140  }
141 }
142 
143 template <bool is_ad>
144 void
146 {
147  switch (_upwinding)
148  {
149  case UpwindingType::none:
151  break;
152  case UpwindingType::full:
153  fullUpwind(JacRes::CALCULATE_JACOBIAN);
154  break;
155  }
156 }
157 
158 template <bool is_ad>
159 void
161 {
162  // The number of nodes in the element
163  const unsigned int num_nodes = _test.size();
164 
165  // Even if we are computing the Jacobian we still need to compute the outflow from each node to
166  // see which nodes are upwind and which are downwind
167  _my_local_re.resize(_var.dofIndices().size());
168 
169  if (!is_ad && (res_or_jac == JacRes::CALCULATE_JACOBIAN))
170  prepareMatrixTag(this->_assembly, _var.number(), _var.number());
171 
172  // Compute the outflux from each node and store in _my_local_re
173  // If _my_local_re is positive at the node, mass (or whatever the Variable represents) is flowing
174  // out of the node
175  _upwind_node.resize(num_nodes);
176  for (_i = 0; _i < num_nodes; ++_i)
177  {
178  for (_qp = 0; _qp < this->_qrule->n_points(); _qp++)
179  _my_local_re(_i) += this->_JxW[_qp] * this->_coord[_qp] * negSpeedQp();
180  _upwind_node[_i] = (MetaPhysicL::raw_value(_my_local_re(_i)) >= 0.0);
181  }
182 
183  // Variables used to ensure mass conservation
184  GenericReal<is_ad> total_mass_out = 0.0;
185  GenericReal<is_ad> total_in = 0.0;
186  if (res_or_jac == JacRes::CALCULATE_JACOBIAN)
187  _dtotal_mass_out.assign(num_nodes, 0.0);
188 
189  for (const auto n : make_range(num_nodes))
190  {
191  if (_upwind_node[n])
192  {
193  if constexpr (!is_ad)
194  if (res_or_jac == JacRes::CALCULATE_JACOBIAN)
195  {
196  if (_test.size() == _phi.size())
197  /* u at node=n depends only on the u at node=n, by construction. For
198  * linear-lagrange variables, this means that Jacobian entries involving the derivative
199  * will only be nonzero for derivatives wrt variable at node=n. Hence the
200  * (n, n) in the line below. The above "if" statement catches other variable types
201  * (eg constant monomials)
202  */
203  _local_ke(n, n) += _my_local_re(n);
204 
205  _dtotal_mass_out[n] += _local_ke(n, n);
206  }
207  _my_local_re(n) *= _u_nodal[n];
208  total_mass_out += _my_local_re(n);
209  }
210  else // downwind node
211  total_in -= _my_local_re(n); // note the -= means the result is positive
212  }
213 
214  // Conserve mass over all phases by proportioning the total_mass_out mass to the inflow nodes,
215  // weighted by their local_re values
216  for (const auto n : make_range(num_nodes))
217  if (!_upwind_node[n]) // downwind node
218  {
219  if constexpr (!is_ad)
220  if (res_or_jac == JacRes::CALCULATE_JACOBIAN)
221  for (_j = 0; _j < _phi.size(); _j++)
222  _local_ke(n, _j) += _my_local_re(n) * _dtotal_mass_out[_j] / total_in;
223  _my_local_re(n) *= total_mass_out / total_in;
224  }
225 
226  // Add the result to the residual and jacobian
227  if (res_or_jac == JacRes::CALCULATE_RESIDUAL)
228  {
229  this->addResiduals(this->_assembly, _my_local_re, _var.dofIndices(), _var.scalingFactor());
230 
231  if (this->_has_save_in)
232  {
233  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
234  for (const auto & var : this->_save_in)
235  var->sys().solution().add_vector(MetaPhysicL::raw_value(_my_local_re), var->dofIndices());
236  }
237  }
238 
239  if (res_or_jac == JacRes::CALCULATE_JACOBIAN)
240  {
241  if constexpr (is_ad)
242  this->addJacobian(this->_assembly, _my_local_re, _var.dofIndices(), _var.scalingFactor());
243  else
244  accumulateTaggedLocalMatrix();
245  }
246 }
247 
248 template class ConservativeAdvectionTempl<false>;
249 template class ConservativeAdvectionTempl<true>;
Moose::GenericType< Real, is_ad > GenericReal
Definition: MooseTypes.h:649
virtual void computeResidual() override
Compute this Kernel&#39;s contribution to the residual.
Definition: Kernel.C:92
virtual void computeResidual() override
Compute this Kernel&#39;s contribution to the residual.
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:435
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:333
static InputParameters validParams()
Definition: GenericKernel.h:19
T * get(const std::unique_ptr< T > &u)
The MooseUtils::get() specializations are used to support making forwards-compatible code changes fro...
Definition: MooseUtils.h:1155
const MooseArray< GenericRealVectorValue< is_ad > > * _velocity
advection velocity
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
enum ConservativeAdvectionTempl::UpwindingType _upwinding
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: KernelBase.h:69
GenericReal< is_ad > negSpeedQp() const
Returns - _grad_test * velocity.
registerMooseObject("MooseApp", ConservativeAdvection)
void deprecateParam(const std::string &old_name, const std::string &new_name, const std::string &removal_date)
static InputParameters generalParams()
virtual GenericReal< is_ad > computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
ConservativeAdvectionTempl(const InputParameters &parameters)
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
Advection of the variable by the velocity provided by the user.
void deprecateCoupledVar(const std::string &old_name, const std::string &new_name, const std::string &removal_date)
void fullUpwind(JacRes res_or_jac)
Calculates the fully-upwind Residual and Jacobian (depending on res_or_jac)
virtual Real computeQpJacobian() override
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
void addCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
virtual void computeJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.
Definition: Kernel.C:115
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
JacRes
enum to make the code clearer
IntRange< T > make_range(T beg, T end)
static InputParameters validParams()
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Definition: MooseBase.h:195
virtual void computeJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.