www.mooseframework.org
Public Member Functions | Protected Attributes | List of all members
TaggingInterface Class Reference

#include <TaggingInterface.h>

Inheritance diagram for TaggingInterface:
[legend]

Public Member Functions

 TaggingInterface (const MooseObject *moose_object)
 
virtual ~TaggingInterface ()
 
void useVectorTag (const TagName &tag_name)
 
void useMatrixTag (const TagName &tag_name)
 
void useVectorTag (TagID tag_id)
 
void useMatrixTag (TagID tag_id)
 
bool isVectorTagged ()
 
bool isMatrixTagged ()
 
const std::set< TagID > & getVectorTags () const
 
const std::set< TagID > & getMatrixTags () const
 
void prepareVectorTag (Assembly &assembly, unsigned int ivar)
 Prepare data for computing element residual the according to active tags. More...
 
void prepareVectorTagNeighbor (Assembly &assembly, unsigned int ivar)
 Prepare data for computing element residual the according to active tags for DG and interface kernels. More...
 
void prepareMatrixTag (Assembly &assembly, unsigned int ivar, unsigned int jvar)
 Prepare data for computing element jacobian according to the ative tags. More...
 
void prepareMatrixTagNeighbor (Assembly &assembly, unsigned int ivar, unsigned int jvar, Moose::DGJacobianType type)
 Prepare data for computing element jacobian according to the ative tags for DG and interface kernels. More...
 
void accumulateTaggedLocalResidual ()
 Local residual blocks will be appended by adding the current local kernel residual. More...
 
void assignTaggedLocalResidual ()
 Local residual blocks will assigned as the current local kernel residual. More...
 
void accumulateTaggedLocalMatrix ()
 Local Jacobian blocks will be appended by adding the current local kernel Jacobian. More...
 
void assignTaggedLocalMatrix ()
 Local Jacobian blocks will assigned as the current local kernel Jacobian. More...
 

Protected Attributes

std::set< TagID_vector_tags
 The vectors this Kernel will contribute to. More...
 
std::set< TagID_matrix_tags
 The matrices this Kernel will contribute to. More...
 
const MooseObject_moose_object
 Moose objct this tag works on. More...
 
const InputParameters_tag_params
 Parameters from moose object. More...
 
SubProblem_subproblem
 SubProblem that contains tag info. More...
 
std::vector< DenseVector< Number > * > _re_blocks
 Residual blocks Vectors For each Tag. More...
 
std::vector< DenseMatrix< Number > * > _ke_blocks
 Kernel blocks Vectors For each Tag. More...
 
DenseVector< Number > _local_re
 Holds residual entries as they are accumulated by this Kernel. More...
 
DenseMatrix< Number > _local_ke
 Holds residual entries as they are accumulated by this Kernel. More...
 

Detailed Description

Definition at line 37 of file TaggingInterface.h.

Constructor & Destructor Documentation

◆ TaggingInterface()

TaggingInterface::TaggingInterface ( const MooseObject moose_object)

Definition at line 49 of file TaggingInterface.C.

