www.mooseframework.org
KernelGrad.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "KernelGrad.h"
11 #include "Assembly.h"
12 #include "SystemBase.h"
13 #include "libmesh/quadrature.h"
14 
16 
19 {
20  return Kernel::validParams();
21 }
22 
23 KernelGrad::KernelGrad(const InputParameters & parameters) : Kernel(parameters) {}
24 
25 void
27 {
28  DenseVector<Number> & re = _assembly.residualBlock(_var.number());
29  _local_re.resize(re.size());
30  _local_re.zero();
31 
32  const unsigned int n_test = _test.size();
33  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
34  {
35  RealGradient value = precomputeQpResidual() * _JxW[_qp] * _coord[_qp];
36  for (_i = 0; _i < n_test; _i++) // target for auto vectorization
38  }
39 
40  re += _local_re;
41 
42  if (_has_save_in)
43  {
44  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
45  for (const auto & var : _save_in)
46  var->sys().solution().add_vector(_local_re, var->dofIndices());
47  }
48 }
49 
50 void
52 {
53  DenseMatrix<Number> & ke = _assembly.jacobianBlock(_var.number(), _var.number());
54  _local_ke.resize(ke.m(), ke.n());
55  _local_ke.zero();
56 
57  const unsigned int n_test = _test.size();
58  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
59  for (_j = 0; _j < _phi.size(); _j++)
60  {
61  RealGradient value = precomputeQpJacobian() * _JxW[_qp] * _coord[_qp];
62  for (_i = 0; _i < n_test; _i++) // target for auto vectorization
63  _local_ke(_i, _j) += value * _grad_test[_i][_qp];
64  }
65 
66  ke += _local_ke;
67 
69  {
70  const unsigned int rows = ke.m();
71  DenseVector<Number> diag(rows);
72  for (unsigned int i = 0; i < rows; i++) // target for auto vectorization
73  diag(i) = _local_ke(i, i);
74 
75  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
76  for (const auto & var : _diag_save_in)
77  var->sys().solution().add_vector(diag, var->dofIndices());
78  }
79 }
80 
81 void
83 {
84  size_t jvar_num = jvar.number();
85  if (jvar_num == _var.number())
87  else
88  {
89  DenseMatrix<Number> & Ke = _assembly.jacobianBlock(_var.number(), jvar_num);
90 
91  // This (undisplaced) jvar could potentially yield the wrong phi size if this object is acting
92  // on the displaced mesh
93  auto phi_size = _sys.getVariable(_tid, jvar.number()).dofIndices().size();
94 
95  for (_j = 0; _j < phi_size; _j++)
96  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
97  for (_i = 0; _i < _test.size(); _i++)
98  Ke(_i, _j) += _JxW[_qp] * _coord[_qp] * computeQpOffDiagJacobian(jvar_num);
99  }
100 }
101 
102 Real
104 {
105  return 0.0;
106 }
107 
108 RealGradient
110 {
111  return RealGradient(0.0);
112 }
KernelBase::_i
unsigned int _i
current index for the test function
Definition: KernelBase.h:159
MooseVariableFEBase
Definition: MooseVariableFEBase.h:27
MooseVariableInterface< Real >::value
virtual const OutputTools< Real >::VariableValue & value()
The value of the variable this object is operating on.
Definition: MooseVariableInterface.C:77
SystemBase::getVariable
MooseVariableFEBase & getVariable(THREAD_ID tid, const std::string &var_name)
Gets a reference to a variable of with specified name.
Definition: SystemBase.C:109
KernelBase::_diag_save_in
std::vector< MooseVariableFEBase * > _diag_save_in
Definition: KernelBase.h:171
Kernel::_phi
const VariablePhiValue & _phi
the current shape functions
Definition: Kernel.h:84
KernelBase::_assembly
Assembly & _assembly
Reference to this Kernel's assembly object.
Definition: KernelBase.h:133
SystemBase.h
KernelGrad::precomputeQpJacobian
virtual RealGradient precomputeQpJacobian()
Called before forming the jacobian for an element.
Definition: KernelGrad.C:109
TaggingInterface::_local_re
DenseVector< Number > _local_re
Holds residual entries as they are accumulated by this Kernel.
Definition: TaggingInterface.h:156
KernelBase::_j
unsigned int _j
current index for the shape function
Definition: KernelBase.h:162
Assembly::residualBlock
DenseVector< Number > & residualBlock(unsigned int var_num, TagID tag_id=0)
Get local residual block for a variable and a tag.
Definition: Assembly.h:817
Assembly.h
Kernel::validParams
static InputParameters validParams()
Definition: Kernel.C:25
Kernel::computeQpOffDiagJacobian
virtual Real computeQpOffDiagJacobian(unsigned int)
This is the virtual that derived classes should override for computing an off-diagonal Jacobian compo...
Definition: Kernel.h:64
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
KernelBase::_qp
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:144
KernelGrad::KernelGrad
KernelGrad(const InputParameters &parameters)
Factory constructor initializes all internal references needed for residual computation.
Definition: KernelGrad.C:23
KernelGrad::computeOffDiagJacobian
virtual void computeOffDiagJacobian(MooseVariableFEBase &jvar) override
Computes d-residual / d-jvar... storing the result in Ke.
Definition: KernelGrad.C:82
KernelGrad
The KernelGrad class is responsible for calculating the residuals in form:
Definition: KernelGrad.h:27
KernelBase::_save_in
std::vector< MooseVariableFEBase * > _save_in
Definition: KernelBase.h:166
KernelGrad::computeResidual
virtual void computeResidual() override
Compute this Kernel's contribution to the residual.
Definition: KernelGrad.C:26
KernelBase::_tid
THREAD_ID _tid
The thread ID for this kernel.
Definition: KernelBase.h:130
KernelBase::_JxW
const MooseArray< Real > & _JxW
The current quadrature point weight value.
Definition: KernelBase.h:153
Kernel::_test
const VariableTestValue & _test
the current test function
Definition: Kernel.h:78
defineLegacyParams
defineLegacyParams(KernelGrad)
Kernel::_grad_test
const VariableTestGradient & _grad_test
gradient of the test function
Definition: Kernel.h:81
Kernel
Definition: Kernel.h:20
KernelBase::_sys
SystemBase & _sys
Reference to the EquationSystem object.
Definition: KernelBase.h:127
KernelBase::_has_diag_save_in
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: KernelBase.h:170
KernelGrad::computeQpResidual
virtual Real computeQpResidual() override
Compute this Kernel's contribution to the residual at the current quadrature point.
Definition: KernelGrad.C:103
KernelGrad::computeJacobian
virtual void computeJacobian() override
Compute this Kernel's contribution to the diagonal Jacobian entries.
Definition: KernelGrad.C:51
TaggingInterface::_local_ke
DenseMatrix< Number > _local_ke
Holds residual entries as they are accumulated by this Kernel.
Definition: TaggingInterface.h:159
KernelBase::_qrule
const QBase *const & _qrule
active quadrature rule
Definition: KernelBase.h:150
KernelGrad.h
KernelGrad::validParams
static InputParameters validParams()
Definition: KernelGrad.C:18
Kernel::_var
MooseVariable & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: Kernel.h:75
Assembly::jacobianBlock
DenseMatrix< Number > & jacobianBlock(unsigned int ivar, unsigned int jvar, TagID tag=0)
Get local Jacobian block for a pair of variables and a tag.
Definition: Assembly.h:841
KernelGrad::precomputeQpResidual
virtual RealGradient precomputeQpResidual()=0
Called before forming the residual for an element.
KernelBase::_has_save_in
bool _has_save_in
The aux variables to save the residual contributions to.
Definition: KernelBase.h:165
KernelBase::_coord
const MooseArray< Real > & _coord
The scaling factor to convert from cartesian to another coordinate system (e.g rz,...
Definition: KernelBase.h:156
MooseVariableBase::number
unsigned int number() const
Get variable number coming from libMesh.
Definition: MooseVariableBase.h:48