Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
ADMortarScalarBase.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 "ADMortarScalarBase.h"
11 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "SystemBase.h"
15 #include "MooseVariable.h"
16 #include "MooseVariableScalar.h"
17 #include "ADUtils.h"
18 
21 {
23  // This parameter can get renamed in derived class to a more relevant variable name
24  params.addCoupledVar("scalar_variable", "Primary coupled scalar variable");
25  params.addParam<bool>("compute_scalar_residuals", true, "Whether to compute scalar residuals");
26  return params;
27 }
28 
30  : ADMortarConstraint(parameters),
31  _use_scalar(isParamValid("scalar_variable") ? true : false),
32  _compute_scalar_residuals(!_use_scalar ? false : getParam<bool>("compute_scalar_residuals")),
33  _kappa_var_ptr(_use_scalar ? getScalarVar("scalar_variable", 0) : nullptr),
34  _kappa_var(_use_scalar ? _kappa_var_ptr->number() : 0),
35  _k_order(_use_scalar ? _kappa_var_ptr->order() : 0),
36  _kappa(_use_scalar ? _kappa_var_ptr->adSln() : _ad_zero)
37 {
38 }
39 
40 void
42 {
44 
46  return;
47 
48  std::vector<Real> scalar_residuals(_k_order);
49  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
50  {
52  for (_h = 0; _h < _k_order; _h++)
53  scalar_residuals[_h] += _JxW_msm[_qp] * _coord[_qp] * raw_value(computeScalarQpResidual());
54  }
57 }
58 
59 void
61 {
62  // d-_var-residual / d-_var and d-_var-residual / d-jvar
64 
66  return;
67 
68  std::vector<ADReal> scalar_residuals;
69  scalar_residuals.resize(_k_order, 0);
70  for (_qp = 0; _qp < _qrule_msm->n_points(); _qp++)
71  {
73  for (_h = 0; _h < _k_order; _h++)
74  scalar_residuals[_h] += _JxW_msm[_qp] * _coord[_qp] * computeScalarQpResidual();
75  }
78 }
static InputParameters validParams()
const MooseVariableScalar *const _kappa_var_ptr
(Pointer to) Scalar variable this kernel operates on
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.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:73
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const bool _compute_scalar_residuals
Whether to compute scalar contributions for this instance.
const std::vector< Real > & _JxW_msm
The element Jacobian times weights.
const unsigned int _k_order
Order of the scalar variable, used in several places.
virtual void computeResidual() override
Method for computing the residual.
const libMesh::QBase *const & _qrule_msm
The quadrature rule on the mortar segment element.
void addResidualsAndJacobianWithoutConstraints(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the provided d...
virtual const std::vector< dof_id_type > & dofIndices() const
Get local DoF indices.
unsigned int n_points() const
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 computeResidual() override
Computes _var-residuals as well as _kappa-residual.
virtual void computeJacobian() override
Method for computing the Jacobian.
static InputParameters validParams()
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...
ADMortarScalarBase(const InputParameters &parameters)
virtual void initScalarQpResidual()
Put necessary evaluations depending on qp but independent of test functions here. ...
const MooseArray< Real > & _coord
Member for handling change of coordinate systems (xyz, rz, spherical)
unsigned int _qp
Definition: Constraint.h:36
void scalingFactor(const std::vector< Real > &factor)
Set the scaling factor for this variable.
unsigned int _h
Used internally to iterate over each scalar component.
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 ADReal computeScalarQpResidual()
Method for computing the scalar part of residual at quadrature points.