https://mooseframework.inl.gov
TaggingInterface.h
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 #pragma once
11 
12 #include "DenseMatrix.h"
13 #include "MooseTypes.h"
14 #include "MultiMooseEnum.h"
15 #include "Assembly.h"
16 #include "MooseVariableFE.h"
17 #include "SystemBase.h"
18 
19 #include "libmesh/dense_vector.h"
20 #include "metaphysicl/raw_type.h"
21 
22 #include <vector>
23 
24 // Forward declarations
25 class InputParameters;
26 class MooseObject;
27 class SubProblem;
28 class Assembly;
29 
30 #ifdef MOOSE_KOKKOS_ENABLED
31 namespace Moose::Kokkos
32 {
33 class ResidualObject;
34 class System;
35 }
36 #endif
37 
38 template <typename T>
40 
50 {
51  const DenseVector<ADReal> & residuals;
52  const std::vector<dof_id_type> & dof_indices;
54 };
55 
57 {
58 public:
59  TaggingInterface(const MooseObject * moose_object);
60 
61 #ifdef MOOSE_KOKKOS_ENABLED
62 
66 #endif
67 
68  virtual ~TaggingInterface();
69 
71 
77  {
78  friend class AttribVectorTags;
79  friend class NonlinearEigenSystem;
81  template <typename>
83 #ifdef MOOSE_KOKKOS_ENABLED
85  friend class Moose::Kokkos::System;
86 #endif
87 
90  };
91 
97  {
98  friend class AttribMatrixTags;
99  friend class NonlinearEigenSystem;
101  template <typename>
103 #ifdef MOOSE_KOKKOS_ENABLED
105  friend class Moose::Kokkos::System;
106 #endif
107 
110  };
111 
115  enum class ResidualTagType
116  {
117  NonReference,
118  Reference
119  };
120 
121  void useVectorTag(const TagName & tag_name, VectorTagsKey);
122 
123  void useMatrixTag(const TagName & tag_name, MatrixTagsKey);
124 
125  void useVectorTag(TagID tag_id, VectorTagsKey);
126 
127  void useMatrixTag(TagID tag_id, MatrixTagsKey);
128 
129  bool isVectorTagged() { return _vector_tags.size() > 0; }
130 
131  bool isMatrixTagged() { return _matrix_tags.size() > 0; }
132 
133  bool hasVectorTags() const { return !_vector_tags.empty(); }
134 
135  const std::set<TagID> & getVectorTags(VectorTagsKey) const { return _vector_tags; }
136 
137  const std::set<TagID> & getMatrixTags(MatrixTagsKey) const { return _matrix_tags; }
138 
139 protected:
146  void prepareVectorTag(Assembly & assembly, unsigned int ivar);
147 
155  void prepareVectorTag(Assembly & assembly, unsigned int ivar, ResidualTagType tag_type);
156 
164  void prepareVectorTagNeighbor(Assembly & assembly, unsigned int ivar);
165 
171  void prepareVectorTagLower(Assembly & assembly, unsigned int ivar);
172 
179  void prepareMatrixTag(Assembly & assembly, unsigned int ivar, unsigned int jvar);
180 
181  void prepareMatrixTag(Assembly & assembly,
182  unsigned int ivar,
183  unsigned int jvar,
184  DenseMatrix<Number> & k) const;
185 
192  void prepareMatrixTagNonlocal(Assembly & assembly, unsigned int ivar, unsigned int jvar);
193 
201  void prepareMatrixTagNeighbor(Assembly & assembly,
202  unsigned int ivar,
203  unsigned int jvar,
204  Moose::DGJacobianType type);
205 
206  void prepareMatrixTagNeighbor(Assembly & assembly,
207  unsigned int ivar,
208  unsigned int jvar,
210  DenseMatrix<Number> & k) const;
211 
217  void prepareMatrixTagLower(Assembly & assembly,
218  unsigned int ivar,
219  unsigned int jvar,
221 
227 
233 
239 
240  void accumulateTaggedLocalMatrix(Assembly & assembly,
241  unsigned int ivar,
242  unsigned int jvar,
243  const DenseMatrix<Number> & k);
244 
245  void accumulateTaggedLocalMatrix(Assembly & assembly,
246  unsigned int ivar,
247  unsigned int jvar,
249  const DenseMatrix<Number> & k);
250 
256 
262 
266  template <typename Residuals, typename Indices>
267  void addResiduals(Assembly & assembly,
268  const Residuals & residuals,
269  const Indices & dof_indices,
270  Real scaling_factor);
271 
272  template <typename Residuals, typename Indices>
273  void addResiduals(Assembly & assembly,
274  const Residuals & residuals,
275  const Indices & dof_indices,
276  const std::vector<Real> & scaling_factors);
277 
281  template <typename T, typename Indices>
282  void addResiduals(Assembly & assembly,
283  const DenseVector<T> & residuals,
284  const Indices & dof_indices,
285  Real scaling_factor);
286 
291  template <typename Residuals, typename Indices>
292  void addResidualsAndJacobian(Assembly & assembly,
293  const Residuals & residuals,
294  const Indices & dof_indices,
295  Real scaling_factor);
296 
300  template <typename Residuals, typename Indices>
301  void addJacobian(Assembly & assembly,
302  const Residuals & residuals,
303  const Indices & dof_indices,
304  Real scaling_factor);
305 
310  template <typename Residuals, typename Indices>
311  void addJacobian(Assembly & assembly,
312  const Residuals & residuals,
313  const Indices & dof_indices,
314  const std::vector<Real> & scaling_factors);
315 
319  void addResiduals(Assembly & assembly, const ADResidualsPacket & packet);
320 
325  void addResidualsAndJacobian(Assembly & assembly, const ADResidualsPacket & packet);
326 
330  void addJacobian(Assembly & assembly, const ADResidualsPacket & packet);
331 
337  template <typename Residuals, typename Indices>
339  const Residuals & residuals,
340  const Indices & dof_indices,
341  Real scaling_factor);
342 
348  template <typename Residuals, typename Indices>
350  const Residuals & residuals,
351  const Indices & dof_indices,
352  Real scaling_factor);
353 
359  template <typename Residuals, typename Indices>
360  void addJacobianWithoutConstraints(Assembly & assembly,
361  const Residuals & residuals,
362  const Indices & dof_indices,
363  Real scaling_factor);
364 
368  void addJacobianElement(Assembly & assembly,
369  Real value,
370  dof_id_type row_index,
371  dof_id_type column_index,
372  Real scaling_factor);
373 
377  void addJacobian(Assembly & assembly,
378  DenseMatrix<Real> & local_k,
379  const std::vector<dof_id_type> & row_indices,
380  const std::vector<dof_id_type> & column_indices,
381  Real scaling_factor);
382 
386  template <typename T>
387  void setResidual(SystemBase & sys, const T & residual, MooseVariableFE<T> & var);
388 
392  void setResidual(SystemBase & sys, Real residual, dof_id_type dof_index);
393 
397  template <typename SetResidualFunctor>
398  void setResidual(SystemBase & sys, SetResidualFunctor set_residual_functor);
399 
402 
405 
408 
411 
412 private:
419  void prepareVectorTagInternal(Assembly & assembly,
420  unsigned int ivar,
421  const std::set<TagID> & vector_tags,
422  const std::set<TagID> & absolute_value_vector_tags);
423 
424 #ifndef NDEBUG
425 
428  void checkForNans() const;
429 #endif
430 
432  std::set<TagID> _vector_tags;
433 
435  std::set<TagID> _abs_vector_tags;
436 
438  std::set<TagID> _matrix_tags;
439 
442  std::set<TagID> _non_ref_vector_tags;
443 
446  std::set<TagID> _non_ref_abs_vector_tags;
447 
451  std::set<TagID> _ref_vector_tags;
452 
456  std::set<TagID> _ref_abs_vector_tags;
457 
460 
463 
465  std::vector<DenseVector<Number> *> _re_blocks;
466 
468  std::vector<DenseVector<Number> *> _absre_blocks;
469 
471  std::vector<DenseMatrix<Number> *> _ke_blocks;
472 
475  std::vector<Real> _absolute_residuals;
476 
477  friend class NonlinearSystemBase;
478 };
479 
480 #define usingTaggingInterfaceMembers \
481  using TaggingInterface::_subproblem; \
482  using TaggingInterface::accumulateTaggedLocalResidual; \
483  using TaggingInterface::accumulateTaggedLocalMatrix; \
484  using TaggingInterface::prepareVectorTag; \
485  using TaggingInterface::prepareMatrixTag; \
486  using TaggingInterface::prepareVectorTagNeighbor; \
487  using TaggingInterface::_local_re; \
488  using TaggingInterface::prepareVectorTagLower; \
489  using TaggingInterface::prepareMatrixTagNeighbor; \
490  using TaggingInterface::prepareMatrixTagLower; \
491  using TaggingInterface::_local_ke
492 
493 template <typename Residuals, typename Indices>
494 void
496  const Residuals & residuals,
497  const Indices & dof_indices,
498  const Real scaling_factor)
499 {
500  assembly.cacheResiduals(
501  residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _vector_tags);
502  if (!_abs_vector_tags.empty())
503  {
504  _absolute_residuals.resize(residuals.size());
505  for (const auto i : index_range(residuals))
507 
509  dof_indices,
510  scaling_factor,
513  }
514 }
515 
516 template <typename Residuals, typename Indices>
517 void
519  const Residuals & residuals,
520  const Indices & dof_indices,
521  const std::vector<Real> & scaling_factors)
522 {
523  const auto count = scaling_factors.size();
524  mooseAssert(dof_indices.size() % count == 0,
525  "The number of dof indices should be divided cleanly by the variable count");
526  const auto nshapes = dof_indices.size() / count;
527 
528  for (const auto j : make_range(count))
529  // The Residuals type may not offer operator[] (e.g. eigen vectors) but more commonly it
530  // should offer data()
531  addResiduals(assembly,
532  Moose::makeSpan(residuals, j * nshapes, nshapes),
533  Moose::makeSpan(dof_indices, j * nshapes, nshapes),
534  scaling_factors[j]);
535 }
536 
537 template <typename T, typename Indices>
538 void
540  const DenseVector<T> & residuals,
541  const Indices & dof_indices,
542  const Real scaling_factor)
543 {
544  addResiduals(assembly, residuals.get_values(), dof_indices, scaling_factor);
545 }
546 
547 template <typename Residuals, typename Indices>
548 void
550  const Residuals & residuals,
551  const Indices & dof_indices,
552  const Real scaling_factor)
553 {
555  residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _vector_tags);
556  if (!_abs_vector_tags.empty())
557  {
558  _absolute_residuals.resize(residuals.size());
559  for (const auto i : index_range(residuals))
561 
563  dof_indices,
564  scaling_factor,
567  }
568 }
569 
570 template <typename Residuals, typename Indices>
571 void
573  const Residuals & residuals,
574  const Indices & dof_indices,
575  Real scaling_factor)
576 {
577  addResiduals(assembly, residuals, dof_indices, scaling_factor);
578  addJacobian(assembly, residuals, dof_indices, scaling_factor);
579 }
580 
581 template <typename Residuals, typename Indices>
582 void
584  const Residuals & residuals,
585  const Indices & dof_indices,
586  Real scaling_factor)
587 {
588  assembly.cacheJacobian(
589  residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _matrix_tags);
590 }
591 
592 template <typename Residuals, typename Indices>
593 void
595  const Residuals & residuals,
596  const Indices & dof_indices,
597  const std::vector<Real> & scaling_factors)
598 {
599  const auto count = scaling_factors.size();
600  mooseAssert(dof_indices.size() % count == 0,
601  "The number of dof indices should be divided cleanly by the variable count");
602  const auto nshapes = dof_indices.size() / count;
603 
604  for (const auto j : make_range(count))
605  // The Residuals type may not offer operator[] (e.g. eigen vectors) but more commonly it
606  // should offer data()
607  addJacobian(assembly,
608  Moose::makeSpan(residuals, j * nshapes, nshapes),
609  Moose::makeSpan(dof_indices, j * nshapes, nshapes),
610  scaling_factors[j]);
611 }
612 
613 template <typename Residuals, typename Indices>
614 void
616  const Residuals & residuals,
617  const Indices & dof_indices,
618  Real scaling_factor)
619 {
620  addResidualsWithoutConstraints(assembly, residuals, dof_indices, scaling_factor);
621  addJacobianWithoutConstraints(assembly, residuals, dof_indices, scaling_factor);
622 }
623 
624 template <typename Residuals, typename Indices>
625 void
627  const Residuals & residuals,
628  const Indices & dof_indices,
629  Real scaling_factor)
630 {
632  residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _matrix_tags);
633 }
634 
635 inline void
637  const Real value,
638  const dof_id_type row_index,
639  const dof_id_type column_index,
640  const Real scaling_factor)
641 {
642  assembly.cacheJacobian(
643  row_index, column_index, value * scaling_factor, Assembly::LocalDataKey{}, _matrix_tags);
644 }
645 
646 inline void
648  DenseMatrix<Real> & local_k,
649  const std::vector<dof_id_type> & row_indices,
650  const std::vector<dof_id_type> & column_indices,
651  const Real scaling_factor)
652 {
653  for (const auto matrix_tag : _matrix_tags)
654  assembly.cacheJacobianBlock(
655  local_k, row_indices, column_indices, scaling_factor, Assembly::LocalDataKey{}, matrix_tag);
656 }
657 
658 template <typename T>
659 void
661 {
662  for (const auto tag_id : _vector_tags)
663  if (sys.hasVector(tag_id))
664  var.insertNodalValue(sys.getVector(tag_id), residual);
665 }
666 
667 inline void
668 TaggingInterface::setResidual(SystemBase & sys, const Real residual, const dof_id_type dof_index)
669 {
670  for (const auto tag_id : _vector_tags)
671  if (sys.hasVector(tag_id))
672  sys.getVector(tag_id).set(dof_index, residual);
673 }
674 
675 template <typename SetResidualFunctor>
676 void
677 TaggingInterface::setResidual(SystemBase & sys, const SetResidualFunctor set_residual_functor)
678 {
679  for (const auto tag_id : _vector_tags)
680  if (sys.hasVector(tag_id))
681  set_residual_functor(sys.getVector(tag_id));
682 }
Nonlinear eigenvalue system to be solved.
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
Definition: EigenADReal.h:50
void accumulateTaggedNonlocalMatrix()
Nonlocal Jacobian blocks will be appended by adding the current nonlocal kernel Jacobian.
std::set< TagID > _ref_abs_vector_tags
A set of either size 1 or 0.
SubProblem & _subproblem
SubProblem that contains tag info.
void addResidualsAndJacobian(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...
const InputParameters & _tag_params
Parameters from moose object.
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
bool hasVector(const std::string &tag_name) const
Check if the named vector exists in the system.
Definition: SystemBase.C:925
bool hasVectorTags() const
Keeps track of stuff related to assembling.
Definition: Assembly.h:109
unsigned int TagID
Definition: MooseTypes.h:238
const std::vector< dof_id_type > & dof_indices
void assignTaggedLocalMatrix()
Local Jacobian blocks will assigned as the current local kernel Jacobian.
Class that is used as a parameter to some of our matrix tag APIs that allows only blessed framework c...
MatrixTagsKey(const MatrixTagsKey &)
VectorTagsKey(const VectorTagsKey &)
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.
void assignTaggedLocalResidual()
Local residual blocks will assigned as the current local kernel residual.
auto raw_value(const Eigen::Map< T > &in)
Definition: EigenADReal.h:100
std::vector< DenseVector< Number > * > _absre_blocks
Residual blocks for absolute value residual tags.
void insertNodalValue(libMesh::NumericVector< libMesh::Number > &residual, const DofValue &v)
Write a nodal value to the passed-in solution vector.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
The base class for Kokkos residual objects.
const Real scaling_factor
ResidualTagType
Enumerate whether a (residual) vector tag is to be of a non-reference or reference tag type...
std::vector< DenseVector< Number > * > _re_blocks
Residual blocks Vectors For each Tag.
const DenseVector< ADReal > & residuals
std::set< TagID > _ref_vector_tags
A set of either size 1 or 0.
Base class for a system (of equations)
Definition: SystemBase.h:85
The Kokkos system class.
Definition: KokkosSystem.h:34
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
std::set< TagID > _abs_vector_tags
The absolute value residual tag ids.
void cacheJacobian(GlobalDataKey)
Takes the values that are currently in _sub_Kee and appends them to the cached values.
Definition: Assembly.C:4046
Nonlinear system to be solved.
void prepareMatrixTagNeighbor(Assembly &assembly, unsigned int ivar, unsigned int jvar, Moose::DGJacobianType type)
Prepare data for computing element jacobian according to the active tags for DG and interface kernels...
TaggingInterface(const MooseObject *moose_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.
InputParameters validParams()
std::set< TagID > _non_ref_abs_vector_tags
A set to hold absolute value vector tags excluding the reference residual tag.
void prepareMatrixTagNonlocal(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing nonlocal element jacobian according to the active tags.
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)
Every object that can be built by the factory should be derived from this class.
Definition: MooseObject.h:28
std::set< TagID > _matrix_tags
The matrices this Kernel will contribute to.
void cacheJacobianWithoutConstraints(const Residuals &residuals, const Indices &row_indices, Real scaling_factor, LocalDataKey, const std::set< TagID > &matrix_tags)
Process the derivatives() data of a vector of ADReals.
Definition: Assembly.h:3178
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...
std::set< TagID > _non_ref_vector_tags
A set to hold vector tags excluding the reference residual tag.
std::set< TagID > _vector_tags
The vector tag ids this Kernel will contribute to.
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
Utility structure for packaging up all of the residual object&#39;s information needed to add into the sy...
A storage container for MooseObjects that inherit from SetupInterface.
void cacheResiduals(const Residuals &residuals, const Indices &row_indices, Real scaling_factor, LocalDataKey, const std::set< TagID > &vector_tags)
Process the supplied residual values.
Definition: Assembly.h:3058
void useMatrixTag(const TagName &tag_name, MatrixTagsKey)
DGJacobianType
Definition: MooseTypes.h:798
void prepareVectorTagNeighbor(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual the according to active tags for DG and interface kernels...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
This is the common base class for objects that give residual contributions.
Generic class for solving transient nonlinear problems.
Definition: SubProblem.h:78
This is the common base class for objects that give contributions to a linear system.
ConstraintJacobianType
Definition: MooseTypes.h:845
void prepareMatrixTagLower(Assembly &assembly, unsigned int ivar, unsigned int jvar, Moose::ConstraintJacobianType type)
Prepare data for computing the jacobian according to the active tags for mortar.
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
const MooseObject & _moose_object
Moose objct this tag works on.
const std::set< TagID > & getMatrixTags(MatrixTagsKey) const
IntRange< T > make_range(T beg, T end)
void addResidualsWithoutConstraints(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
Add the provided incoming residuals corresponding to the provided dof indices.
virtual ~TaggingInterface()
virtual void set(const numeric_index_type i, const Number value)=0
void checkForNans() const
Checks _local_re for NaNs/Infs and returns an error if found.
void cacheJacobianBlock(DenseMatrix< Number > &jac_block, const std::vector< dof_id_type > &idof_indices, const std::vector< dof_id_type > &jdof_indices, Real scaling_factor, LocalDataKey, TagID tag)
Cache a local Jacobian block with the provided rows (idof_indices) and columns (jdof_indices) for eve...
std::vector< DenseMatrix< Number > * > _ke_blocks
Kernel blocks Vectors For each Tag.
void prepareVectorTagLower(Assembly &assembly, unsigned int ivar)
Prepare data for computing the residual according to active tags for mortar constraints.
void prepareVectorTagInternal(Assembly &assembly, unsigned int ivar, const std::set< TagID > &vector_tags, const std::set< TagID > &absolute_value_vector_tags)
Prepare data for computing element residual according to the specified tags Residual blocks for diffe...
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual according to active tags.
void addJacobianWithoutConstraints(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.
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing element jacobian according to the active tags.
void addJacobianElement(Assembly &assembly, Real value, dof_id_type row_index, dof_id_type column_index, Real scaling_factor)
Add into a single Jacobian element.
DenseMatrix< Number > _nonlocal_ke
Holds nonlocal Jacobian entries as they are accumulated by this Kernel.
static InputParameters validParams()
void setResidual(SystemBase &sys, const T &residual, MooseVariableFE< T > &var)
Set residual using the variables&#39; insertion API.
auto makeSpan(C &container, std::size_t offset, std::size_t n)
Helper function for creating a span from a given container.
Definition: MooseTypes.h:1090
virtual NumericVector< Number > & getVector(const std::string &name)
Get a raw NumericVector by name.
Definition: SystemBase.C:934
auto index_range(const T &sizable)
std::vector< Real > _absolute_residuals
A container to hold absolute values of residuals passed into addResiduals.
Key structure for APIs adding/caching local element residuals/Jacobians.
Definition: Assembly.h:862
void useVectorTag(const TagName &tag_name, VectorTagsKey)
Class that is used as a parameter to some of our vector tag APIs that allows only blessed framework c...
void cacheResidualsWithoutConstraints(const Residuals &residuals, const Indices &row_indices, Real scaling_factor, LocalDataKey, const std::set< TagID > &vector_tags)
Process the supplied residual values.
Definition: Assembly.h:3098
uint8_t dof_id_type
const std::set< TagID > & getVectorTags(VectorTagsKey) const