https://mooseframework.inl.gov
LowerDIntegratedBC.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 "LowerDIntegratedBC.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  params.addRequiredCoupledVar("lowerd_variable", "Lagrange multiplier");
25  return params;
26 }
27 
29  : IntegratedBC(parameters),
30  _lowerd_var(*getVar("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 {
36  const auto & lower_domains = _lowerd_var.activeSubdomains();
37  for (const auto & id : _mesh.boundaryLowerDBlocks())
38  if (lower_domains.count(id) == 0)
40  "moose",
41  29151,
42  "'lowerd_variable' must be defined on the boundary lower-dimensional subdomain '" +
44  "' that is added by Mesh/build_all_side_lowerd_mesh=true.\nThe check could be overly "
45  "restrictive.");
46 
47  for (const auto & id : _var.activeSubdomains())
48  if (_mesh.boundaryLowerDBlocks().count(id) > 0)
49  paramError("variable",
50  "Must not be defined on the boundary lower-dimensional subdomain '" +
51  _mesh.getSubdomainName(id) + "'");
52 
53  // Note: the above two conditions also ensure that the variable and lower-d variable are
54  // different.
55 }
56 
57 void
59 {
61 
63 
64  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
65  {
67  for (_i = 0; _i < _test_lambda.size(); _i++)
69  }
70 
72 }
73 
74 void
76 {
78 
80 }
81 
82 void
84 {
85  mooseAssert(type == Moose::LowerLower || type == Moose::LowerPrimary ||
87  "Jacobian types must have lower in computeLowerDJacobian");
88 
89  const auto & test_space =
91  unsigned int ivar = (type == Moose::LowerLower || type == Moose::LowerPrimary)
93  : _var.number();
94 
95  const auto & loc_phi =
97  const auto jvar = (type == Moose::LowerLower || type == Moose::PrimaryLower)
99  : _var.number();
100 
101  // need to transform the type for assembling Jacobian on boundary to be consistent with
102  // Assembly::addJacobianLowerD() and Assembly::prepareLowerD().
105  ? type
107  prepareMatrixTagLower(_assembly, ivar, jvar, type_tr);
108 
109  if (_local_ke.n() == 0 || _local_ke.m() == 0)
110  return;
111 
112  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
113  {
115  for (_i = 0; _i < test_space.size(); _i++)
116  for (_j = 0; _j < loc_phi.size(); _j++)
118  }
119 
121 }
122 
123 void
124 LowerDIntegratedBC::computeOffDiagJacobian(const unsigned int jvar_num)
125 {
126  if (jvar_num == variable().number())
127  {
130  }
131  else if (jvar_num == _lowerd_var.number())
132  {
135  }
136  else
137  {
142  }
143 }
144 
145 void
147  const unsigned int jvar_num)
148 {
149  mooseAssert(type == Moose::LowerLower || type == Moose::LowerPrimary ||
151  "Jacobian types must have lower in computeLowerDJacobian");
152 
153  const auto & test_space =
155  const auto ivar = (type == Moose::LowerLower || type == Moose::LowerPrimary)
156  ? _lowerd_var.number()
157  : _var.number();
158 
159  const auto & jvar = getVariable(jvar_num);
160 
161  prepareMatrixTagLower(_assembly, ivar, jvar_num, type);
162  if (_local_ke.n() == 0 || _local_ke.m() == 0)
163  return;
164 
165  if (jvar.fieldType() == Moose::VarFieldType::VAR_FIELD_STANDARD)
166  {
167  const auto & jv0 = static_cast<const MooseVariable &>(jvar);
168  const auto & loc_phi =
169  (type == Moose::LowerLower || type == Moose::PrimaryLower) ? jv0.phiLower() : jv0.phiFace();
170 
171  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
172  {
174  for (_i = 0; _i < test_space.size(); _i++)
175  for (_j = 0; _j < loc_phi.size(); _j++)
177  }
178  }
179  else if (jvar.fieldType() == Moose::VarFieldType::VAR_FIELD_ARRAY)
180  mooseError("Array variable cannot be coupled into integrated BC currently");
181  else
182  mooseError("Vector variable cannot be coupled into integrated BC currently");
183 
185 }
MooseMesh & _mesh
Reference to this Kernel&#39;s mesh object.
const VariableTestValue & _test
test function values (in QPs)
Definition: IntegratedBC.h:97
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-ivar-residual / d-jvar...
Definition: IntegratedBC.C:139
virtual void initLowerDQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
virtual Real computeLowerDQpResidual()=0
Method for computing the Lower part of residual at quadrature points.
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
unsigned int number() const
Get variable number coming from libMesh.
static InputParameters validParams()
Definition: IntegratedBC.C:23
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void initLowerDQpJacobian(Moose::ConstraintJacobianType)
Put necessary evaluations depending on qp but independent on test and shape functions here...
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
const VariableTestValue & _test_lambda
test functions
virtual const MooseVariable & variable() const override
Returns the variable that this object operates on.
const VariablePhiValue & _phi
shape function values (in QPs)
Definition: IntegratedBC.h:90
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.
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 computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: IntegratedBC.C:109
LowerDIntegratedBC(const InputParameters &parameters)
virtual void computeLowerDJacobian(Moose::ConstraintJacobianType type)
Method for computing the LowerLower, PrimaryLower and LowerPrimary parts of Jacobian.
virtual Real computeLowerDQpJacobian(Moose::ConstraintJacobianType)=0
Method for computing the LowerLower, PrimaryLower and LowerPrimary parts of Jacobian at quadrature po...
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
virtual Real computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Method for computing an off-diagonal jacobian component at quadrature points.
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:89
Base class for deriving any boundary condition of a integrated type.
Definition: IntegratedBC.h:18
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-ivar-residual / d-jvar...
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
static InputParameters validParams()
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.
MooseVariable & _var
Definition: IntegratedBC.h:82
const std::set< SubdomainID > & activeSubdomains() const
The subdomains the variable is active on.
void computeLowerDOffDiagJacobian(Moose::ConstraintJacobianType type, const unsigned int jvar_num)
Method for computing an off-diagonal jacobian component.
const QBase *const & _qrule
active quadrature rule
const VariablePhiValue & _phi_lambda
Shape functions.
ConstraintJacobianType
Definition: MooseTypes.h:797
virtual const FieldVariablePhiValue & phiLower() const override
Return the variable&#39;s shape functions on a lower-dimensional element.
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.
const MooseVariable & _lowerd_var
Variable this kernel operates on.
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 void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: IntegratedBC.C:85
virtual void initLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
void prepareVectorTagLower(Assembly &assembly, unsigned int ivar)
Prepare data for computing the residual according to active tags for mortar constraints.
unsigned int n() const
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
const MooseArray< Real > & _JxW
transformed Jacobian weights