50  : _moose_object(*moose_object),
53 {
54  auto & vector_tag_names = _tag_params.get<MultiMooseEnum>("vector_tags");
55 
56  if (!vector_tag_names.isValid())
57  mooseError("MUST provide at least one vector_tag for Kernel: ", _moose_object.name());
58 
59  for (auto & vector_tag_name : vector_tag_names)
60  _vector_tags.insert(_subproblem.getVectorTagID(vector_tag_name.name()));
61 
62  // Add extra vector tags. These tags should be created in the System already, otherwise
63  // we can not add the extra tags
64  auto & extra_vector_tags = _tag_params.get<std::vector<TagName>>("extra_vector_tags");
65 
66  for (auto & vector_tag_name : extra_vector_tags)
67  _vector_tags.insert(_subproblem.getVectorTagID(vector_tag_name));
68 
69  auto & matrix_tag_names = _tag_params.get<MultiMooseEnum>("matrix_tags");
70 
71  if (!matrix_tag_names.isValid())
72  mooseError("MUST provide at least one matrix_tag for Kernel: ", _moose_object.name());
73 
74  for (auto & matrix_tag_name : matrix_tag_names)
75  _matrix_tags.insert(_subproblem.getMatrixTagID(matrix_tag_name.name()));
76 
77  auto & extra_matrix_tags = _tag_params.get<std::vector<TagName>>("extra_matrix_tags");
78 
79  for (auto & matrix_tag_name : extra_matrix_tags)
80  _matrix_tags.insert(_subproblem.getMatrixTagID(matrix_tag_name));
81 
82  _re_blocks.resize(_vector_tags.size());
83  _ke_blocks.resize(_matrix_tags.size());
84 }
virtual TagID getVectorTagID(const TagName &tag_name)
Get a TagID from a TagName.
Definition: SubProblem.C:80
SubProblem & _subproblem
SubProblem that contains tag info.
const InputParameters & _tag_params
Parameters from moose object.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:208
T getCheckedPointerParam(const std::string &name, const std::string &error_string="") const
Verifies that the requested parameter exists and is not NULL and returns it to the caller...
std::vector< DenseVector< Number > * > _re_blocks
Residual blocks Vectors For each Tag.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseObject.h:57
std::set< TagID > _matrix_tags
The matrices this Kernel will contribute to.
std::set< TagID > _vector_tags
The vectors this Kernel will contribute to.
Generic class for solving transient nonlinear problems.
Definition: SubProblem.h:53
virtual TagID getMatrixTagID(const TagName &tag_name)
Get a TagID from a TagName.
Definition: SubProblem.C:132
const MooseObject & _moose_object
Moose objct this tag works on.
const std::string & name() const
Get the name of the object.
Definition: MooseObject.h:51
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
std::vector< DenseMatrix< Number > * > _ke_blocks
Kernel blocks Vectors For each Tag.

◆ ~TaggingInterface()

TaggingInterface::~TaggingInterface ( )
virtual

Definition at line 205 of file TaggingInterface.C.

205 {}

Member Function Documentation

◆ accumulateTaggedLocalMatrix()

void TaggingInterface::accumulateTaggedLocalMatrix ( )

◆ accumulateTaggedLocalResidual()

void TaggingInterface::accumulateTaggedLocalResidual ( )

Local residual blocks will be appended by adding the current local kernel residual.

It should be called after the local element vector has been computed.

Definition at line 178 of file TaggingInterface.C.

Referenced by InterfaceKernel::computeElemNeighResidual(), DGKernel::computeElemNeighResidual(), Kernel::computeResidual(), ODEKernel::computeResidual(), ODETimeKernel::computeResidual(), TimeKernel::computeResidual(), IntegratedBC::computeResidual(), VectorIntegratedBC::computeResidual(), DiracKernel::computeResidual(), and ConservativeAdvection::fullUpwind().

179 {
180  for (auto & re : _re_blocks)
181  *re += _local_re;
182 }
std::vector< DenseVector< Number > * > _re_blocks
Residual blocks Vectors For each Tag.
DenseVector< Number > _local_re
Holds residual entries as they are accumulated by this Kernel.

◆ assignTaggedLocalMatrix()

void TaggingInterface::assignTaggedLocalMatrix ( )

Local Jacobian blocks will assigned as the current local kernel Jacobian.

It should be called after the local element matrix has been computed.

Definition at line 199 of file TaggingInterface.C.

Referenced by NodalEqualValueConstraint::computeJacobian().

200 {
201  for (auto & ke : _ke_blocks)
202  *ke = _local_ke;
203 }
DenseMatrix< Number > _local_ke
Holds residual entries as they are accumulated by this Kernel.
std::vector< DenseMatrix< Number > * > _ke_blocks
Kernel blocks Vectors For each Tag.

◆ assignTaggedLocalResidual()

void TaggingInterface::assignTaggedLocalResidual ( )

Local residual blocks will assigned as the current local kernel residual.

It should be called after the local element vector has been computed.

Definition at line 185 of file TaggingInterface.C.

Referenced by NodalEqualValueConstraint::computeResidual().

186 {
187  for (auto & re : _re_blocks)
188  *re = _local_re;
189 }
std::vector< DenseVector< Number > * > _re_blocks
Residual blocks Vectors For each Tag.
DenseVector< Number > _local_re
Holds residual entries as they are accumulated by this Kernel.

◆ getMatrixTags()

const std::set<TagID>& TaggingInterface::getMatrixTags ( ) const
inline

