Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
ArrayDGLowerDKernel.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 "ArrayDGLowerDKernel.h"
11 
12 #include "Assembly.h"
13 #include "MooseVariable.h"
14 #include "Problem.h"
15 #include "SubProblem.h"
16 #include "SystemBase.h"
17 #include "MaterialData.h"
18 #include "ParallelUniqueId.h"
19 
20 #include "libmesh/dof_map.h"
21 #include "libmesh/dense_vector.h"
22 #include "libmesh/numeric_vector.h"
23 #include "libmesh/dense_subvector.h"
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/quadrature.h"
26 
29 {
31  params.addRequiredCoupledVar("lowerd_variable", "Lagrange multiplier");
32  return params;
33 }
34 
36  : ArrayDGKernel(parameters),
37  _lowerd_var(*getArrayVar("lowerd_variable", 0)),
38  _lambda(_is_implicit ? _lowerd_var.slnLower() : _lowerd_var.slnLowerOld()),
39 
40  _phi_lambda(_lowerd_var.phiLower()),
41  _test_lambda(_lowerd_var.phiLower()),
42  _work_vector(_count)
43 {
44  const auto & lower_domains = _lowerd_var.activeSubdomains();
45  for (const auto & id : _mesh.interiorLowerDBlocks())
46  if (lower_domains.count(id) == 0)
48  "moose",
49  29151,
50  "'lowerd_variable' must be defined on the interior lower-dimensional subdomain '" +
52  "' that is added by Mesh/build_all_side_lowerd_mesh=true.\nThe check could be overly "
53  "restrictive.");
54 
55  for (const auto & id : _var.activeSubdomains())
56  if (_mesh.interiorLowerDBlocks().count(id) > 0)
57  paramError("variable",
58  "Must not be defined on the interior lower-dimensional subdomain '" +
59  _mesh.getSubdomainName(id) + "'");
60 
61  if (_lowerd_var.count() != _count)
62  paramError("lowerd_variable",
63  "The number of components must be equal to the number of "
64  "components of 'variable'");
65 }
66 
67 void
69 {
70  if (!excludeBoundary())
71  {
73 
74  // Compute the residual for this element
76 
77  // Compute the residual for the neighbor
79 
81  }
82 }
83 
84 void
86 {
88 
89  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
90  {
92  for (_i = 0; _i < _test_lambda.size(); _i++)
93  {
94  _work_vector.setZero();
96  mooseAssert(_work_vector.size() == _count,
97  "Size of local residual is not equal to the number of array variable compoments");
100  }
101  }
102 
104 }
105 
106 void
108 {
109  if (!excludeBoundary())
110  {
112 
113  // Compute element-element Jacobian
115 
116  // Compute element-neighbor Jacobian
118 
119  // Compute neighbor-element Jacobian
121 
122  // Compute neighbor-neighbor Jacobian
124 
125  // Compute the other five pieces of Jacobian related with lower-d variable
127  }
128 }
129 
130 void
132 {
135  "Jacobian types without lower should be handled in computeElemNeighJacobian");
136 
137  const ArrayVariableTestValue & test_space =
139  ? _test_lambda
141  unsigned int ivar =
143  ? _lowerd_var.number()
144  : _var.number();
145 
146  const ArrayVariableTestValue & loc_phi =
148  ? _phi_lambda
150  unsigned int jvar =
152  ? _lowerd_var.number()
153  : _var.number();
154 
155  prepareMatrixTagLower(_assembly, ivar, jvar, type);
156 
157  if (_local_ke.n() == 0 || _local_ke.m() == 0)
158  return;
159 
160  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
161  {
163  for (_i = 0; _i < test_space.size(); _i++)
164  for (_j = 0; _j < loc_phi.size(); _j++)
165  {
166  // vector or matrix?
169  _local_ke, _i, test_space.size(), _j, loc_phi.size(), ivar, v);
170  }
171  }
172 
174 }
175 
176 void
178 {
179  if (!excludeBoundary())
180  {
181  precalculateOffDiagJacobian(jvar_num);
182 
183  /*
184  * Note: we cannot call compute jacobian functions like in DGLowerDKernel
185  * because we could have cross component couplings for the array variables
186  */
187 
188  const auto & jvar = getVariable(jvar_num);
189 
190  // Compute element-element Jacobian
192 
193  // Compute element-neighbor Jacobian
195 
196  // Compute neighbor-element Jacobian
198 
199  // Compute neighbor-neighbor Jacobian
201 
202  // Compute the other five pieces of Jacobian related with lower-d variable
208  }
209 }
210 
211 void
213  const MooseVariableFEBase & jvar)
214 {
217  "Jacobian types without lower should be handled in computeOffDiagElemNeighJacobian");
218 
219  const ArrayVariableTestValue & test_space =
221  ? _test_lambda
223  unsigned int ivar =
225  ? _lowerd_var.number()
226  : _var.number();
227 
228  prepareMatrixTagLower(_assembly, ivar, jvar.number(), type);
229  if (_local_ke.n() == 0 || _local_ke.m() == 0)
230  return;
231 
233  {
234  const auto & jv0 = static_cast<const MooseVariable &>(jvar);
235  const VariableTestValue & loc_phi =
237  ? jv0.phiLower()
238  : (type == Moose::LowerPrimary ? jv0.phiFace() : jv0.phiFaceNeighbor());
239 
240  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
241  {
243  for (_i = 0; _i < test_space.size(); _i++)
244  for (_j = 0; _j < loc_phi.size(); _j++)
245  {
248  _local_ke, _i, test_space.size(), _j, loc_phi.size(), ivar, jvar.number(), v);
249  }
250  }
251  }
253  {
254  const auto & jv1 = static_cast<const ArrayMooseVariable &>(jvar);
255  const ArrayVariableTestValue & loc_phi =
257  ? jv1.phiLower()
258  : (type == Moose::LowerPrimary ? jv1.phiFace() : jv1.phiFaceNeighbor());
259 
260  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
261  {
263  for (_i = 0; _i < test_space.size(); _i++)
264  for (_j = 0; _j < loc_phi.size(); _j++)
265  {
268  _local_ke, _i, test_space.size(), _j, loc_phi.size(), ivar, jvar.number(), v);
269  }
270  }
271  }
272  else
273  mooseError("Vector variable cannot be coupled into array DG kernel currently");
274 
276 }
277 
280  const MooseVariableFEBase & jvar)
281 {
282  if (jvar.number() == _var.number())
283  {
285  return computeLowerDQpJacobian(type).asDiagonal();
286  }
287  else if (jvar.number() == _lowerd_var.number())
288  {
290  return computeLowerDQpJacobian(type).asDiagonal();
291  }
292 
293  return RealEigenMatrix::Zero(_var.count(), jvar.count());
294 }
MooseMesh & _mesh
Reference to this Kernel&#39;s mesh object.
const QBase *const & _qrule
Quadrature rule.
Definition: DGKernelBase.h:107
virtual void computeLowerDJacobian(Moose::ConstraintJacobianType jacobian_type)
Computes one of the five pieces of Jacobian involving lower-d.
const std::set< SubdomainID > & interiorLowerDBlocks() const
Definition: MooseMesh.h:1388
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
virtual void computeLowerDResidual()
Computes the Lower part of residual for the variable on the lower-d element.
virtual void computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, const MooseVariableFEBase &jvar) override
Computes the element-element off-diagonal Jacobian.
virtual void precalculateOffDiagJacobian(unsigned int)
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:1795
unsigned int count() const
Get the number of components Note: For standard and vector variables, the number is one...
virtual void precalculateResidual()
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
RealEigenVector _work_vector
Work vector for residual computation.
unsigned int m() const
This class provides an interface for common operations on field variables of both FE and FV types wit...
unsigned int _i
Definition: DGKernelBase.h:136
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
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:1776
virtual void computeJacobian() override
Computes the nine pieces of element/neighbor/lower-d - element/neighbor/lower-d Jacobian.
const ArrayVariableTestValue & _test_lambda
test functions
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:1732
ArrayDGLowerDKernel(const InputParameters &parameters)
virtual void computeElemNeighJacobian(Moose::DGJacobianType type) override
Computes the element/neighbor-element/neighbor Jacobian.
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
Definition: ArrayDGKernel.C:27
virtual void computeElemNeighResidual(Moose::DGResidualType type) override
Computes the residual for this element or the neighbor.
const ArrayVariableTestValue & _test
test functions
virtual void initLowerDQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
unsigned int _qp
Definition: DGKernelBase.h:134
void mooseDocumentedError(const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
Emits a documented error with object name and type.
const unsigned int _count
Number of components of the array variable.
unsigned int _j
Definition: DGKernelBase.h:136
virtual void precalculateJacobian()
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:320
bool excludeBoundary() const
Check current element if it contains broken boundary.
Definition: DGKernelBase.C:187
const std::string & type() const
Get the type of this class.
Definition: MooseBase.h:51
const MooseArray< Real > & _coord
Coordinate transform mainly for curvilinear coordinates.
Definition: DGKernelBase.h:111
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 ...
Eigen::Matrix< Real, Eigen::Dynamic, Eigen::Dynamic > RealEigenMatrix
Definition: MooseTypes.h:145
virtual void computeLowerDQpResidual(RealEigenVector &residual)=0
Method for computing the Lower part of residual at quadrature points, to be filled in residual...
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
Definition: ArrayDGKernel.h:20
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.
const std::set< SubdomainID > & activeSubdomains() const
The subdomains the variable is active on.
virtual void initLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
virtual void computeResidual() override
Computes the residual for this element, the neighbor and the lower-d element.
virtual void computeOffDiagLowerDJacobian(Moose::ConstraintJacobianType type, const MooseVariableFEBase &jvar)
Computes one of the five pieces of off-diagonal Jacobian involving lower-d.
virtual RealEigenVector computeLowerDQpJacobian(Moose::ConstraintJacobianType jacobian_type)=0
Computes one of the five pieces of Jacobian involving lower-d at quadrature points.
ConstraintJacobianType
Definition: MooseTypes.h:790
const ArrayMooseVariable & _lowerd_var
Variable this kernel operates on.
virtual void initLowerDQpJacobian(Moose::ConstraintJacobianType)
Put necessary evaluations depending on qp but independent on test and shape functions here...
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 MooseArray< Real > & _JxW
Jacobian det times quadrature weighting on quadrature points.
Definition: DGKernelBase.h:109
virtual Moose::VarFieldType fieldType() const =0
Filed type of this variable.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
const ArrayVariablePhiValue & _phi_neighbor
Side shape function.
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:142
void prepareVectorTagLower(Assembly &assembly, unsigned int ivar)
Prepare data for computing the residual according to active tags for mortar constraints.
unsigned int n() const
const ArrayVariablePhiValue & _phi_lambda
Shape functions.
OutputTools< RealEigenVector >::VariableTestValue ArrayVariableTestValue
Definition: MooseTypes.h:355
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:1820
const ArrayVariableTestValue & _test_neighbor
Side test function.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar...
virtual RealEigenMatrix computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType type, const MooseVariableFEBase &jvar)
Computes one of the five pieces of off-diagonal Jacobian involving lower-d at quadrature points...
ArrayMooseVariable & _var
Variable this kernel operates on.
Definition: ArrayDGKernel.h:96
const ArrayVariablePhiValue & _phi
Shape functions.