Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
KernelScalarBase.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 "Kernel.h"
13 
23 class KernelScalarBase : public Kernel
24 {
25 public:
27 
29 
34  {
35  mooseAssert(_kappa_var_ptr, "kappa pointer should have been set in the constructor");
36  return *_kappa_var_ptr;
37  }
38 
39  virtual void computeResidual() override;
40  virtual void computeJacobian() override;
44  virtual void computeOffDiagJacobian(unsigned int jvar_num) override;
49  void computeOffDiagJacobianScalar(unsigned int svar_num) override;
50 
51 protected:
57  virtual void precalculateOffDiagJacobianScalar(unsigned int /* svar_num */) {}
58 
62  virtual void computeScalarResidual();
63 
67  virtual Real computeScalarQpResidual();
68 
72  virtual void computeScalarJacobian();
73 
78  virtual Real computeScalarQpJacobian() { return 0; }
79 
83  virtual void computeScalarOffDiagJacobian(const unsigned int jvar_num);
84 
88  virtual Real computeScalarQpOffDiagJacobian(const unsigned int /*jvar_num*/) { return 0; }
89 
94  virtual void computeOffDiagJacobianScalarLocal(const unsigned int svar_num);
95 
99  virtual void computeScalarOffDiagJacobianScalar(const unsigned int svar_num);
100 
104  virtual Real computeScalarQpOffDiagJacobianScalar(const unsigned int /*svar_num*/) { return 0; }
105 
109  virtual void initScalarQpResidual() {}
110 
114  virtual void initScalarQpJacobian(const unsigned int /*svar_num*/) {}
115 
121 
123  const bool _use_scalar;
124 
127 
130 
133 
135  const unsigned int _kappa_var;
136 
138  const unsigned int _k_order;
139 
142 
144  unsigned int _h;
145  unsigned int _l;
146 };
147 
148 inline Real
150 {
151  mooseError(
152  "A scalar_variable has been set and compute_scalar_residuals=true, ",
153  "but the computeScalarQpResidual method was not overridden. Accidental call of base class?");
154  return 0;
155 }
virtual void computeScalarResidual()
Method for computing the scalar part of residual.
virtual void computeOffDiagJacobianScalarLocal(const unsigned int svar_num)
Method for computing an off-diagonal jacobian component d-_var-residual / d-scalar Revised version of...
unsigned int _h
Used internally to iterate over each scalar component.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const bool _compute_scalar_residuals
Whether to compute scalar contributions for this instance.
This class provides an interface for common operations on field variables of both FE and FV types wit...
const MooseVariableScalar & scalarVariable() const
The scalar variable that this kernel operates on.
const unsigned int _k_order
Order of the scalar variable, used in several places.
void computeOffDiagJacobianScalar(unsigned int svar_num) override
Computes jacobian block with respect to a scalar variable.
virtual void precalculateOffDiagJacobianScalar(unsigned int)
Precalculate method that is executed prior to scalar coupled variable loop within computeOffDiagJacob...
virtual void computeOffDiagJacobian(unsigned int jvar_num) override
Computes d-_var-residual / d-jvar as well as d-_kappa-residual / d-jvar.
virtual void computeScalarOffDiagJacobianScalar(const unsigned int svar_num)
Method for computing an off-diagonal jacobian component d-_kappa-residual / d-scalar.
virtual void computeScalarJacobian()
Method for computing the scalar variable part of Jacobian.
virtual Real computeScalarQpJacobian()
Method for computing the scalar variable part of Jacobian at quadrature points.
const VariableValue & _kappa
Reference to the current solution at the current quadrature point.
virtual void computeJacobian() override
Compute this Kernel's contribution to the diagonal Jacobian entries.
virtual Real computeScalarQpOffDiagJacobian(const unsigned int)
Method for computing an off-diagonal jacobian component at quadrature points.
KernelScalarBase(const InputParameters &parameters)
virtual void initScalarQpJacobian(const unsigned int)
Put necessary evaluations depending on qp but independent of test and shape functions here...
virtual void computeScalarOffDiagJacobian(const unsigned int jvar_num)
Method for computing an off-diagonal jacobian component d-_kappa-residual / d-jvar.
This Kernel adds standardized methods for assembling to a primary scalar variable associated with the...
virtual Real computeScalarQpOffDiagJacobianScalar(const unsigned int)
Method for computing an off-diagonal jacobian component at quadrature points.
forward declarations
virtual Real computeScalarQpResidual()
Method for computing the scalar part of residual at quadrature points.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static InputParameters validParams()
Class for scalar variables (they are different).
Definition: Kernel.h:15
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
virtual void initScalarQpResidual()
Put necessary evaluations depending on qp but independent of test functions here. ...
const InputParameters & parameters() const
Get the parameters of the object.
const MooseVariableScalar *const _kappa_var_ptr
(Pointer to) Scalar variable this kernel operates on
const unsigned int _kappa_var
The unknown scalar variable ID.
const bool _compute_field_residuals
Whether to compute field contributions for this instance.
virtual void computeResidual() override
Compute this Kernel's contribution to the residual.
const bool _use_scalar
Whether a scalar variable is declared for this kernel.
virtual void initScalarQpOffDiagJacobian(const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent of test and shape functions here for off-di...