www.mooseframework.org
ArrayIntegratedBC.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 "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 
21 
24 {
26  return params;
27 }
28 
30  : IntegratedBCBase(parameters),
32  false,
33  "variable",
36  _var(*mooseVariable()),
37  _normals(_assembly.normals()),
38  _phi(_assembly.phiFace(_var)),
39  _test(_var.phiFace()),
40  _u(_is_implicit ? _var.sln() : _var.slnOld()),
41  _count(_var.count())
42 {
44 
45  _save_in.resize(_save_in_strings.size());
46  _diag_save_in.resize(_diag_save_in_strings.size());
47 
48  for (unsigned int i = 0; i < _save_in_strings.size(); i++)
49  {
51 
52  if (var->feType() != _var.feType())
53  paramError(
54  "save_in",
55  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
57  _save_in[i] = var;
60  }
61 
62  _has_save_in = _save_in.size() > 0;
63 
64  for (unsigned int i = 0; i < _diag_save_in_strings.size(); i++)
65  {
67 
68  if (var->feType() != _var.feType())
69  paramError(
70  "diag_save_in",
71  "saved-in auxiliary variable is incompatible with the object's nonlinear variable: ",
73 
74  _diag_save_in[i] = var;
77  }
78 
79  _has_diag_save_in = _diag_save_in.size() > 0;
80 }
81 
82 void
84 {
86 
87  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
88  {
90  for (_i = 0; _i < _test.size(); _i++)
91  {
93  mooseAssert(residual.size() == _count,
94  "Size of local residual is not equal to the number of array variable compoments");
96  }
97  }
98 
100 
101  if (_has_save_in)
102  {
103  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
104  for (const auto & var : _save_in)
105  {
106  auto * avar = dynamic_cast<ArrayMooseVariable *>(var);
107  if (avar)
108  avar->addSolution(_local_re);
109  else
110  mooseError("Save-in variable for an array kernel must be an array variable");
111  }
112  }
113 }
114 
115 void
117 {
119 
120  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
121  {
122  initQpJacobian();
123  for (_i = 0; _i < _test.size(); _i++)
124  for (_j = 0; _j < _phi.size(); _j++)
125  {
128  _local_ke, _i, _test.size(), _j, _phi.size(), _var.number(), v);
129  }
130  }
131 
133 
134  if (_has_diag_save_in)
135  {
136  DenseVector<Number> diag = _assembly.getJacobianDiagonal(_local_ke);
137  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
138  for (const auto & var : _diag_save_in)
139  {
140  auto * avar = dynamic_cast<ArrayMooseVariable *>(var);
141  if (avar)
142  avar->addSolution(diag);
143  else
144  mooseError("Save-in variable for an array kernel must be an array variable");
145  }
146  }
147 }
148 
149 void
151 {
152  size_t jvar_num = jvar.number();
153  bool same_var = jvar_num == _var.number();
154 
155  prepareMatrixTag(_assembly, _var.number(), jvar_num);
156 
157  // This (undisplaced) jvar could potentially yield the wrong phi size if this object is acting on
158  // the displaced mesh
159  auto phi_size = _sys.getVariable(_tid, jvar_num).dofIndices().size();
160 
161  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
162  {
163  initQpOffDiagJacobian(jvar);
164  for (_i = 0; _i < _test.size(); _i++)
165  for (_j = 0; _j < phi_size; _j++)
166  {
169  _local_ke, _i, _test.size(), _j, jvar.phiSize(), _var.number(), jvar.number(), v);
170  }
171  }
172 
174 
175  if (_has_diag_save_in && same_var)
176  {
177  DenseVector<Number> diag = _assembly.getJacobianDiagonal(_local_ke);
178  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
179  for (const auto & var : _diag_save_in)
180  {
181  auto * avar = dynamic_cast<ArrayMooseVariable *>(var);
182  if (avar)
183  avar->addSolution(diag);
184  else
185  mooseError("Save-in variable for an array kernel must be an array variable");
186  }
187  }
188 }
189 
190 void
192 {
194 
196  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
197  for (_i = 0; _i < _test.size(); _i++)
198  {
201  _local_ke, _i, _test.size(), 0, 1, _var.number(), jvar, v);
202  }
203 
205 }
IntegratedBCBase::_i
unsigned int _i
i-th, j-th index for enumerating test and shape functions
Definition: IntegratedBCBase.h:91
Moose::VarFieldType
VarFieldType
Definition: MooseTypes.h:613
MooseVariableFEBase
Definition: MooseVariableFEBase.h:27
Moose
Definition: Moose.h:116
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
TaggingInterface::accumulateTaggedLocalMatrix
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
Definition: TaggingInterface.C:213
MooseVariableInterface< RealEigenVector >::mooseVariable
MooseVariableFE< RealEigenVector > * mooseVariable() const
Get the variable that this object is using.
Definition: MooseVariableInterface.C:70
SystemBase.h
MooseVariableScalar
Class for scalar variables (they are different).
Definition: MooseVariableScalar.h:30
MooseVariableBase::feType
const FEType & feType() const
Get the type of finite element object.
Definition: MooseVariableBase.h:53
MooseObject::mooseError
void mooseError(Args &&... args) const
Definition: MooseObject.h:141
MooseVariableFEBase::phiSize
virtual size_t phiSize() const =0
Return phi size.
Assembly::saveLocalArrayResidual
void saveLocalArrayResidual(DenseVector< Number > &re, unsigned int i, unsigned int ntest, const RealEigenVector &v)
Helper function for assembling residual contriubutions on local quadrature points for an array variab...
Definition: Assembly.h:1276
ArrayIntegratedBC::computeJacobianBlockScalar
void computeJacobianBlockScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: ArrayIntegratedBC.C:191
TaggingInterface::_local_re
DenseVector< Number > _local_re
Holds residual entries as they are accumulated by this Kernel.
Definition: TaggingInterface.h:156
ArrayIntegratedBC
Base class for deriving any boundary condition of a integrated type.
Definition: ArrayIntegratedBC.h:25
MooseVariableDependencyInterface::addMooseVariableDependency
void addMooseVariableDependency(MooseVariableFEBase *var)
Call this function to add the passed in MooseVariableFEBase as a variable that this object depends on...
Definition: MooseVariableDependencyInterface.h:36
IntegratedBCBase::_has_diag_save_in
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: IntegratedBCBase.h:99
SystemBase::addVariableToZeroOnJacobian
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:183
ArrayIntegratedBC::_phi
const ArrayVariablePhiValue & _phi
shape function values (in QPs)
Definition: ArrayIntegratedBC.h:104
ArrayIntegratedBC::computeJacobianBlock
virtual void computeJacobianBlock(MooseVariableFEBase &jvar) override
Computes d-ivar-residual / d-jvar...
Definition: ArrayIntegratedBC.C:150
ArrayIntegratedBC::computeQpResidual
virtual RealEigenVector computeQpResidual()=0
Method for computing the residual at quadrature points.
IntegratedBCBase::_coord
const MooseArray< Real > & _coord
coordinate transformation
Definition: IntegratedBCBase.h:89
TaggingInterface::prepareMatrixTag
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing element jacobian according to the ative tags.
Definition: TaggingInterface.C:157
Assembly.h
ArrayIntegratedBC::validParams
static InputParameters validParams()
Definition: ArrayIntegratedBC.C:23
SystemBase::getScalarVariable
virtual MooseVariableScalar & getScalarVariable(THREAD_ID tid, const std::string &var_name)
Gets a reference to a scalar variable with specified number.
Definition: SystemBase.C:142
SubProblem::getArrayVariable
virtual ArrayMooseVariable & getArrayVariable(THREAD_ID tid, const std::string &var_name)=0
Returns the variable reference for requested ArrayMooseVariable which may be in any system.
ArrayIntegratedBC::computeJacobian
virtual void computeJacobian() override
Definition: ArrayIntegratedBC.C:116
IntegratedBCBase::_save_in
std::vector< MooseVariableFEBase * > _save_in
Definition: IntegratedBCBase.h:95
MooseObject::paramError
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: MooseObject.h:215
defineLegacyParams
defineLegacyParams(ArrayIntegratedBC)
IntegratedBCBase::_qrule
const QBase *const & _qrule
active quadrature rule
Definition: IntegratedBCBase.h:83
ArrayIntegratedBC::initQpOffDiagJacobian
virtual void initQpOffDiagJacobian(MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
Definition: ArrayIntegratedBC.h:96
InputParameters
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system.
Definition: InputParameters.h:53
ArrayIntegratedBC::computeQpOffDiagJacobianScalar
virtual RealEigenMatrix computeQpOffDiagJacobianScalar(MooseVariableScalar &jvar)
This is the virtual that derived classes should override for computing a full Jacobian component.
Definition: ArrayIntegratedBC.h:77
MooseVariableScalar.h
IntegratedBCBase
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBCBase.h:30
Assembly::saveFullLocalArrayJacobian
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)
Helper function for assembling full Jacobian contriubutions on local quadrature points for an array v...
Definition: Assembly.h:1320
ArrayIntegratedBC::_count
const unsigned int _count
Number of components of the array variable.
Definition: ArrayIntegratedBC.h:113
TaggingInterface::accumulateTaggedLocalResidual
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
Definition: TaggingInterface.C:199
IntegratedBCBase::_diag_save_in_strings
std::vector< AuxVariableName > _diag_save_in_strings
Definition: IntegratedBCBase.h:101
SystemBase::addVariableToZeroOnResidual
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:171
ArrayIntegratedBC::ArrayIntegratedBC
ArrayIntegratedBC(const InputParameters &parameters)
Definition: ArrayIntegratedBC.C:29
ArrayIntegratedBC::initQpResidual
virtual void initQpResidual()
Put necessary evaluations depending on qp but independent on test functions here.
Definition: ArrayIntegratedBC.h:85
IntegratedBCBase::validParams
static InputParameters validParams()
Definition: IntegratedBCBase.C:21
TaggingInterface::prepareVectorTag
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual the according to active tags.
Definition: TaggingInterface.C:121
ArrayIntegratedBC::computeQpOffDiagJacobian
virtual RealEigenMatrix computeQpOffDiagJacobian(MooseVariableFEBase &jvar)
Method for computing an off-diagonal jacobian component at quadrature points.
Definition: ArrayIntegratedBC.h:60
moose::internal::incompatVarMsg
std::string incompatVarMsg(MooseVariableFEBase &var1, MooseVariableFEBase &var2)
Builds and returns a string of the form:
Definition: MooseError.C:23
BoundaryCondition::_assembly
Assembly & _assembly
Reference to assembly.
Definition: BoundaryCondition.h:109
BoundaryCondition::_sys
SystemBase & _sys
Reference to SystemBase.
Definition: BoundaryCondition.h:103
Moose::VAR_FIELD_ARRAY
Definition: MooseTypes.h:618
Assembly::getJacobianDiagonal
DenseVector< Real > getJacobianDiagonal(DenseMatrix< Number > &ke)
Definition: Assembly.h:1339
ArrayIntegratedBC::computeQpJacobian
virtual RealEigenVector computeQpJacobian()
Method for computing the diagonal Jacobian at quadrature points.
Definition: ArrayIntegratedBC.h:55
MooseVariableFE.h
ArrayIntegratedBC.h
libMesh::RealEigenMatrix
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:136
SubProblem.h
ArrayIntegratedBC::initQpJacobian
virtual void initQpJacobian()
Put necessary evaluations depending on qp but independent on test and shape functions here.
Definition: ArrayIntegratedBC.h:90
Moose::VAR_NONLINEAR
Definition: MooseTypes.h:608
ArrayIntegratedBC::_var
ArrayMooseVariable & _var
Definition: ArrayIntegratedBC.h:98
Assembly::saveDiagLocalArrayJacobian
void saveDiagLocalArrayJacobian(DenseMatrix< Number > &ke, unsigned int i, unsigned int ntest, unsigned int j, unsigned int nphi, unsigned int ivar, const RealEigenVector &v)
Helper function for assembling diagonal Jacobian contriubutions on local quadrature points for an arr...
Definition: Assembly.h:1295
MooseVariableBase::dofIndices
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
Definition: MooseVariableBase.h:114
IntegratedBCBase::_qp
unsigned int _qp
quadrature point index
Definition: IntegratedBCBase.h:81
IntegratedBCBase::_JxW
const MooseArray< Real > & _JxW
transformed Jacobian weights
Definition: IntegratedBCBase.h:87
Moose::VarKindType
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:606
MooseVariableBase::sys
SystemBase & sys()
Get the system this variable is part of.
Definition: MooseVariableBase.h:58
IntegratedBCBase::_save_in_strings
std::vector< AuxVariableName > _save_in_strings
Definition: IntegratedBCBase.h:96
IntegratedBCBase::_j
unsigned int _j
Definition: IntegratedBCBase.h:91
TaggingInterface::_local_ke
DenseMatrix< Number > _local_ke
Holds residual entries as they are accumulated by this Kernel.
Definition: TaggingInterface.h:159
ArrayIntegratedBC::computeResidual
virtual void computeResidual() override
Definition: ArrayIntegratedBC.C:83
IntegratedBCBase::_diag_save_in
std::vector< MooseVariableFEBase * > _diag_save_in
Definition: IntegratedBCBase.h:100
BoundaryCondition::_tid
THREAD_ID _tid
Thread id.
Definition: BoundaryCondition.h:106
MooseArray::size
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
IntegratedBCBase::_has_save_in
bool _has_save_in
The aux variables to save the residual contributions to.
Definition: IntegratedBCBase.h:94
libMesh::RealEigenVector
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:133
MooseVariableInterface
Interface for objects that need to get values of MooseVariables.
Definition: MooseVariableInterface.h:24
ArrayIntegratedBC::_test
const ArrayVariableTestValue & _test
test function values (in QPs)
Definition: ArrayIntegratedBC.h:107
MooseVariableFE< RealEigenVector >
BoundaryCondition::_subproblem
SubProblem & _subproblem
Reference to SubProblem.
Definition: BoundaryCondition.h:97
MooseVariableBase::number
unsigned int number() const
Get variable number coming from libMesh.
Definition: MooseVariableBase.h:48