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 
19 template <>
22 {
24  params += validParams<RandomInterface>();
26  return params;
27 }
28 
30  : IntegratedBCBase(parameters),
32  false,
33  "variable",
36  _var(*mooseVariable()),
37  _normals(_assembly.normals()),
38  _phi(_assembly.phiFace(_var)),
39  _test(_var.phiFace()),
40  _u(_is_implicit ? _var.sln() : _var.slnOld())
41 {
43 }
44 
45 void
47 {
49 
50  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
51  for (_i = 0; _i < _test.size(); _i++)
52  {
53  Real residual = _JxW[_qp] * _coord[_qp] * computeQpResidual();
54  _local_re(_i) += residual;
55  }
56 
58 }
59 
60 void
62 {
64 
65  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
66  for (_i = 0; _i < _test.size(); _i++)
67  for (_j = 0; _j < _phi.size(); _j++)
69 
71 }
72 
73 void
75 {
76  size_t jvar_num = jvar.number();
77  prepareMatrixTag(_assembly, _var.number(), jvar_num);
78 
79  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
80  for (_i = 0; _i < _test.size(); _i++)
81  for (_j = 0; _j < jvar.phiFaceSize(); _j++)
82  {
83  if (_var.number() == jvar_num)
85  else
87  }
88 
90 }
91 
92 void
94 {
96 
98  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
99  for (_i = 0; _i < _test.size(); _i++)
100  for (_j = 0; _j < jv.order(); _j++)
102 
104 }
InputParameters validParams< MaterialPropertyInterface >()
VarFieldType
Definition: MooseTypes.h:488
InputParameters validParams< IntegratedBCBase >()
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
void computeJacobianBlockScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
virtual Real computeQpJacobian()
virtual void computeResidual() override
VectorValue< Real > RealVectorValue
Definition: Assembly.h:31
unsigned int number() const
Get variable number coming from libMesh.
virtual void computeJacobian() override
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
DenseMatrix< Number > _local_ke
Holds residual entries as they are accumulated by this Kernel.
unsigned int _qp
quadrature point index
void addMooseVariableDependency(MooseVariableFEBase *var)
Call this function to add the passed in MooseVariableFEBase as a variable that this object depends on...
MooseVariableFE< RealVectorValue > * mooseVariable() const
Get the variable that this object is using.
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:259
virtual void computeJacobianBlock(MooseVariableFEBase &jvar) override
Computes d-ivar-residual / d-jvar...
InputParameters validParams< RandomInterface >()
VarKindType
Framework-wide stuff.
Definition: MooseTypes.h:481
Assembly & _assembly
Reference to assembly.
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
virtual size_t phiFaceSize() const =0
Return phiFace size.
const QBase *const & _qrule
active quadrature rule
InputParameters validParams< VectorIntegratedBC >()
virtual Real computeQpResidual()=0
method for computing the residual at quadrature points
SystemBase & _sys
Reference to SystemBase.
virtual Real computeQpOffDiagJacobian(unsigned int)
This is the virtual that derived classes should override for computing an off-diagonal jacobian compo...
DenseVector< Number > _local_re
Holds residual entries as they are accumulated by this Kernel.
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.
Definition: Moose.h:112
VectorIntegratedBC(const InputParameters &parameters)
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual the according to active tags.
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing element jacobian according to the ative tags.
VectorMooseVariable & _var
THREAD_ID _tid
Thread id.
virtual MooseVariableScalar & getScalarVariable(THREAD_ID tid, const std::string &var_name)
Gets a reference to a scalar variable with specified number.
Definition: SystemBase.C:149
const MooseArray< Real > & _JxW
transformed Jacobian weights
const VectorVariablePhiValue & _phi
shape function values (in QPs)