https://mooseframework.inl.gov
ADKernelScalarBase.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 "ADKernelScalarBase.h"
11 
12 #include "Assembly.h"
13 #include "SubProblem.h"
14 #include "SystemBase.h"
15 #include "MooseVariableFE.h"
16 #include "MooseVariableScalar.h"
17 #include "ADUtils.h"
18 
19 #include "libmesh/quadrature.h"
20 
23 {
25  // This parameter can get renamed in derived class to a more relevant variable name
26  params.addCoupledVar("scalar_variable", "Primary coupled scalar variable");
27  params.addParam<bool>("compute_scalar_residuals", true, "Whether to compute scalar residuals");
28  params.addParam<bool>(
29  "compute_field_residuals", true, "Whether to compute residuals for the field variable.");
30  return params;
31 }
32 
34  : ADKernel(parameters),
35  _use_scalar(isParamValid("scalar_variable") ? true : false),
36  _compute_scalar_residuals(!_use_scalar ? false : getParam<bool>("compute_scalar_residuals")),
37  _compute_field_residuals(getParam<bool>("compute_field_residuals")),
38  _kappa_var_ptr(_use_scalar ? getScalarVar("scalar_variable", 0) : nullptr),
39  _kappa_var(_use_scalar ? _kappa_var_ptr->number() : 0),
40  _k_order(_use_scalar ? _kappa_var_ptr->order() : 0),
41  _kappa(_use_scalar ? _kappa_var_ptr->adSln() : _ad_zero)
42 {
43 }
44 
45 void
47 {
49  ADKernel::computeResidual(); // compute and assemble regular variable contributions
50 
52  {
53  std::vector<Real> scalar_residuals(_k_order);
54  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
55  {
57  for (_h = 0; _h < _k_order; _h++)
58  scalar_residuals[_h] += _JxW[_qp] * _coord[_qp] * raw_value(computeScalarQpResidual());
59  }
62  }
63 }
64 
65 void
67 {
70 
72  {
78  }
79 }
80 
81 void
82 ADKernelScalarBase::computeOffDiagJacobian(const unsigned int jvar_num)
83 {
84  // Only need to do this once because AD does all the derivatives at once
85  if (jvar_num == _var.number())
87 }
88 
89 void
90 ADKernelScalarBase::computeOffDiagJacobianScalar(const unsigned int /*jvar_num*/)
91 {
92 }
93 
94 void
96 {
99 
101  {
107  }
108 }
109 
110 void
112 {
113  if (_scalar_residuals.size() != _k_order)
114  _scalar_residuals.resize(_k_order, 0);
115  for (auto & sr : _scalar_residuals)
116  sr = 0;
117 
118  // precalculateResidual was already run for the field variable
119  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
120  for (_h = 0; _h < _k_order; _h++)
122 }
void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: ADKernel.C:110
void addResidualsAndJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the provided d...
const MooseArray< Real > & _JxW
The current quadrature point weight value.
Definition: KernelBase.h:52
unsigned int number() const
Get variable number coming from libMesh.
std::vector< ADReal > _scalar_residuals
virtual void initScalarQpResidual()
Put necessary evaluations depending on qp but independent of test functions here. ...
void addResiduals(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided incoming residuals corresponding to the provided dof indices.
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...
void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: ADKernel.C:162
unsigned int _h
Used internally to iterate over each scalar component.
const unsigned int _k_order
Order of the scalar variable, used in several places.
static InputParameters validParams()
void computeResidualAndJacobian() override
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
Definition: ADKernel.C:205
const MooseVariableScalar *const _kappa_var_ptr
(Pointer to) Scalar variable this kernel operates on
virtual void computeScalarResidualsForJacobian()
compute the _scalar_residuals member for filling the Jacobian.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-_var-residual / d-jvar as well as d-_kappa-residual / d-jvar.
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
void computeOffDiagJacobianScalar(unsigned int) override
Computes jacobian block with respect to a scalar variable.
const QBase *const & _qrule
active quadrature rule
Definition: KernelBase.h:49
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
static InputParameters validParams()
Definition: ADKernel.C:24
ADKernelScalarBase(const InputParameters &parameters)
const MooseArray< Real > & _coord
The scaling factor to convert from cartesian to another coordinate system (e.g rz, spherical, etc.)
Definition: KernelBase.h:55
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
void addCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
void computeResidualAndJacobian() override
Computes residual and jacobian block for field and scalar variables.
const bool _compute_field_residuals
Whether to compute field contributions for this instance.
MooseVariableFE< T > & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: ADKernel.h:80
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...
const bool _compute_scalar_residuals
Whether to compute scalar contributions for this instance.
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
virtual ADReal computeScalarQpResidual()
Method for computing the scalar part of residual at quadrature points.
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:43