www.mooseframework.org
VectorIntegratedBC.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "VectorIntegratedBC.h"
11 #include "Assembly.h"
12 #include "SubProblem.h"
13 #include "SystemBase.h"
14 #include "MooseVariableFE.h"
15 #include "MooseVariableScalar.h"
16 
17 #include "libmesh/quadrature.h"
18 
21 {
23  return params;
24 }
25 
27  : IntegratedBCBase(parameters),
29  false,
30  "variable",
33  _var(*mooseVariable()),
34  _normals(_assembly.normals()),
35  _phi(_assembly.phiFace(_var)),
36  _test(_var.phiFace()),
37  _u(_is_implicit ? _var.sln() : _var.slnOld())
38 {
40 }
41 
42 void
44 {
46 
47  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
48  for (_i = 0; _i < _test.size(); _i++)
49  {
50  Real residual = _JxW[_qp] * _coord[_qp] * computeQpResidual();
51  _local_re(_i) += residual;
52  }
53 
55 }
56 
57 void
59 {
61 
62  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
63  for (_i = 0; _i < _test.size(); _i++)
64  for (_j = 0; _j < _phi.size(); _j++)
66 
68 }
69 
70 void
71 VectorIntegratedBC::computeOffDiagJacobian(const unsigned int jvar_num)
72 {
73  const auto & jvar = getVariable(jvar_num);
74 
75  prepareMatrixTag(_assembly, _var.number(), jvar_num);
76 
77  // This (undisplaced) jvar could potentially yield the wrong phi size if this object is acting on
78  // the displaced mesh
79  auto phi_size = jvar.dofIndices().size();
80 
81  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
82  for (_i = 0; _i < _test.size(); _i++)
83  for (_j = 0; _j < phi_size; _j++)
84  {
85  if (_var.number() == jvar_num)
87  else
89  }
90 
92 }
93 
94 void
96 {
98 
100  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
101  for (_i = 0; _i < _test.size(); _i++)
102  for (_j = 0; _j < jv.order(); _j++)
104 
106 }
VarFieldType
Definition: MooseTypes.h:634
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
virtual Real computeQpOffDiagJacobianScalar(unsigned int)
Method for computing an off-diagonal jacobian component from a scalar var.
const VectorMooseVariable & _var
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
unsigned int number() const
Get variable number coming from libMesh.
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
unsigned int _i
i-th, j-th index for enumerating test and shape functions
THREAD_ID _tid
The thread ID for this kernel.
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.
unsigned int _qp
quadrature point index
virtual Real computeQpOffDiagJacobian(unsigned int)
Method for computing an off-diagonal jacobian component at quadrature points.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-ivar-residual / d-jvar...
SystemBase & _sys
Reference to the EquationSystem object.
MooseVariableFE< RealVectorValue > * mooseVariable() const
Return the MooseVariableFE object that this interface acts on.
virtual Real computeQpJacobian()
Method for computing the diagonal Jacobian at quadrature points.
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:627
const VectorVariableTestValue & _test
test function values (in QPs)
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
const MooseArray< Real > & _coord
coordinate transformation
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
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:141
const QBase *const & _qrule
active quadrature rule
void addMooseVariableDependency(MooseVariableFieldBase *var)
Call this function to add the passed in MooseVariableFieldBase as a variable that this object depends...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real computeQpResidual()=0
Method for computing the residual at quadrature points.
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
static InputParameters validParams()
Interface for objects that need to get values of MooseVariables.
Class for scalar variables (they are different).
Base class for deriving any boundary condition of a integrated type.
static InputParameters validParams()
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
VectorIntegratedBC(const InputParameters &parameters)
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual according to active tags.
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing element jacobian according to the active tags.
const MooseArray< Real > & _JxW
transformed Jacobian weights
const VectorVariablePhiValue & _phi
shape function values (in QPs)