https://mooseframework.inl.gov
ArrayIntegratedBC.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 "ArrayIntegratedBC.h"
11 
12 #include "Assembly.h"
13 #include "SubProblem.h"
14 #include "SystemBase.h"
15 #include "MooseVariableFE.h"
16 #include "MooseVariableScalar.h"
17 
18 #include "libmesh/quadrature.h"
19 
22 {
24  return params;
25 }
26 
28  : IntegratedBCBase(parameters),
30  false,
31  "variable",
34  _var(*mooseVariable()),
35  _normals(_assembly.normals()),
36  _phi(_assembly.phiFace(_var)),
37  _test(_var.phiFace()),
38  _u(_is_implicit ? _var.sln() : _var.slnOld()),
39  _count(_var.count()),
40  _work_vector(_count)
41 {
43 
44  _save_in.resize(_save_in_strings.size());
45  _diag_save_in.resize(_diag_save_in_strings.size());
46 
47  for (unsigned int i = 0; i < _save_in_strings.size(); i++)
48  {
50 
51  if (var->feType() != _var.feType())
52  paramError(
53  "save_in",
54  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
56  _save_in[i] = var;
59  }
60 
61  _has_save_in = _save_in.size() > 0;
62 
63  for (unsigned int i = 0; i < _diag_save_in_strings.size(); i++)
64  {
66 
67  if (var->feType() != _var.feType())
68  paramError(
69  "diag_save_in",
70  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
72 
73  _diag_save_in[i] = var;
76  }
77 
78  _has_diag_save_in = _diag_save_in.size() > 0;
79 }
80 
81 void
83 {
85 
86  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
87  {
89  for (_i = 0; _i < _test.size(); _i++)
90  {
91  _work_vector.setZero();
93  mooseAssert(_work_vector.size() == _count,
94  "Size of local residual is not equal to the number of array variable compoments");
97  }
98  }
99 
101 
102  if (_has_save_in)
103  for (const auto & var : _save_in)
104  {
105  auto * avar = dynamic_cast<ArrayMooseVariable *>(var);
106  if (avar)
107  avar->addSolution(_local_re);
108  else
109  mooseError("Save-in variable for an array kernel must be an array variable");
110  }
111 }
112 
113 void
115 {
117 
118  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
119  {
120  initQpJacobian();
121  for (_i = 0; _i < _test.size(); _i++)
122  for (_j = 0; _j < _phi.size(); _j++)
123  {
126  _local_ke, _i, _test.size(), _j, _phi.size(), _var.number(), v);
127  }
128  }
129 
131 
132  if (_has_diag_save_in)
133  {
135  for (const auto & var : _diag_save_in)
136  {
137  auto * avar = dynamic_cast<ArrayMooseVariable *>(var);
138  if (avar)
139  avar->addSolution(diag);
140  else
141  mooseError("Save-in variable for an array kernel must be an array variable");
142  }
143  }
144 }
145 
148 {
149  return RealEigenVector::Zero(_var.count());
150 }
151 
152 void
153 ArrayIntegratedBC::computeOffDiagJacobian(const unsigned int jvar_num)
154 {
155  const auto & jvar = getVariable(jvar_num);
156  if (!jvar.dofIndices().size())
157  // If we have no dof indices then data like phiSize() can't be trusted. For instance a variable
158  // of constant monomial finite element type that lives only on lower-dimensional subdomains will
159  // have zero dof indices here (assuming this BC is being applied on a higher-D face) but will
160  // have a phiSize() of 1 corresponding to the number of shapes for that finite element type
161  // on the higher-dimensional element.
162  return;
163 
164  const auto phi_size = jvar.phiSize();
165 
166  bool same_var = jvar_num == _var.number();
167 
168  prepareMatrixTag(_assembly, _var.number(), jvar_num);
169 
170  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
171  {
172  initQpOffDiagJacobian(jvar);
173  for (_i = 0; _i < _test.size(); _i++)
174  for (_j = 0; _j < phi_size; _j++)
175  {
178  _local_ke, _i, _test.size(), _j, phi_size, _var.number(), jvar_num, v);
179  }
180  }
181 
183 
184  if (_has_diag_save_in && same_var)
185  {
187  for (const auto & var : _diag_save_in)
188  {
189  auto * avar = dynamic_cast<ArrayMooseVariable *>(var);
190  if (avar)
191  avar->addSolution(diag);
192  else
193  mooseError("Save-in variable for an array kernel must be an array variable");
194  }
195  }
196 }
197 
200 {
201  if (jvar.number() == _var.number())
202  return computeQpJacobian().asDiagonal();
203  else
204  return RealEigenMatrix::Zero(_var.count(), jvar.count());
205 }
206 
207 void
209 {
211 
212  const MooseVariableScalar & jv = _sys.getScalarVariable(_tid, jvar);
213  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
214  for (_i = 0; _i < _test.size(); _i++)
215  {
218  _local_ke, _i, _test.size(), 0, 1, _var.number(), jvar, v);
219  }
220 
222 }
223 
226 {
227  return RealEigenMatrix::Zero(_var.count(), (unsigned int)jvar.order() + 1);
228 }
VarFieldType
Definition: MooseTypes.h:770
const libMesh::FEType & feType() const
Get the type of finite element object.
RealEigenVector _work_vector
Work vector for residual.
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
std::vector< MooseVariableFEBase * > _diag_save_in
Class for stuff related to variables.
Definition: Adaptivity.h:31
const ArrayVariablePhiValue & _phi
shape function values (in QPs)
virtual void computeQpResidual(RealEigenVector &residual)=0
Method for computing the residual at quadrature points, to be filled in residual. ...
void paramError(const std::string &param, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
Definition: MooseBase.h:439
std::string incompatVarMsg(MooseVariableFieldBase &var1, MooseVariableFieldBase &var2)
Builds and returns a string of the form:
Definition: MooseError.C:26
std::vector< MooseVariableFEBase * > _save_in
unsigned int number() const
Get variable number coming from libMesh.
void saveDiagLocalArrayJacobian(DenseMatrix< Number > &ke, unsigned int i, unsigned int ntest, unsigned int j, unsigned int nphi, unsigned int ivar, const RealEigenVector &v) const
Helper function for assembling diagonal Jacobian contriubutions on local quadrature points for an arr...
Definition: Assembly.h:1844
void addSolution(const DenseVector< libMesh::Number > &v)
Add passed in local DOF values onto the current solution.
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one...
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
unsigned int _i
i-th, j-th index for enumerating test and shape functions
This class provides an interface for common operations on field variables of both FE and FV types wit...
const unsigned int _count
Number of components of the array variable.
static InputParameters validParams()
void saveLocalArrayResidual(DenseVector< Number > &re, unsigned int i, unsigned int ntest, const RealEigenVector &v) const
Helper function for assembling residual contriubutions on local quadrature points for an array kernel...
Definition: Assembly.h:1816
virtual RealEigenMatrix computeQpOffDiagJacobianScalar(const MooseVariableScalar &jvar)
This is the virtual that derived classes should override for computing a full Jacobian component...
THREAD_ID _tid
The thread ID for this kernel.
DenseVector< Real > getJacobianDiagonal(DenseMatrix< Number > &ke)
Definition: Assembly.h:1895
ArrayIntegratedBC(const InputParameters &parameters)
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
const MooseVariableFieldBase & getVariable(unsigned int jvar_num) const
Retrieve the variable object from our system associated with jvar_num.
unsigned int _qp
quadrature point index
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
virtual RealEigenMatrix computeQpOffDiagJacobian(const MooseVariableFEBase &jvar)
Method for computing an off-diagonal jacobian component at quadrature points.
SystemBase & _sys
Reference to the EquationSystem object.
MooseVariableFE< RealEigenVector > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
std::vector< AuxVariableName > _save_in_strings
SubProblem & _subproblem
Reference to this kernel&#39;s SubProblem.
virtual void initQpJacobian()
Put necessary evaluations depending on qp but independent on test and shape functions here...
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:763
ArrayMooseVariable & _var
virtual void initQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
virtual ArrayMooseVariable & getArrayVariable(const THREAD_ID tid, const std::string &var_name)=0
Returns the variable reference for requested ArrayMooseVariable which may be in any system...
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:151
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
const MooseArray< Real > & _coord
coordinate transformation
bool _has_save_in
The aux variables to save the residual contributions to.
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-ivar-residual / d-jvar...
libMesh::Order order() const
Get the order of this variable Note: Order enum can be implicitly converted to unsigned int...
virtual MooseVariableScalar & getScalarVariable(THREAD_ID tid, const std::string &var_name) const
Gets a reference to a scalar variable with specified number.
Definition: SystemBase.C:145
const QBase *const & _qrule
active quadrature rule
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
virtual void initQpOffDiagJacobian(const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
std::vector< AuxVariableName > _diag_save_in_strings
virtual void addVariableToZeroOnResidual(std::string var_name)
Adds this variable to the list of variables to be zeroed during each residual evaluation.
Definition: SystemBase.C:174
const ArrayVariableTestValue & _test
test function values (in QPs)
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
static InputParameters validParams()
Interface for objects that need to get values of MooseVariables.
Class for scalar variables (they are different).
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type and optionally a file path to the top-level block p...
Definition: MooseBase.h:271
Base class for deriving any boundary condition of a integrated type.
virtual void addVariableToZeroOnJacobian(std::string var_name)
Adds this variable to the list of variables to be zeroed during each Jacobian evaluation.
Definition: SystemBase.C:180
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:147
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual according to active tags.
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing element jacobian according to the active tags.
void saveFullLocalArrayJacobian(DenseMatrix< Number > &ke, unsigned int i, unsigned int ntest, unsigned int j, unsigned int nphi, unsigned int ivar, unsigned int jvar, const RealEigenMatrix &v) const
Helper function for assembling full Jacobian contriubutions on local quadrature points for an array k...
Definition: Assembly.h:1869
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
SystemBase & sys()
Get the system this variable is part of.
void ErrorVector unsigned int
const MooseArray< Real > & _JxW
transformed Jacobian weights
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
virtual RealEigenVector computeQpJacobian()
Method for computing the diagonal Jacobian at quadrature points.