Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
DGLowerDKernel.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 "DGLowerDKernel.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  : DGKernel(parameters),
37  _lowerd_var(*getVar("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 {
43  const auto & lower_domains = _lowerd_var.activeSubdomains();
44  for (const auto & id : _mesh.interiorLowerDBlocks())
45  if (lower_domains.count(id) == 0)
47  "moose",
48  29151,
49  "'lowerd_variable' must be defined on the interior lower-dimensional subdomain '" +
51  "' that is added by Mesh/build_all_side_lowerd_mesh=true.\nThe check could be overly "
52  "restrictive.");
53 
54  for (const auto & id : _var.activeSubdomains())
55  if (_mesh.interiorLowerDBlocks().count(id) > 0)
56  paramError("variable",
57  "Must not be defined on the interior lower-dimensional subdomain'" +
58  _mesh.getSubdomainName(id) + "'");
59 
60  // Note: the above two conditions also ensure that the variable and lower-d variable are
61  // different.
62 }
63 
64 void
66 {
67  if (!excludeBoundary())
68  {
70 
71  // Compute the residual for this element
73 
74  // Compute the residual for the neighbor
76 
78  }
79 }
80 
81 void
83 {
85 
86  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
87  {
89  for (_i = 0; _i < _test_lambda.size(); _i++)
91  }
92 
94 }
95 
96 void
98 {
99  if (!excludeBoundary())
100  {
102 
103  // Compute element-element Jacobian
105 
106  // Compute element-neighbor Jacobian
108 
109  // Compute neighbor-element Jacobian
111 
112  // Compute neighbor-neighbor Jacobian
114 
115  // Compute the other five pieces of Jacobian related with lower-d variable
117  }
118 }
119 
120 void
122 {
125  "Jacobian types without lower should be handled in computeElemNeighJacobian");
126 
127  const auto & test_space =
129  ? _test_lambda
131  const auto ivar =
133  ? _lowerd_var.number()
134  : _var.number();
135 
136  const auto & loc_phi =
138  ? _phi_lambda
140  const auto jvar =
142  ? _lowerd_var.number()
143  : _var.number();
144 
145  prepareMatrixTagLower(_assembly, ivar, jvar, type);
146 
147  if (_local_ke.n() == 0 || _local_ke.m() == 0)
148  return;
149 
150  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
151  {
153  for (_i = 0; _i < test_space.size(); _i++)
154  for (_j = 0; _j < loc_phi.size(); _j++)
156  }
157 
159 }
160 
161 void
162 DGLowerDKernel::computeOffDiagJacobian(const unsigned int jvar_num)
163 {
164  if (!excludeBoundary())
165  {
166  precalculateOffDiagJacobian(jvar_num);
167 
168  if (jvar_num == variable().number())
169  {
176  }
177  else if (jvar_num == _lowerd_var.number())
178  {
182  }
183  else
184  {
185  const auto & jvar = getVariable(jvar_num);
186 
187  // Compute element-element Jacobian
189 
190  // Compute element-neighbor Jacobian
192 
193  // Compute neighbor-element Jacobian
195 
196  // Compute neighbor-neighbor Jacobian
198 
199  // Compute the other five pieces of Jacobian related with lower-d variable
205  }
206  }
207 }
208 
209 void
211  const MooseVariableFEBase & jvar)
212 {
215  "Jacobian types without lower should be handled in computeOffDiagElemNeighJacobian");
216 
217  const auto & test_space =
219  ? _test_lambda
221  const auto ivar =
223  ? _lowerd_var.number()
224  : _var.number();
225 
226  prepareMatrixTagLower(_assembly, ivar, jvar.number(), type);
227  if (_local_ke.n() == 0 || _local_ke.m() == 0)
228  return;
229 
231  {
232  const auto & jv0 = static_cast<const MooseVariable &>(jvar);
233  const VariableTestValue & loc_phi =
235  ? jv0.phiLower()
236  : (type == Moose::LowerPrimary ? jv0.phiFace() : jv0.phiFaceNeighbor());
237 
238  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
239  {
241  for (_i = 0; _i < test_space.size(); _i++)
242  for (_j = 0; _j < loc_phi.size(); _j++)
244  }
245  }
247  mooseError("Array variable cannot be coupled into DG kernel currently");
248  else
249  mooseError("Vector variable cannot be coupled into DG kernel currently");
250 
252 }
MooseMesh & _mesh
Reference to this Kernel&#39;s mesh object.
const QBase *const & _qrule
Quadrature rule.
Definition: DGKernelBase.h:107
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.
const MooseVariable & _lowerd_var
Variable this kernel operates on.
virtual void precalculateOffDiagJacobian(unsigned int)
unsigned int number() const
Get variable number coming from libMesh.
virtual void computeOffDiagLowerDJacobian(Moose::ConstraintJacobianType type, const MooseVariableFEBase &jvar)
Computes one of the five pieces of off-diagonal Jacobian involving lower-d.
virtual void computeElemNeighResidual(Moose::DGResidualType type) override
Computes the residual for this element or the neighbor.
Definition: DGKernel.C:111
virtual Real computeLowerDQpResidual()=0
Method for computing the Lower part of residual at quadrature points, to be filled in residual...
virtual void precalculateResidual()
const VariablePhiValue & _phi_neighbor
Side shape function.
Definition: DGKernel.h:106
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
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
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
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:1752
const MooseVariableFEBase & variable() const override final
The variable that this kernel operates on.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar...
MooseVariable & _var
Variable this kernel operates on.
Definition: DGKernel.h:92
static InputParameters validParams()
Factory constructor initializes all internal references needed for residual computation.
Definition: DGKernel.C:27
const VariableTestValue & _test_neighbor
Side test function.
Definition: DGKernel.h:110
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 VariablePhiValue & _phi_lambda
Shape functions.
unsigned int _j
Definition: DGKernelBase.h:136
const VariableTestValue & _test_lambda
test functions
The DGKernel class is responsible for calculating the residuals for various physics on internal sides...
Definition: DGKernel.h:18
virtual void precalculateJacobian()
virtual void initLowerDQpResidual()
Put necessary evaluations depending on qp but independent on test functions here. ...
OutputTools< Real >::VariableTestValue VariableTestValue
Definition: MooseTypes.h:324
virtual void computeResidual() override
Computes the residual for this element, the neighbor and the lower-d element.
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
virtual void computeJacobian() override
Computes the nine pieces of element/neighbor/lower-d - element/neighbor/lower-d Jacobian.
virtual void computeLowerDResidual()
Computes the Lower part of residual for the variable on the lower-d element.
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 ...
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
virtual Real computeLowerDQpJacobian(Moose::ConstraintJacobianType jacobian_type)=0
Computes one of the five pieces of Jacobian involving lower-d at quadrature points.
virtual void computeLowerDJacobian(Moose::ConstraintJacobianType jacobian_type)
Computes one of the five pieces of Jacobian involving lower-d.
virtual void computeElemNeighJacobian(Moose::DGJacobianType type) override
Computes the element/neighbor-element/neighbor Jacobian.
Definition: DGKernel.C:148
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 initLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Put necessary evaluations depending on qp but independent on test and shape functions here for off-di...
const VariableTestValue & _test
test functions
Definition: DGKernel.h:102
const std::set< SubdomainID > & activeSubdomains() const
The subdomains the variable is active on.
ConstraintJacobianType
Definition: MooseTypes.h:796
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 Real computeLowerDQpOffDiagJacobian(Moose::ConstraintJacobianType, const MooseVariableFEBase &)
Computes one of the five pieces of off-diagonal Jacobian involving lower-d at quadrature points...
DGLowerDKernel(const InputParameters &parameters)
virtual Moose::VarFieldType fieldType() const =0
Field type of this variable.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
const VariablePhiValue & _phi
Shape functions.
Definition: DGKernel.h:98
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 computeOffDiagElemNeighJacobian(Moose::DGJacobianType type, const MooseVariableFEBase &jvar) override
Computes the element-element off-diagonal Jacobian.
Definition: DGKernel.C:189
virtual void initLowerDQpJacobian(Moose::ConstraintJacobianType)
Put necessary evaluations depending on qp but independent on test and shape functions here...