Definition at line 57 of file TaggingInterface.h.

57 { return _matrix_tags; }
std::set< TagID > _matrix_tags
The matrices this Kernel will contribute to.

◆ getVectorTags()

const std::set<TagID>& TaggingInterface::getVectorTags ( ) const
inline

Definition at line 55 of file TaggingInterface.h.

55 { return _vector_tags; }
std::set< TagID > _vector_tags
The vectors this Kernel will contribute to.

◆ isMatrixTagged()

bool TaggingInterface::isMatrixTagged ( )
inline

Definition at line 53 of file TaggingInterface.h.

53 { return _matrix_tags.size() > 0; }
std::set< TagID > _matrix_tags
The matrices this Kernel will contribute to.

◆ isVectorTagged()

bool TaggingInterface::isVectorTagged ( )
inline

Definition at line 51 of file TaggingInterface.h.

51 { return _vector_tags.size() > 0; }
std::set< TagID > _vector_tags
The vectors this Kernel will contribute to.

◆ prepareMatrixTag()

void TaggingInterface::prepareMatrixTag ( Assembly assembly,
unsigned int  ivar,
unsigned int  jvar 
)

Prepare data for computing element jacobian according to the ative tags.

Jacobian blocks for different tags will be extracted from Assembly. A local Jacobian will be zeroed. It should be called right before the local element matrix is computed.

Definition at line 149 of file TaggingInterface.C.

Referenced by InterfaceKernel::computeElemNeighJacobian(), DGKernel::computeElemNeighJacobian(), TimeDerivative::computeJacobian(), Kernel::computeJacobian(), ODEKernel::computeJacobian(), IntegratedBC::computeJacobian(), NodalEqualValueConstraint::computeJacobian(), VectorIntegratedBC::computeJacobian(), DiracKernel::computeJacobian(), IntegratedBC::computeJacobianBlock(), VectorIntegratedBC::computeJacobianBlock(), IntegratedBC::computeJacobianBlockScalar(), VectorIntegratedBC::computeJacobianBlockScalar(), InterfaceKernel::computeOffDiagElemNeighJacobian(), DGKernel::computeOffDiagElemNeighJacobian(), ODEKernel::computeOffDiagJacobian(), Kernel::computeOffDiagJacobian(), DiracKernel::computeOffDiagJacobian(), Kernel::computeOffDiagJacobianScalar(), and ConservativeAdvection::fullUpwind().

150 {
151  _ke_blocks.resize(_matrix_tags.size());
152  mooseAssert(_matrix_tags.size() >= 1, "we need at least one active tag");
153  auto mat_vector = _matrix_tags.begin();
154  for (auto i = beginIndex(_matrix_tags); i < _matrix_tags.size(); i++, ++mat_vector)
155  _ke_blocks[i] = &assembly.jacobianBlock(ivar, jvar, *mat_vector);
156 
157  _local_ke.resize(_ke_blocks[0]->m(), _ke_blocks[0]->n());
158  _local_ke.zero();
159 }
DenseMatrix< Number > & jacobianBlock(unsigned int ivar, unsigned int jvar, TagID tag=0)
Definition: Assembly.C:934
DenseMatrix< Number > _local_ke
Holds residual entries as they are accumulated by this Kernel.
PetscInt m
std::set< TagID > _matrix_tags
The matrices this Kernel will contribute to.
PetscInt n
std::vector< DenseMatrix< Number > * > _ke_blocks
Kernel blocks Vectors For each Tag.

◆ prepareMatrixTagNeighbor()

void TaggingInterface::prepareMatrixTagNeighbor ( Assembly assembly,
unsigned int  ivar,
unsigned int  jvar,
Moose::DGJacobianType  type 
)

Prepare data for computing element jacobian according to the ative tags for DG and interface kernels.

Jacobian blocks for different tags will be extracted from Assembly. A local Jacobian will be zeroed. It should be called right before the local element matrix is computed.

Definition at line 162 of file TaggingInterface.C.

Referenced by InterfaceKernel::computeElemNeighJacobian(), DGKernel::computeElemNeighJacobian(), InterfaceKernel::computeOffDiagElemNeighJacobian(), and DGKernel::computeOffDiagElemNeighJacobian().

