https://mooseframework.inl.gov
ArrayLowerDIntegratedBC.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 
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  params.addRequiredCoupledVar("lowerd_variable", "Lagrange multiplier");
25  return params;
26 }
27 
29  : ArrayIntegratedBC(parameters),
30  _lowerd_var(*getArrayVar("lowerd_variable", 0)),
31  _lambda(_is_implicit ? _lowerd_var.slnLower() : _lowerd_var.slnLowerOld()),
32 
33  _phi_lambda(_lowerd_var.phiLower()),
34  _test_lambda(_lowerd_var.phiLower()),
35  _work_vector(_count)
36 {
37  const auto & lower_domains = _lowerd_var.activeSubdomains();
38  for (const auto & id : _mesh.boundaryLowerDBlocks())
39  if (lower_domains.count(id) == 0)
41  "moose",
42  29151,
43  "'lowerd_variable' must be defined on the boundary lower-dimensional subdomain '" +
45  "' that is added by Mesh/build_all_side_lowerd_mesh=true.\nThe check could be overly "
46  "restrictive.");
47 
48  for (const auto & id : _var.activeSubdomains())
49  if (_mesh.boundaryLowerDBlocks().count(id) > 0)
50  paramError("variable",
51  "Must not be defined on the boundary lower-dimensional subdomain '" +
52  _mesh.getSubdomainName(id) + "'");
53 
54  if (_lowerd_var.count() != _count)
55  paramError("lowerd_variable",
56  "The number of components must be equal to the number of "
57  "components of 'variable'");
58 }
59 
60 void
62 {
64 
66 
67  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
68  {
70  for (_i = 0; _i < _test_lambda.size(); _i++)
71  {
72  _work_vector.setZero();
74  mooseAssert(_work_vector.size() == _count,
75  "Size of local residual is not equal to the number of array variable compoments");
78  }
79  }
80 
82 }
83 
84 void
86 {
88 
90 }
91 
92 void
94 {
95  mooseAssert(type == Moose::LowerLower || type == Moose::LowerPrimary ||
97  "Jacobian types must have lower in computeLowerDJacobian");
98 
99  const ArrayVariableTestValue & test_space =
101  unsigned int ivar = (type == Moose::LowerLower || type == Moose::LowerPrimary)
102  ? _lowerd_var.number()
103  : _var.number();
104 
105  const ArrayVariableTestValue & loc_phi =
107  unsigned int jvar = (type == Moose::LowerLower || type == Moose::PrimaryLower)
108  ? _lowerd_var.number()
109  : _var.number();
110 
111  // need to transform the type for assembling Jacobian on boundary to be consistent with
112  // Assembly::addJacobianLowerD() and Assembly::prepareLowerD().
113  Moose::ConstraintJacobianType type_transformed =
115  ? type
117  prepareMatrixTagLower(_assembly, ivar, jvar, type_transformed);
118 
119  if (_local_ke.n() == 0 || _local_ke.m() == 0)
120  return;
121 
122  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
123  {
125  for (_i = 0; _i < test_space.size(); _i++)
126  for (_j = 0; _j < loc_phi.size(); _j++)
127  {
130  _local_ke, _i, test_space.size(), _j, loc_phi.size(), ivar, v);
131  }
132  }
133 
135 }
136 
137 void
139 {
141 
145 }
146 
147 void
149  const unsigned int jvar_num)
150 {
151  mooseAssert(type == Moose::LowerLower || type == Moose::LowerPrimary ||
153  "Jacobian types must have lower in computeLowerDJacobian");
154 
155  const ArrayVariableTestValue & test_space =
157  unsigned int ivar = (type == Moose::LowerLower || type == Moose::LowerPrimary)
158  ? _lowerd_var.number()
159  : _var.number();
160 
161  const auto & jvar = getVariable(jvar_num);
162 
163  // need to transform the type for assembling Jacobian on boundary to be consistent with
164  // Assembly::addJacobianLowerD() and Assembly::prepareLowerD().
165  Moose::ConstraintJacobianType type_transformed =
167  ? type
169  prepareMatrixTagLower(_assembly, ivar, jvar_num, type_transformed);
170  if (_local_ke.n() == 0 || _local_ke.m() == 0)
171  return;
172 
173  if (jvar.fieldType() == Moose::VarFieldType::VAR_FIELD_STANDARD)
174  {
175  const auto & jv0 = static_cast<const MooseVariable &>(jvar);
176  const VariableTestValue & loc_phi =
177  (type == Moose::LowerLower || type == Moose::PrimaryLower) ? jv0.phiLower() : jv0.phiFace();
178 
179  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
180  {
182  for (_i = 0; _i < test_space.size(); _i++)
183  for (_j = 0; _j < loc_phi.size(); _j++)
184  {
187  _local_ke, _i, test_space.size(), _j, loc_phi.size(), ivar, jvar_num, v);
188  }
189  }
190  }
191  else if (jvar.fieldType() == Moose::VarFieldType::VAR_FIELD_ARRAY)
192  {
193  const auto & jv1 = static_cast<const ArrayMooseVariable &>(jvar);
194  const ArrayVariableTestValue & loc_phi =
195  (type == Moose::LowerLower || type == Moose::PrimaryLower) ? jv1.phiLower() : jv1.phiFace();
196 
197  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
198  {
200  for (_i = 0; _i < test_space.size(); _i++)
201  for (_j = 0; _j < loc_phi.size(); _j++)
202  {
205  _local_ke, _i, test_space.size(), _j, loc_phi.size(), ivar, jvar_num, v);
206  }
207  }
208  }
209  else
210  mooseError("Vector variable cannot be coupled into array DG kernel currently");
211 
213 }
214 
217  const MooseVariableFEBase & jvar)
218 {
219  if (jvar.number() == _var.number())
220  {
221  if (type == Moose::LowerPrimary)
222  return computeLowerDQpJacobian(type).asDiagonal();
223  }
224  else if (jvar.number() == _lowerd_var.number())
225  {
227  return computeLowerDQpJacobian(type).asDiagonal();
228  }
229 
230  return RealEigenMatrix::Zero(_count, jvar.count());
231 }
MooseMesh & _mesh
Reference to this Kernel&#39;s mesh object.
ArrayLowerDIntegratedBC(const InputParameters &parameters)
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
const ArrayVariablePhiValue & _phi
shape function values (in QPs)
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:435
virtual void computeLowerDJacobian(Moose::ConstraintJacobianType type)
Method for computing the LowerLower, PrimaryLower and LowerPrimary parts of Jacobian.
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
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:1796
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...
static InputParameters validParams()
unsigned int _i
i-th, j-th index for enumerating test and shape functions
void mooseDocumentedError(const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
Definition: MooseBase.h:273
unsigned int m() const
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:1777
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
const ArrayVariableTestValue & _test_lambda
test functions
const MooseVariableFieldBase & getVariable(unsigned int jvar_num) const
Retrieve the variable object from our system associated with jvar_num.
const std::string & getSubdomainName(SubdomainID subdomain_id) const
Return the name of a block given an id.
Definition: MooseMesh.C:1763
unsigned int _qp
quadrature point index
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
const ArrayMooseVariable & _lowerd_var
Variable this kernel operates on.
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:325
ArrayMooseVariable & _var
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-ivar-residual / d-jvar...
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:89
virtual void initLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:149
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
const MooseArray< Real > & _coord
coordinate transformation
const std::set< SubdomainID > & boundaryLowerDBlocks() const
Definition: MooseMesh.h:1407
void computeLowerDOffDiagJacobian(Moose::ConstraintJacobianType type, const unsigned int jvar_num)
Method for computing an off-diagonal jacobian component.
const ArrayVariablePhiValue & _phi_lambda
Shape functions.
virtual void initLowerDQpJacobian(Moose::ConstraintJacobianType)
Put necessary evaluations depending on qp but independent on test and shape functions here...
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-ivar-residual / d-jvar...
const std::set< SubdomainID > & activeSubdomains() const
The subdomains the variable is active on.
const QBase *const & _qrule
active quadrature rule
forward declarations
const ArrayVariableTestValue & _test
test function values (in QPs)
ConstraintJacobianType
Definition: MooseTypes.h:797
virtual const FieldVariablePhiValue & phiLower() const override
Return the variable&#39;s shape functions on a lower-dimensional element.
virtual RealEigenMatrix computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType type, const MooseVariableFEBase &jvar)
Method for computing an off-diagonal jacobian component at quadrature points.
void prepareMatrixTagLower(Assembly &assembly, unsigned int ivar, unsigned int jvar, Moose::ConstraintJacobianType type)
Prepare data for computing the jacobian according to the active tags for mortar.
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
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:267
virtual RealEigenVector computeLowerDQpJacobian(Moose::ConstraintJacobianType)=0
Method for computing the LowerLower, PrimaryLower and LowerPrimary parts of Jacobian at quadrature po...
Base class for deriving any boundary condition of a integrated type.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:146
void prepareVectorTagLower(Assembly &assembly, unsigned int ivar)
Prepare data for computing the residual according to active tags for mortar constraints.
virtual void computeLowerDQpResidual(RealEigenVector &residual)=0
Method for computing the Lower part of residual at quadrature points, to be filled in residual...
unsigned int n() const
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:1821
virtual void initLowerDQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
RealEigenVector _work_vector
Work vector for residual.
const MooseArray< Real > & _JxW
transformed Jacobian weights
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.