https://mooseframework.inl.gov
MortarScalarBase.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 "MortarScalarBase.h"
11 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "SystemBase.h"
15 #include "MooseVariable.h"
16 #include "MooseVariableScalar.h"
17 
20 {
22  // This parameter can get renamed in derived class to a more relevant variable name
23  params.addCoupledVar("scalar_variable", "Primary coupled scalar variable");
24  params.addParam<bool>("compute_scalar_residuals", true, "Whether to compute scalar residuals");
25  return params;
26 }
27 
29  : MortarConstraint(parameters),
30  _use_scalar(isParamValid("scalar_variable") ? true : false),
31  _compute_scalar_residuals(!_use_scalar ? false : getParam<bool>("compute_scalar_residuals")),
32  _kappa_var_ptr(_use_scalar ? getScalarVar("scalar_variable", 0) : nullptr),
33  _kappa_var(_use_scalar ? _kappa_var_ptr->number() : 0),
34  _k_order(_use_scalar ? _kappa_var_ptr->order() : 0),
35  _kappa(_use_scalar ? (_is_implicit ? _kappa_var_ptr->sln() : _kappa_var_ptr->slnOld()) : _zero)
36 {
37 }
38 
39 void
41 {
43 
45  return;
46 
47  std::vector<Real> scalar_residuals(_k_order);
48  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
49  {
51  for (_h = 0; _h < _k_order; _h++)
52  scalar_residuals[_h] += _JxW_msm[_qp] * _coord[_qp] * computeScalarQpResidual();
53  }
54 
57 }
58 
59 void
61 {
62  // d-_var-residual / d-_var and d-_var-residual / d-jvar
64 
65  if (!_use_scalar)
66  return;
67 
68  // Get the list of coupled scalar vars and compute their off-diag jacobians
69  const auto & coupled_scalar_vars = getCoupledMooseScalarVars();
70 
71  // Handle ALL d-_var-residual / d-scalar columns like computeOffDiagJacobianScalar
73  // Do: dvar / dscalar_var, only want to process only nl-variables (not aux ones)
74  for (const auto & svariable : coupled_scalar_vars)
75  if (_sys.hasScalarVariable(svariable->name()))
76  {
77  // Compute the jacobian for the secondary interior primal dofs
79  // Compute the jacobian for the primary interior primal dofs.
81  }
82 
84  // Do: dvar / dscalar_var, only want to process only nl-variables (not aux ones)
85  for (const auto & svariable : coupled_scalar_vars)
86  if (_sys.hasScalarVariable(svariable->name()))
87  // Compute the jacobian for the lower dimensional LM dofs (if we even have an LM variable)
89 
91  {
92  // Handle ALL d-_kappa-residual / d-_var and d-_kappa-residual / d-jvar columns
94 
95  // Do: d-_kappa-residual / d-_kappa and d-_kappa-residual / d-jvar,
96  // only want to process only nl-variables (not aux ones)
97  for (const auto & svariable : coupled_scalar_vars)
98  {
99  if (_sys.hasScalarVariable(svariable->name()))
100  {
101  const unsigned int svar_num = svariable->number();
102  if (svar_num == _kappa_var)
103  computeScalarJacobian(); // d-_kappa-residual / d-_kappa
104  else
105  computeScalarOffDiagJacobianScalar(svar_num); // d-_kappa-residual / d-svar
106  }
107  }
108  }
109 }
110 
111 void
113 {
115  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
116  {
118  for (_h = 0; _h < _k_order; _h++)
119  for (_l = 0; _l < _k_order; _l++)
121  }
122 
124  _local_ke,
128 }
129 
130 void
132 {
133  typedef Moose::MortarType MType;
134  std::array<MType, 3> mortar_types = {{MType::Secondary, MType::Primary, MType::Lower}};
135 
137  for (const auto & it : ce)
138  {
139  MooseVariableScalar & ivariable = *(it.first);
140  MooseVariableFEBase & jvariable = *(it.second);
141 
142  const unsigned int ivar_num = ivariable.number();
143  const unsigned int jvar_num = jvariable.number();
144 
145  if (ivar_num != _kappa_var) // only do the row for _kappa_var in this object
146  continue;
147 
148  // Load shape functions of different types for easy access; identical to MortarConstraint.C
149  std::array<size_t, 3> shape_space_sizes{{jvariable.dofIndices().size(),
150  jvariable.dofIndicesNeighbor().size(),
151  jvariable.dofIndicesLower().size()}};
152  std::array<const VariablePhiValue *, 3> phis;
153  std::array<const VariablePhiGradient *, 3> grad_phis;
154  std::array<const VectorVariablePhiValue *, 3> vector_phis;
155  std::array<const VectorVariablePhiGradient *, 3> vector_grad_phis;
156  if (jvariable.isVector())
157  {
158  const auto & temp_var = static_cast<MooseVariableFE<RealVectorValue> &>(jvariable);
159  vector_phis = {{&temp_var.phiFace(), &temp_var.phiFaceNeighbor(), &temp_var.phiLower()}};
160  vector_grad_phis = {
161  {&temp_var.gradPhiFace(), &temp_var.gradPhiFaceNeighbor(), &temp_var.gradPhiLower()}};
162  }
163  else
164  {
165  const auto & temp_var = static_cast<MooseVariableFE<Real> &>(jvariable);
166  phis = {{&temp_var.phiFace(), &temp_var.phiFaceNeighbor(), &temp_var.phiLower()}};
167  grad_phis = {
168  {&temp_var.gradPhiFace(), &temp_var.gradPhiFaceNeighbor(), &temp_var.gradPhiLower()}};
169  }
170 
171  // Loop over 3 types of spatial variables, find out what jvar_num is
172  for (MooseIndex(3) type_index = 0; type_index < 3; ++type_index)
173  {
174  const auto mortar_type = mortar_types[type_index];
175  const auto shape_space_size = shape_space_sizes[type_index];
176  std::vector<dof_id_type> dof_indices;
177  switch (mortar_type)
178  {
179  case MType::Secondary:
180  dof_indices = jvariable.dofIndices();
181  break;
182 
183  case MType::Primary:
184  dof_indices = jvariable.dofIndicesNeighbor();
185  break;
186 
187  case MType::Lower:
188  dof_indices = jvariable.dofIndicesLower();
189  break;
190  }
191 
193  if (jvariable.isVector())
194  {
195  _vector_phi = vector_phis[type_index];
196  _vector_grad_phi = vector_grad_phis[type_index];
197  }
198  else
199  {
200  _phi = phis[type_index];
201  _grad_phi = grad_phis[type_index];
202  }
203 
204  _local_ke.resize(_k_order, shape_space_size);
205 
206  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
207  {
208  initScalarQpOffDiagJacobian(mortar_type, jvar_num);
209  const Real dV = _JxW_msm[_qp] * _coord[_qp];
210  for (_h = 0; _h < _k_order; _h++)
211  {
212  for (_j = 0; _j < shape_space_size; _j++)
213  {
214  _local_ke(_h, _j) += computeScalarQpOffDiagJacobian(mortar_type, jvar_num) * dV;
215  }
216  }
217  }
218 
220  _local_ke,
222  dof_indices,
224  }
225  }
226 }
227 
228 void
230 {
231  unsigned int test_space_size = 0;
232  std::vector<dof_id_type> dof_indices;
233  Real scaling_factor = 1;
234  switch (mortar_type)
235  {
237  test_space_size = _test_secondary.size();
238  dof_indices = _secondary_var.dofIndices();
239  scaling_factor = _secondary_var.scalingFactor();
240  break;
241 
243  test_space_size = _test_primary.size();
244  dof_indices = _primary_var.dofIndicesNeighbor();
245  scaling_factor = _primary_var.scalingFactor();
246  break;
247 
249  mooseAssert(_var, "LM variable is null");
250  test_space_size = _test.size();
251  dof_indices = _var->dofIndicesLower();
252  scaling_factor = _var->scalingFactor();
253  break;
254  }
255 
256  // Get dofs and order of this scalar; at least one will be _kappa_var
257  const auto & svar = _sys.getScalarVariable(_tid, svar_num);
258  const unsigned int s_order = svar.order();
259 
260  _local_ke.resize(test_space_size, s_order);
261 
262  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
263  {
264  initScalarQpOffDiagJacobian(mortar_type, svar_num);
265  const Real dV = _JxW_msm[_qp] * _coord[_qp];
266  for (_h = 0; _h < s_order; _h++)
267  for (_i = 0; _i < test_space_size; _i++)
268  { // This assumes Galerkin, i.e. the test and trial functions are the
269  // same
270  _j = _i;
271  _local_ke(_i, _h) += computeQpOffDiagJacobianScalar(mortar_type, svar_num) * dV;
272  }
273  }
274 
275  addJacobian(_assembly, _local_ke, dof_indices, svar.dofIndices(), scaling_factor);
276 }
277 
278 void
280 {
281  // Get dofs and order of this scalar; will NOT be _kappa_var
282  const auto & svar = _sys.getScalarVariable(_tid, svar_num);
283  const unsigned int s_order = svar.order();
284 
285  _local_ke.resize(_k_order, s_order);
286  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
287  {
288  initScalarQpJacobian(svar_num);
289  for (_h = 0; _h < _k_order; _h++)
290  for (_l = 0; _l < s_order; _l++)
291  _local_ke(_h, _l) +=
293  }
294 
296  _local_ke,
298  svar.dofIndices(),
300 }
void computeScalarOffDiagJacobian()
Method for computing an off-diagonal jacobian component d-_kappa-residual / d-jvar.
MooseVariableField< Real > & _secondary_var
Reference to the secondary variable.
static InputParameters validParams()
const VectorVariablePhiValue * _vector_phi
The current shape functions for vector variables.
virtual void computeJacobian() override
Computes d-_var-residual / d-_var and d-_var-residual / d-jvar, as well as d-_kappa-residual / d-_var...
virtual Real computeScalarQpOffDiagJacobian(const Moose::MortarType, const unsigned int)
Method for computing an off-diagonal jacobian component at quadrature points.
const MooseVariableScalar *const _kappa_var_ptr
(Pointer to) Scalar variable this kernel operates on
const VariableTestValue & _test_secondary
The shape functions corresponding to the secondary interior primal variable.
static InputParameters validParams()
const VariablePhiValue * _phi
The current shape functions.
unsigned int number() const
Get variable number coming from libMesh.
MortarType
Definition: MooseTypes.h:771
void computeOffDiagJacobianScalar(unsigned int) override final
Computes jacobian block with respect to a scalar variable.
const VariablePhiGradient * _grad_phi
The current shape function gradients.
void addResiduals(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided incoming residuals corresponding to the provided dof indices.
const std::vector< dof_id_type > & dofIndicesLower() const final
Get dof indices for the current lower dimensional element (this is meaningful when performing mortar ...
const VectorVariablePhiGradient * _vector_grad_phi
The current shape function gradients for vector variables.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
virtual void computeResidual() override
Computes _var-residuals as well as _kappa-residual.
const bool _compute_primal_residuals
Whether to compute primal residuals.
const std::vector< Real > & _JxW_msm
The element Jacobian times weights.
This class provides an interface for common operations on field variables of both FE and FV types wit...
unsigned int _i
Definition: Constraint.h:35
const bool _compute_scalar_residuals
Whether to compute scalar contributions for this instance.
THREAD_ID _tid
The thread ID for this kernel.
unsigned int _h
Used internally to iterate over each scalar component.
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
void computeScalarOffDiagJacobianScalar(const unsigned int svar_num)
Method for computing an off-diagonal jacobian component d-_kappa-residual / d-scalar.
virtual void computeResidual() override
Method for computing the residual.
const libMesh::QBase *const & _qrule_msm
The quadrature rule on the mortar segment element.
const VariableTestValue & _test
The shape functions corresponding to the lagrange multiplier variable.
SystemBase & _sys
Reference to the EquationSystem object.
void addJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided residual derivatives into the Jacobian for the provided dof indices.
unsigned int _j
Definition: Constraint.h:35
virtual bool isVector() const =0
virtual const std::vector< dof_id_type > & dofIndicesNeighbor() const =0
Get neighbor DOF indices for currently selected element.
const std::vector< MooseVariableScalar * > & getCoupledMooseScalarVars()
Get the list of coupled scalar variables.
const bool _use_scalar
Whether a scalar variable is declared for this constraint.
virtual void initScalarQpResidual()
Put necessary evaluations depending on qp but independent of test functions here. ...
const unsigned int _kappa_var
The unknown scalar variable ID.
virtual Real computeScalarQpResidual()
Method for computing the scalar part of residual at quadrature points.
virtual Real computeQpOffDiagJacobianScalar(const Moose::MortarType, unsigned int)
For coupling scalar variables.
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
virtual void initScalarQpOffDiagJacobian(const Moose::MortarType, const unsigned int)
Put necessary evaluations depending on qp but independent of test and shape functions here for off-di...
unsigned int n_points() const
virtual Real computeScalarQpJacobian()
Method for computing the scalar variable part of Jacobian at quadrature points.
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
void addCoupledVar(const std::string &name, const std::string &doc_string)
This method adds a coupled variable name pair.
virtual void initScalarQpJacobian(const unsigned int)
Put necessary evaluations depending on qp but independent of test and shape functions here...
libMesh::Order order() const
Get the order of this variable Note: Order enum can be implicitly converted to unsigned int...
virtual MooseVariableScalar & getScalarVariable(THREAD_ID tid, const std::string &var_name) const
Gets a reference to a scalar variable with specified number.
Definition: SystemBase.C:145
MooseVariable *const _var
Pointer to the lagrange multipler variable. nullptr if none.
const std::vector< std::pair< MooseVariableScalar *, MooseVariableFieldBase * > > & scalarFieldCouplingEntries() const
Definition: Assembly.h:1274
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const unsigned int _k_order
Order of the scalar variable, used in several places.
virtual void computeJacobian() override
Method for computing the Jacobian.
void resize(const unsigned int new_m, const unsigned int new_n)
Class for scalar variables (they are different).
const FieldVariablePhiValue & phiFace() const override final
Return the variable&#39;s shape functions on an element face.
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an optional parameter and a documentation string to the InputParameters object...
virtual const std::vector< dof_id_type > & dofIndicesLower() const =0
Get dof indices for the current lower dimensional element (this is meaningful when performing mortar ...
MortarScalarBase(const InputParameters &parameters)
const MooseArray< Real > & _coord
Member for handling change of coordinate systems (xyz, rz, spherical)
virtual bool hasScalarVariable(const std::string &var_name) const
Definition: SystemBase.C:868
virtual Real computeScalarQpOffDiagJacobianScalar(const unsigned int)
Method for computing an off-diagonal jacobian component at quadrature points.
const VariableTestValue & _test_primary
The shape functions corresponding to the primary interior primal variable.
virtual void computeScalarJacobian()
Method for computing the scalar variable part of Jacobian.
unsigned int _qp
Definition: Constraint.h:36
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
const bool _compute_lm_residuals
Whether to compute lagrange multiplier residuals.
MooseVariableField< Real > & _primary_var
Reference to the primary variable.