166 {
167  _ke_blocks.resize(_matrix_tags.size());
168  mooseAssert(_matrix_tags.size() >= 1, "we need at least one active tag");
169  auto mat_vector = _matrix_tags.begin();
170  for (auto i = beginIndex(_matrix_tags); i < _matrix_tags.size(); i++, ++mat_vector)
171  _ke_blocks[i] = &assembly.jacobianBlockNeighbor(type, ivar, jvar, *mat_vector);
172 
173  _local_ke.resize(_ke_blocks[0]->m(), _ke_blocks[0]->n());
174  _local_ke.zero();
175 }
DenseMatrix< Number > _local_ke
Holds residual entries as they are accumulated by this Kernel.
PetscInt m
std::set< TagID > _matrix_tags
The matrices this Kernel will contribute to.
MatType type
PetscInt n
DenseMatrix< Number > & jacobianBlockNeighbor(Moose::DGJacobianType type, unsigned int ivar, unsigned int jvar, TagID tag=0)
Definition: Assembly.C:948
std::vector< DenseMatrix< Number > * > _ke_blocks
Kernel blocks Vectors For each Tag.

◆ prepareVectorTag()

void TaggingInterface::prepareVectorTag ( Assembly assembly,
unsigned int  ivar 
)

Prepare data for computing element residual the according to active tags.

Residual blocks for different tags will be extracted from Assembly. A local residual will be zeroed. It should be called right before the local element vector is computed.

Definition at line 123 of file TaggingInterface.C.

Referenced by InterfaceKernel::computeElemNeighResidual(), DGKernel::computeElemNeighResidual(), Kernel::computeResidual(), TimeKernel::computeResidual(), ODEKernel::computeResidual(), ODETimeKernel::computeResidual(), IntegratedBC::computeResidual(), NodalEqualValueConstraint::computeResidual(), VectorIntegratedBC::computeResidual(), DiracKernel::computeResidual(), and ConservativeAdvection::fullUpwind().

124 {
125  _re_blocks.resize(_vector_tags.size());
126  mooseAssert(_vector_tags.size() >= 1, "we need at least one active tag");
127  auto vector_tag = _vector_tags.begin();
128  for (auto i = beginIndex(_vector_tags); i < _vector_tags.size(); i++, ++vector_tag)
129  _re_blocks[i] = &assembly.residualBlock(ivar, *vector_tag);
130 
131  _local_re.resize(_re_blocks[0]->size());
132  _local_re.zero();
133 }
std::vector< DenseVector< Number > * > _re_blocks
Residual blocks Vectors For each Tag.
std::set< TagID > _vector_tags
The vectors this Kernel will contribute to.
DenseVector< Number > _local_re
Holds residual entries as they are accumulated by this Kernel.
DenseVector< Number > & residualBlock(unsigned int var_num, TagID tag_id=0)
Definition: Assembly.h:610

◆ prepareVectorTagNeighbor()

void TaggingInterface::prepareVectorTagNeighbor ( Assembly assembly,
unsigned int  ivar 
)

Prepare data for computing element residual the according to active tags for DG and interface kernels.

Residual blocks for different tags will be extracted from Assembly. A local residual will be zeroed. It should be called right before the local element vector is computed.

Definition at line 136 of file TaggingInterface.C.

Referenced by InterfaceKernel::computeElemNeighResidual(), and DGKernel::computeElemNeighResidual().

137 {
138  _re_blocks.resize(_vector_tags.size());
139  mooseAssert(_vector_tags.size() >= 1, "we need at least one active tag");
140  auto vector_tag = _vector_tags.begin();
141  for (auto i = beginIndex(_vector_tags); i < _vector_tags.size(); i++, ++vector_tag)
142  _re_blocks[i] = &assembly.residualBlockNeighbor(ivar, *vector_tag);
143 
144  _local_re.resize(_re_blocks[0]->size());
145  _local_re.zero();
146 }
DenseVector< Number > & residualBlockNeighbor(unsigned int var_num, TagID tag_id=0)
Definition: Assembly.h:615
std::vector< DenseVector< Number > * > _re_blocks
Residual blocks Vectors For each Tag.
std::set< TagID > _vector_tags
The vectors this Kernel will contribute to.
DenseVector< Number > _local_re
Holds residual entries as they are accumulated by this Kernel.

