https://mooseframework.inl.gov
ConservativeAdvection.h
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 #pragma once
11 
12 #include "GenericKernel.h"
13 #include "libmesh/dense_vector.h"
14 
19 template <bool is_ad>
21 {
22 public:
25 
27 
28 protected:
29  virtual GenericReal<is_ad> computeQpResidual() override;
30  virtual Real computeQpJacobian() override;
31  virtual void computeResidual() override;
32  virtual void computeJacobian() override;
33 
36 
39 
41  enum class JacRes
42  {
45  };
46 
48  const enum class UpwindingType { none, full } _upwinding;
49 
52 
54  std::vector<bool> _upwind_node;
55 
57  std::vector<GenericReal<is_ad>> _dtotal_mass_out;
58 
61 
63  void fullUpwind(JacRes res_or_jac);
64 
66 
67 private:
70 };
71 
Moose::GenericType< Real, is_ad > GenericReal
Definition: MooseTypes.h:648
virtual void computeResidual() override
Compute this Kernel&#39;s contribution to the residual.
const MooseArray< GenericRealVectorValue< is_ad > > * _velocity
advection velocity
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
enum ConservativeAdvectionTempl::UpwindingType _upwinding
GenericReal< is_ad > negSpeedQp() const
Returns - _grad_test * velocity.
std::vector< GenericReal< is_ad > > _dtotal_mass_out
In the full-upwind scheme d(total_mass_out)/d(variable_at_node_i)
const MooseArray< GenericReal< is_ad > > & _adv_quant
advected quantity
ConservativeAdvectionTempl< true > ADConservativeAdvection
Moose::GenericType< VariableValue, is_ad > GenericVariableValue
Definition: MooseTypes.h:662
static InputParameters generalParams()
virtual GenericReal< is_ad > computeQpResidual() override
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
std::vector< bool > _upwind_node
In the full-upwind scheme, whether a node is an upwind node.
ConservativeAdvectionTempl(const InputParameters &parameters)
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&#39;s contribution to the Jacobian at the current quadrature point.
forward declarations
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
JacRes
enum to make the code clearer
ConservativeAdvectionTempl< false > ConservativeAdvection
static InputParameters validParams()
const InputParameters & parameters() const
Get the parameters of the object.
libMesh::DenseVector< GenericReal< is_ad > > _my_local_re
A container for holding the local residuals.
const GenericVariableValue< is_ad > & _u_nodal
Nodal value of u, used for full upwinding.
virtual void computeJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.