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");
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.");
44 params.
addParam<MaterialPropertyName>(
"velocity_material",
"Velocity vector given as a material");
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");
62 _velocity(this->isParamValid(
"velocity_variable")
63 ? &this->template coupledGenericVectorValue<is_ad>(
"velocity_variable")
64 : (this->isParamValid(
"velocity_material")
70 isParamValid(
"advected_quantity")
71 ? this->template getGenericMaterialProperty<
Real, is_ad>(
"advected_quantity").
get()
75 _u_nodal(_var.template genericDofValues<is_ad>()),
82 "Upwinding is not compatable with an advected quantity that is not the primary variable.");
87 "Either velocity_variable or velocity_material must be specified");
91 "_local_ke not computed for global AD indexing. Save-in is deprecated anyway. Use " 92 "the tagging system instead.");
99 return -_grad_test[_i][_qp] * (*_velocity)[_qp];
102 template <
bool is_ad>
108 return negSpeedQp() * _adv_quant[_qp];
117 return negSpeedQp() * _phi[_j][_qp];
124 mooseError(
"Internal error, should never get here when using AD");
128 template <
bool is_ad>
134 case UpwindingType::none:
137 case UpwindingType::full:
138 fullUpwind(JacRes::CALCULATE_RESIDUAL);
143 template <
bool is_ad>
149 case UpwindingType::none:
152 case UpwindingType::full:
153 fullUpwind(JacRes::CALCULATE_JACOBIAN);
158 template <
bool is_ad>
163 const unsigned int num_nodes = _test.size();
167 _my_local_re.resize(_var.dofIndices().size());
169 if (!is_ad && (res_or_jac == JacRes::CALCULATE_JACOBIAN))
170 prepareMatrixTag(this->_assembly, _var.number(), _var.number());
175 _upwind_node.resize(num_nodes);
176 for (_i = 0; _i < num_nodes; ++_i)
178 for (_qp = 0; _qp < this->_qrule->n_points(); _qp++)
179 _my_local_re(_i) += this->_JxW[_qp] * this->_coord[_qp] * negSpeedQp();
186 if (res_or_jac == JacRes::CALCULATE_JACOBIAN)
187 _dtotal_mass_out.assign(num_nodes, 0.0);
193 if constexpr (!is_ad)
194 if (res_or_jac == JacRes::CALCULATE_JACOBIAN)
196 if (_test.size() == _phi.size())
203 _local_ke(n, n) += _my_local_re(n);
205 _dtotal_mass_out[n] += _local_ke(n, n);
207 _my_local_re(n) *= _u_nodal[n];
208 total_mass_out += _my_local_re(n);
211 total_in -= _my_local_re(n);
217 if (!_upwind_node[n])
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;
227 if (res_or_jac == JacRes::CALCULATE_RESIDUAL)
229 this->addResiduals(this->_assembly, _my_local_re, _var.dofIndices(), _var.scalingFactor());
231 if (this->_has_save_in)
233 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
234 for (
const auto & var : this->_save_in)
239 if (res_or_jac == JacRes::CALCULATE_JACOBIAN)
242 this->addJacobian(this->_assembly, _my_local_re, _var.dofIndices(), _var.scalingFactor());
244 accumulateTaggedLocalMatrix();
Moose::GenericType< Real, is_ad > GenericReal
virtual void computeResidual() override
Compute this Kernel's contribution to the residual.
virtual void computeResidual() override
Compute this Kernel's contribution to the residual.
void paramError(const std::string ¶m, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
static InputParameters validParams()
T * get(const std::unique_ptr< T > &u)
The MooseUtils::get() specializations are used to support making forwards-compatible code changes fro...
const MooseArray< GenericRealVectorValue< is_ad > > * _velocity
advection velocity
UpwindingType
Type of upwinding.
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.
GenericReal< is_ad > negSpeedQp() const
Returns - _grad_test * velocity.
registerMooseObject("MooseApp", ConservativeAdvection)
static InputParameters generalParams()
virtual GenericReal< is_ad > computeQpResidual() override
Compute this Kernel's contribution to the residual at the current quadrature point.
ConservativeAdvectionTempl(const InputParameters ¶meters)
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Advection of the variable by the velocity provided by the user.
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's contribution to the Jacobian at the current quadrature point.
virtual void computeJacobian() override
Compute this Kernel's contribution to the diagonal Jacobian entries.
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()
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
virtual void computeJacobian() override
Compute this Kernel's contribution to the diagonal Jacobian entries.