◆ useMatrixTag() [1/2]

void TaggingInterface::useMatrixTag ( const TagName &  tag_name)

Definition at line 96 of file TaggingInterface.C.

97 {
98  if (!_subproblem.matrixTagExists(tag_name))
99  mooseError("Matrix tag ", tag_name, " does not exsit in system");
100 
101  _matrix_tags.insert(_subproblem.getMatrixTagID(tag_name));
102 }
SubProblem & _subproblem
SubProblem that contains tag info.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:208
virtual bool matrixTagExists(const TagName &tag_name)
Check to see if a particular Tag exists.
Definition: SubProblem.C:118
std::set< TagID > _matrix_tags
The matrices this Kernel will contribute to.
virtual TagID getMatrixTagID(const TagName &tag_name)
Get a TagID from a TagName.
Definition: SubProblem.C:132

◆ useMatrixTag() [2/2]

void TaggingInterface::useMatrixTag ( TagID  tag_id)

Definition at line 114 of file TaggingInterface.C.

115 {
116  if (!_subproblem.matrixTagExists(tag_id))
117  mooseError("Matrix tag ", tag_id, " does not exsit in system");
118 
119  _matrix_tags.insert(tag_id);
120 }
SubProblem & _subproblem
SubProblem that contains tag info.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:208
virtual bool matrixTagExists(const TagName &tag_name)
Check to see if a particular Tag exists.
Definition: SubProblem.C:118
std::set< TagID > _matrix_tags
The matrices this Kernel will contribute to.

◆ useVectorTag() [1/2]

void TaggingInterface::useVectorTag ( const TagName &  tag_name)

Definition at line 87 of file TaggingInterface.C.

88 {
89  if (!_subproblem.vectorTagExists(tag_name))
90  mooseError("Vector tag ", tag_name, " does not exsit in system");
91 
92  _vector_tags.insert(_subproblem.getVectorTagID(tag_name));
93 }
virtual TagID getVectorTagID(const TagName &tag_name)
Get a TagID from a TagName.
Definition: SubProblem.C:80
SubProblem & _subproblem
SubProblem that contains tag info.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:208
std::set< TagID > _vector_tags
The vectors this Kernel will contribute to.
virtual bool vectorTagExists(TagID tag)
Check to see if a particular Tag exists.
Definition: SubProblem.h:100

◆ useVectorTag() [2/2]

void TaggingInterface::useVectorTag ( TagID  tag_id)

Definition at line 105 of file TaggingInterface.C.

106 {
107  if (!_subproblem.vectorTagExists(tag_id))
108  mooseError("Vector tag ", tag_id, " does not exsit in system");
109 
110  _vector_tags.insert(tag_id);
111 }
SubProblem & _subproblem
SubProblem that contains tag info.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:208
std::set< TagID > _vector_tags
The vectors this Kernel will contribute to.
virtual bool vectorTagExists(TagID tag)
Check to see if a particular Tag exists.
Definition: SubProblem.h:100

Member Data Documentation

◆ _ke_blocks

std::vector<DenseMatrix<Number> *> TaggingInterface::_ke_blocks
protected

◆ _local_ke

DenseMatrix<Number> TaggingInterface::_local_ke
protected

◆ _local_re

DenseVector<Number> TaggingInterface::_local_re
protected

◆ _matrix_tags

std::set<TagID> TaggingInterface::_matrix_tags
protected

◆ _moose_object

const MooseObject& TaggingInterface::_moose_object
protected

Moose objct this tag works on.

Definition at line 128 of file TaggingInterface.h.

Referenced by TaggingInterface().

◆ _re_blocks

std::vector<DenseVector<Number> *> TaggingInterface::_re_blocks
protected

◆ _subproblem

SubProblem& TaggingInterface::_subproblem
protected

SubProblem that contains tag info.

Definition at line 134 of file TaggingInterface.h.

Referenced by TaggingInterface(), useMatrixTag(), and useVectorTag().

◆ _tag_params

const InputParameters& TaggingInterface::_tag_params
protected

Parameters from moose object.

Definition at line 131 of file TaggingInterface.h.

Referenced by TaggingInterface().

◆ _vector_tags

std::set<TagID> TaggingInterface::_vector_tags
protected

The documentation for this class was generated from the following files: