LCOV - code coverage report
Current view: top level - include/interfaces - TaggingInterface.h (source / functions) Hit Total Coverage
Test: idaholab/moose framework: 7323e9 Lines: 79 82 96.3 %
Date: 2025-11-05 20:01:15 Functions: 41 44 93.2 %
Legend: Lines: hit not hit

          Line data    Source code
       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             : }
      35             : #endif
      36             : 
      37             : template <typename T>
      38             : InputParameters validParams();
      39             : 
      40             : /**
      41             :  * Utility structure for packaging up all of the residual object's information needed to add into
      42             :  * the system residual and Jacobian in the context of automatic differentiation. The necessary
      43             :  * information includes the vector of residuals and Jacobians represented by dual numbers, the
      44             :  * vector of degrees of freedom (this should be the same length as the vector of dual numbers), and
      45             :  * a scaling factor that will multiply the dual numbers before their addition into the global
      46             :  * residual and Jacobian
      47             :  */
      48             : struct ADResidualsPacket
      49             : {
      50             :   const DenseVector<ADReal> & residuals;
      51             :   const std::vector<dof_id_type> & dof_indices;
      52             :   const Real scaling_factor;
      53             : };
      54             : 
      55             : class TaggingInterface
      56             : {
      57             : public:
      58             :   TaggingInterface(const MooseObject * moose_object);
      59             : 
      60             : #ifdef MOOSE_KOKKOS_ENABLED
      61             :   /**
      62             :    * Special constructor used for Kokkos functor copy during parallel dispatch
      63             :    */
      64             :   TaggingInterface(const TaggingInterface & object, const Moose::Kokkos::FunctorCopy & key);
      65             : #endif
      66             : 
      67             :   virtual ~TaggingInterface();
      68             : 
      69             :   static InputParameters validParams();
      70             : 
      71             :   /**
      72             :    * Class that is used as a parameter to some of our vector tag APIs that allows only
      73             :    * blessed framework classes to call them
      74             :    */
      75             :   class VectorTagsKey
      76             :   {
      77             :     friend class AttribVectorTags;
      78             :     friend class NonlinearEigenSystem;
      79             :     friend class LinearSystemContributionObject;
      80             :     template <typename>
      81             :     friend class MooseObjectTagWarehouse;
      82             : #ifdef MOOSE_KOKKOS_ENABLED
      83             :     friend class Moose::Kokkos::ResidualObject;
      84             : #endif
      85             : 
      86      287245 :     VectorTagsKey() {}
      87             :     VectorTagsKey(const VectorTagsKey &) {}
      88             :   };
      89             : 
      90             :   /**
      91             :    * Class that is used as a parameter to some of our matrix tag APIs that allows only
      92             :    * blessed framework classes to call them
      93             :    */
      94             :   class MatrixTagsKey
      95             :   {
      96             :     friend class AttribMatrixTags;
      97             :     friend class NonlinearEigenSystem;
      98             :     friend class LinearSystemContributionObject;
      99             :     template <typename>
     100             :     friend class MooseObjectTagWarehouse;
     101             : #ifdef MOOSE_KOKKOS_ENABLED
     102             :     friend class Moose::Kokkos::ResidualObject;
     103             : #endif
     104             : 
     105      249046 :     MatrixTagsKey() {}
     106             :     MatrixTagsKey(const MatrixTagsKey &) {}
     107             :   };
     108             : 
     109             :   /**
     110             :    * Enumerate whether a (residual) vector tag is to be of a non-reference or reference tag type
     111             :    */
     112             :   enum class ResidualTagType
     113             :   {
     114             :     NonReference,
     115             :     Reference
     116             :   };
     117             : 
     118             :   void useVectorTag(const TagName & tag_name, VectorTagsKey);
     119             : 
     120             :   void useMatrixTag(const TagName & tag_name, MatrixTagsKey);
     121             : 
     122             :   void useVectorTag(TagID tag_id, VectorTagsKey);
     123             : 
     124             :   void useMatrixTag(TagID tag_id, MatrixTagsKey);
     125             : 
     126             :   bool isVectorTagged() { return _vector_tags.size() > 0; }
     127             : 
     128             :   bool isMatrixTagged() { return _matrix_tags.size() > 0; }
     129             : 
     130   487210808 :   bool hasVectorTags() const { return !_vector_tags.empty(); }
     131             : 
     132      283912 :   const std::set<TagID> & getVectorTags(VectorTagsKey) const { return _vector_tags; }
     133             : 
     134      243653 :   const std::set<TagID> & getMatrixTags(MatrixTagsKey) const { return _matrix_tags; }
     135             : 
     136             : protected:
     137             :   /**
     138             :    * Prepare data for computing element residual according to active tags.
     139             :    * Residual blocks for different tags will be extracted from Assembly.
     140             :    * A local residual will be zeroed. It should be called
     141             :    * right before the local element vector is computed.
     142             :    */
     143             :   void prepareVectorTag(Assembly & assembly, unsigned int ivar);
     144             : 
     145             :   /**
     146             :    * Prepare vector tags in a reference residual problem context
     147             :    * @param Assembly The assembly object that we obtain the local residual blocks from
     148             :    * @param ivar The variable which we are retrieving the local residual blocks for
     149             :    * @param ref_problem A pointer to a reference residual problem. This can be a nullptr
     150             :    * @param tag_type What type of tags to prepare
     151             :    */
     152             :   void prepareVectorTag(Assembly & assembly, unsigned int ivar, ResidualTagType tag_type);
     153             : 
     154             :   /**
     155             :    * Prepare data for computing element residual the according to active tags
     156             :    * for DG and interface kernels.
     157             :    * Residual blocks for different tags will be extracted from Assembly.
     158             :    * A local residual will be zeroed. It should be called
     159             :    * right before the local element vector is computed.
     160             :    */
     161             :   void prepareVectorTagNeighbor(Assembly & assembly, unsigned int ivar);
     162             : 
     163             :   /**
     164             :    * Prepare data for computing the residual according to active tags for mortar constraints.
     165             :    * Residual blocks for different tags will be extracted from Assembly.  A local residual will be
     166             :    * zeroed. It should be called right before the local element vector is computed.
     167             :    */
     168             :   void prepareVectorTagLower(Assembly & assembly, unsigned int ivar);
     169             : 
     170             :   /**
     171             :    * Prepare data for computing element jacobian according to the active tags.
     172             :    * Jacobian blocks for different tags will be extracted from Assembly.
     173             :    * A local Jacobian will be zeroed. It should be called
     174             :    * right before the local element matrix is computed.
     175             :    */
     176             :   void prepareMatrixTag(Assembly & assembly, unsigned int ivar, unsigned int jvar);
     177             : 
     178             :   void prepareMatrixTag(Assembly & assembly,
     179             :                         unsigned int ivar,
     180             :                         unsigned int jvar,
     181             :                         DenseMatrix<Number> & k) const;
     182             : 
     183             :   /**
     184             :    * Prepare data for computing nonlocal element jacobian according to the active tags.
     185             :    * Jacobian blocks for different tags will be extracted from Assembly.
     186             :    * A nonlocal Jacobian will be zeroed. It should be called
     187             :    * right before the nonlocal element matrix is computed.
     188             :    */
     189             :   void prepareMatrixTagNonlocal(Assembly & assembly, unsigned int ivar, unsigned int jvar);
     190             : 
     191             :   /**
     192             :    * Prepare data for computing element jacobian according to the active tags
     193             :    * for DG and interface kernels.
     194             :    * Jacobian blocks for different tags will be extracted from Assembly.
     195             :    * A local Jacobian will be zeroed. It should be called
     196             :    * right before the local element matrix is computed.
     197             :    */
     198             :   void prepareMatrixTagNeighbor(Assembly & assembly,
     199             :                                 unsigned int ivar,
     200             :                                 unsigned int jvar,
     201             :                                 Moose::DGJacobianType type);
     202             : 
     203             :   void prepareMatrixTagNeighbor(Assembly & assembly,
     204             :                                 unsigned int ivar,
     205             :                                 unsigned int jvar,
     206             :                                 Moose::DGJacobianType type,
     207             :                                 DenseMatrix<Number> & k) const;
     208             : 
     209             :   /**
     210             :    * Prepare data for computing the jacobian according to the active tags for mortar.  Jacobian
     211             :    * blocks for different tags will be extracted from Assembly.  A local Jacobian will be zeroed. It
     212             :    * should be called right before the local element matrix is computed.
     213             :    */
     214             :   void prepareMatrixTagLower(Assembly & assembly,
     215             :                              unsigned int ivar,
     216             :                              unsigned int jvar,
     217             :                              Moose::ConstraintJacobianType type);
     218             : 
     219             :   /**
     220             :    * Local residual blocks  will be appended by adding the current local kernel residual.
     221             :    * It should be called after the local element vector has been computed.
     222             :    */
     223             :   void accumulateTaggedLocalResidual();
     224             : 
     225             :   /**
     226             :    * Local residual blocks will assigned as the current local kernel residual.
     227             :    * It should be called after the local element vector has been computed.
     228             :    */
     229             :   void assignTaggedLocalResidual();
     230             : 
     231             :   /**
     232             :    * Local Jacobian blocks  will be appended by adding the current local kernel Jacobian.
     233             :    * It should be called after the local element matrix has been computed.
     234             :    */
     235             :   void accumulateTaggedLocalMatrix();
     236             : 
     237             :   void accumulateTaggedLocalMatrix(Assembly & assembly,
     238             :                                    unsigned int ivar,
     239             :                                    unsigned int jvar,
     240             :                                    const DenseMatrix<Number> & k);
     241             : 
     242             :   void accumulateTaggedLocalMatrix(Assembly & assembly,
     243             :                                    unsigned int ivar,
     244             :                                    unsigned int jvar,
     245             :                                    Moose::DGJacobianType type,
     246             :                                    const DenseMatrix<Number> & k);
     247             : 
     248             :   /**
     249             :    * Nonlocal Jacobian blocks  will be appended by adding the current nonlocal kernel Jacobian.
     250             :    * It should be called after the nonlocal element matrix has been computed.
     251             :    */
     252             :   void accumulateTaggedNonlocalMatrix();
     253             : 
     254             :   /**
     255             :    * Local Jacobian blocks will assigned as the current local kernel Jacobian.
     256             :    * It should be called after the local element matrix has been computed.
     257             :    */
     258             :   void assignTaggedLocalMatrix();
     259             : 
     260             :   /**
     261             :    * Add the provided incoming residuals corresponding to the provided dof indices
     262             :    */
     263             :   template <typename Residuals, typename Indices>
     264             :   void addResiduals(Assembly & assembly,
     265             :                     const Residuals & residuals,
     266             :                     const Indices & dof_indices,
     267             :                     Real scaling_factor);
     268             : 
     269             :   template <typename Residuals, typename Indices>
     270             :   void addResiduals(Assembly & assembly,
     271             :                     const Residuals & residuals,
     272             :                     const Indices & dof_indices,
     273             :                     const std::vector<Real> & scaling_factors);
     274             : 
     275             :   /**
     276             :    * Add the provided incoming residuals corresponding to the provided dof indices
     277             :    */
     278             :   template <typename T, typename Indices>
     279             :   void addResiduals(Assembly & assembly,
     280             :                     const DenseVector<T> & residuals,
     281             :                     const Indices & dof_indices,
     282             :                     Real scaling_factor);
     283             : 
     284             :   /**
     285             :    * Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the
     286             :    * provided dof indices
     287             :    */
     288             :   template <typename Residuals, typename Indices>
     289             :   void addResidualsAndJacobian(Assembly & assembly,
     290             :                                const Residuals & residuals,
     291             :                                const Indices & dof_indices,
     292             :                                Real scaling_factor);
     293             : 
     294             :   /**
     295             :    * Add the provided residual derivatives into the Jacobian for the provided dof indices
     296             :    */
     297             :   template <typename Residuals, typename Indices>
     298             :   void addJacobian(Assembly & assembly,
     299             :                    const Residuals & residuals,
     300             :                    const Indices & dof_indices,
     301             :                    Real scaling_factor);
     302             : 
     303             :   /**
     304             :    * Add the provided incoming residuals corresponding to the provided dof indices
     305             :    */
     306             :   void addResiduals(Assembly & assembly, const ADResidualsPacket & packet);
     307             : 
     308             :   /**
     309             :    * Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the
     310             :    * provided dof indices
     311             :    */
     312             :   void addResidualsAndJacobian(Assembly & assembly, const ADResidualsPacket & packet);
     313             : 
     314             :   /**
     315             :    * Add the provided residual derivatives into the Jacobian for the provided dof indices
     316             :    */
     317             :   void addJacobian(Assembly & assembly, const ADResidualsPacket & packet);
     318             : 
     319             :   /**
     320             :    * Add the provided incoming residuals corresponding to the provided dof indices. This API should
     321             :    * only be used if the caller knows that no libMesh-level constraints (hanging nodes or periodic
     322             :    * boundary conditions) apply to the provided dof indices
     323             :    */
     324             :   template <typename Residuals, typename Indices>
     325             :   void addResidualsWithoutConstraints(Assembly & assembly,
     326             :                                       const Residuals & residuals,
     327             :                                       const Indices & dof_indices,
     328             :                                       Real scaling_factor);
     329             : 
     330             :   /**
     331             :    * Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the
     332             :    * provided dof indices.  This API should only be used if the caller knows that no libMesh-level
     333             :    * constraints (hanging nodes or periodic boundary conditions) apply to the provided dof indices
     334             :    */
     335             :   template <typename Residuals, typename Indices>
     336             :   void addResidualsAndJacobianWithoutConstraints(Assembly & assembly,
     337             :                                                  const Residuals & residuals,
     338             :                                                  const Indices & dof_indices,
     339             :                                                  Real scaling_factor);
     340             : 
     341             :   /**
     342             :    * Add the provided residual derivatives into the Jacobian for the provided dof indices. This API
     343             :    * should only be used if the caller knows that no libMesh-level constraints (hanging nodes or
     344             :    * periodic boundary conditions) apply to the provided dof indices
     345             :    */
     346             :   template <typename Residuals, typename Indices>
     347             :   void addJacobianWithoutConstraints(Assembly & assembly,
     348             :                                      const Residuals & residuals,
     349             :                                      const Indices & dof_indices,
     350             :                                      Real scaling_factor);
     351             : 
     352             :   /**
     353             :    * Add into a single Jacobian element
     354             :    */
     355             :   void addJacobianElement(Assembly & assembly,
     356             :                           Real value,
     357             :                           dof_id_type row_index,
     358             :                           dof_id_type column_index,
     359             :                           Real scaling_factor);
     360             : 
     361             :   /**
     362             :    * Add a local Jacobian matrix
     363             :    */
     364             :   void addJacobian(Assembly & assembly,
     365             :                    DenseMatrix<Real> & local_k,
     366             :                    const std::vector<dof_id_type> & row_indices,
     367             :                    const std::vector<dof_id_type> & column_indices,
     368             :                    Real scaling_factor);
     369             : 
     370             :   /**
     371             :    * Set residual using the variables' insertion API
     372             :    */
     373             :   template <typename T>
     374             :   void setResidual(SystemBase & sys, const T & residual, MooseVariableFE<T> & var);
     375             : 
     376             :   /**
     377             :    * Set residual at a specified degree of freedom index
     378             :    */
     379             :   void setResidual(SystemBase & sys, Real residual, dof_id_type dof_index);
     380             : 
     381             :   /**
     382             :    * Set residuals using the provided functor
     383             :    */
     384             :   template <typename SetResidualFunctor>
     385             :   void setResidual(SystemBase & sys, SetResidualFunctor set_residual_functor);
     386             : 
     387             :   /// SubProblem that contains tag info
     388             :   SubProblem & _subproblem;
     389             : 
     390             :   /// Holds local residual entries as they are accumulated by this Kernel
     391             :   DenseVector<Number> _local_re;
     392             : 
     393             :   /// Holds local Jacobian entries as they are accumulated by this Kernel
     394             :   DenseMatrix<Number> _local_ke;
     395             : 
     396             :   /// Holds nonlocal Jacobian entries as they are accumulated by this Kernel
     397             :   DenseMatrix<Number> _nonlocal_ke;
     398             : 
     399             : private:
     400             :   /**
     401             :    * Prepare data for computing element residual according to the specified tags
     402             :    * Residual blocks for different tags will be extracted from Assembly.
     403             :    * A local residual will be zeroed. It should be called
     404             :    * right before the local element vector is computed.
     405             :    */
     406             :   void prepareVectorTagInternal(Assembly & assembly,
     407             :                                 unsigned int ivar,
     408             :                                 const std::set<TagID> & vector_tags,
     409             :                                 const std::set<TagID> & absolute_value_vector_tags);
     410             : 
     411             :   /// The vector tag ids this Kernel will contribute to
     412             :   std::set<TagID> _vector_tags;
     413             : 
     414             :   /// The absolute value residual tag ids
     415             :   std::set<TagID> _abs_vector_tags;
     416             : 
     417             :   /// The matrices this Kernel will contribute to
     418             :   std::set<TagID> _matrix_tags;
     419             : 
     420             :   /// A set to hold vector tags excluding the reference residual tag. If there is no reference
     421             :   /// residual problem, this container is the same as \p _vector_tags;
     422             :   std::set<TagID> _non_ref_vector_tags;
     423             : 
     424             :   /// A set to hold absolute value vector tags excluding the reference residual tag. If there is no
     425             :   /// reference residual problem, this container is the same as \p _abs_vector_tags;
     426             :   std::set<TagID> _non_ref_abs_vector_tags;
     427             : 
     428             :   /// A set of either size 1 or 0. If we have a reference residual problem and \p _vector_tags holds
     429             :   /// the reference vector tag, then this set holds the reference vector tags, otherwise it holds
     430             :   /// nothing
     431             :   std::set<TagID> _ref_vector_tags;
     432             : 
     433             :   /// A set of either size 1 or 0. If we have a reference residual problem and \p _abs_vector_tags
     434             :   /// holds the reference vector tag, then this set holds the reference vector tags, otherwise it
     435             :   /// holds nothing
     436             :   std::set<TagID> _ref_abs_vector_tags;
     437             : 
     438             :   /// Moose objct this tag works on
     439             :   const MooseObject & _moose_object;
     440             : 
     441             :   /// Parameters from moose object
     442             :   const InputParameters & _tag_params;
     443             : 
     444             :   /// Residual blocks Vectors For each Tag
     445             :   std::vector<DenseVector<Number> *> _re_blocks;
     446             : 
     447             :   /// Residual blocks for absolute value residual tags
     448             :   std::vector<DenseVector<Number> *> _absre_blocks;
     449             : 
     450             :   /// Kernel blocks Vectors For each Tag
     451             :   std::vector<DenseMatrix<Number> *> _ke_blocks;
     452             : 
     453             :   /// A container to hold absolute values of residuals passed into \p addResiduals. We maintain
     454             :   /// this data member to avoid constant dynamic heap allocations
     455             :   std::vector<Real> _absolute_residuals;
     456             : 
     457             :   friend class NonlinearSystemBase;
     458             : };
     459             : 
     460             : #define usingTaggingInterfaceMembers                                                               \
     461             :   using TaggingInterface::_subproblem;                                                             \
     462             :   using TaggingInterface::accumulateTaggedLocalResidual;                                           \
     463             :   using TaggingInterface::accumulateTaggedLocalMatrix;                                             \
     464             :   using TaggingInterface::prepareVectorTag;                                                        \
     465             :   using TaggingInterface::prepareMatrixTag;                                                        \
     466             :   using TaggingInterface::prepareVectorTagNeighbor;                                                \
     467             :   using TaggingInterface::_local_re;                                                               \
     468             :   using TaggingInterface::prepareVectorTagLower;                                                   \
     469             :   using TaggingInterface::prepareMatrixTagNeighbor;                                                \
     470             :   using TaggingInterface::prepareMatrixTagLower;                                                   \
     471             :   using TaggingInterface::_local_ke
     472             : 
     473             : template <typename Residuals, typename Indices>
     474             : void
     475    43501432 : TaggingInterface::addResiduals(Assembly & assembly,
     476             :                                const Residuals & residuals,
     477             :                                const Indices & dof_indices,
     478             :                                const Real scaling_factor)
     479             : {
     480    43501432 :   assembly.cacheResiduals(
     481    43501432 :       residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _vector_tags);
     482    43501432 :   if (!_abs_vector_tags.empty())
     483             :   {
     484       16380 :     _absolute_residuals.resize(residuals.size());
     485       49140 :     for (const auto i : index_range(residuals))
     486       32760 :       _absolute_residuals[i] = std::abs(MetaPhysicL::raw_value(residuals[i]));
     487             : 
     488       16380 :     assembly.cacheResiduals(_absolute_residuals,
     489             :                             dof_indices,
     490             :                             scaling_factor,
     491       16380 :                             Assembly::LocalDataKey{},
     492       16380 :                             _abs_vector_tags);
     493             :   }
     494    43501432 : }
     495             : 
     496             : template <typename Residuals, typename Indices>
     497             : void
     498        1464 : TaggingInterface::addResiduals(Assembly & assembly,
     499             :                                const Residuals & residuals,
     500             :                                const Indices & dof_indices,
     501             :                                const std::vector<Real> & scaling_factors)
     502             : {
     503        1464 :   const auto n = libMesh::cast_int<std::size_t>(std::size(residuals));
     504             : 
     505             :   mooseAssert((n == libMesh::cast_int<std::size_t>(std::size(dof_indices))) &&
     506             :                   (n == scaling_factors.size()),
     507             :               "Residuals, dof indices, and scaling factors must all be the same size");
     508             : 
     509             :   using R = Moose::ContainerElement<Residuals>;
     510             :   using I = Moose::ContainerElement<Indices>;
     511             : 
     512        4392 :   for (const auto i : make_range(n))
     513             :     // The Residuals type may not offer operator[] (e.g. eigen vectors) but more commonly it should
     514             :     // offer data()
     515        5856 :     addResiduals(assembly,
     516        2928 :                  std::array<R, 1>{std::data(residuals)[i]},
     517        2928 :                  std::array<I, 1>{std::data(dof_indices)[i]},
     518        2928 :                  scaling_factors[i]);
     519        1464 : }
     520             : 
     521             : template <typename T, typename Indices>
     522             : void
     523    11631784 : TaggingInterface::addResiduals(Assembly & assembly,
     524             :                                const DenseVector<T> & residuals,
     525             :                                const Indices & dof_indices,
     526             :                                const Real scaling_factor)
     527             : {
     528    11631784 :   addResiduals(assembly, residuals.get_values(), dof_indices, scaling_factor);
     529    11631784 : }
     530             : 
     531             : template <typename Residuals, typename Indices>
     532             : void
     533      175175 : TaggingInterface::addResidualsWithoutConstraints(Assembly & assembly,
     534             :                                                  const Residuals & residuals,
     535             :                                                  const Indices & dof_indices,
     536             :                                                  const Real scaling_factor)
     537             : {
     538      175175 :   assembly.cacheResidualsWithoutConstraints(
     539      175175 :       residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _vector_tags);
     540      175175 :   if (!_abs_vector_tags.empty())
     541             :   {
     542         216 :     _absolute_residuals.resize(residuals.size());
     543         936 :     for (const auto i : index_range(residuals))
     544         720 :       _absolute_residuals[i] = std::abs(MetaPhysicL::raw_value(residuals[i]));
     545             : 
     546         216 :     assembly.cacheResidualsWithoutConstraints(_absolute_residuals,
     547             :                                               dof_indices,
     548             :                                               scaling_factor,
     549         216 :                                               Assembly::LocalDataKey{},
     550         216 :                                               _abs_vector_tags);
     551             :   }
     552      175175 : }
     553             : 
     554             : template <typename Residuals, typename Indices>
     555             : void
     556    15034901 : TaggingInterface::addResidualsAndJacobian(Assembly & assembly,
     557             :                                           const Residuals & residuals,
     558             :                                           const Indices & dof_indices,
     559             :                                           Real scaling_factor)
     560             : {
     561    15034901 :   addResiduals(assembly, residuals, dof_indices, scaling_factor);
     562    15034901 :   addJacobian(assembly, residuals, dof_indices, scaling_factor);
     563    15034901 : }
     564             : 
     565             : template <typename Residuals, typename Indices>
     566             : void
     567    25781968 : TaggingInterface::addJacobian(Assembly & assembly,
     568             :                               const Residuals & residuals,
     569             :                               const Indices & dof_indices,
     570             :                               Real scaling_factor)
     571             : {
     572    25781968 :   assembly.cacheJacobian(
     573    25781968 :       residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _matrix_tags);
     574    25781968 : }
     575             : 
     576             : template <typename Residuals, typename Indices>
     577             : void
     578      175175 : TaggingInterface::addResidualsAndJacobianWithoutConstraints(Assembly & assembly,
     579             :                                                             const Residuals & residuals,
     580             :                                                             const Indices & dof_indices,
     581             :                                                             Real scaling_factor)
     582             : {
     583      175175 :   addResidualsWithoutConstraints(assembly, residuals, dof_indices, scaling_factor);
     584      175175 :   addJacobianWithoutConstraints(assembly, residuals, dof_indices, scaling_factor);
     585      175175 : }
     586             : 
     587             : template <typename Residuals, typename Indices>
     588             : void
     589      175175 : TaggingInterface::addJacobianWithoutConstraints(Assembly & assembly,
     590             :                                                 const Residuals & residuals,
     591             :                                                 const Indices & dof_indices,
     592             :                                                 Real scaling_factor)
     593             : {
     594      175175 :   assembly.cacheJacobianWithoutConstraints(
     595      175175 :       residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _matrix_tags);
     596      175175 : }
     597             : 
     598             : inline void
     599    10046675 : TaggingInterface::addJacobianElement(Assembly & assembly,
     600             :                                      const Real value,
     601             :                                      const dof_id_type row_index,
     602             :                                      const dof_id_type column_index,
     603             :                                      const Real scaling_factor)
     604             : {
     605    10046675 :   assembly.cacheJacobian(
     606    10046675 :       row_index, column_index, value * scaling_factor, Assembly::LocalDataKey{}, _matrix_tags);
     607    10046675 : }
     608             : 
     609             : inline void
     610     6719819 : TaggingInterface::addJacobian(Assembly & assembly,
     611             :                               DenseMatrix<Real> & local_k,
     612             :                               const std::vector<dof_id_type> & row_indices,
     613             :                               const std::vector<dof_id_type> & column_indices,
     614             :                               const Real scaling_factor)
     615             : {
     616    13439638 :   for (const auto matrix_tag : _matrix_tags)
     617     6719819 :     assembly.cacheJacobianBlock(
     618    13439638 :         local_k, row_indices, column_indices, scaling_factor, Assembly::LocalDataKey{}, matrix_tag);
     619     6719819 : }
     620             : 
     621             : template <typename T>
     622             : void
     623    61671475 : TaggingInterface::setResidual(SystemBase & sys, const T & residual, MooseVariableFE<T> & var)
     624             : {
     625   124184201 :   for (const auto tag_id : _vector_tags)
     626    62512726 :     if (sys.hasVector(tag_id))
     627    61719333 :       var.insertNodalValue(sys.getVector(tag_id), residual);
     628    61671475 : }
     629             : 
     630             : inline void
     631      921934 : TaggingInterface::setResidual(SystemBase & sys, const Real residual, const dof_id_type dof_index)
     632             : {
     633     1843868 :   for (const auto tag_id : _vector_tags)
     634      921934 :     if (sys.hasVector(tag_id))
     635      921934 :       sys.getVector(tag_id).set(dof_index, residual);
     636      921934 : }
     637             : 
     638             : template <typename SetResidualFunctor>
     639             : void
     640             : TaggingInterface::setResidual(SystemBase & sys, const SetResidualFunctor set_residual_functor)
     641             : {
     642             :   for (const auto tag_id : _vector_tags)
     643             :     if (sys.hasVector(tag_id))
     644             :       set_residual_functor(sys.getVector(tag_id));
     645             : }
     646             : 
     647             : inline void
     648     6861552 : TaggingInterface::addResiduals(Assembly & assembly, const ADResidualsPacket & packet)
     649             : {
     650     6861552 :   addResiduals(assembly, packet.residuals, packet.dof_indices, packet.scaling_factor);
     651     6861552 : }
     652             : 
     653             : inline void
     654           0 : TaggingInterface::addResidualsAndJacobian(Assembly & assembly, const ADResidualsPacket & packet)
     655             : {
     656           0 :   addResidualsAndJacobian(assembly, packet.residuals, packet.dof_indices, packet.scaling_factor);
     657           0 : }
     658             : 
     659             : inline void
     660     3372130 : TaggingInterface::addJacobian(Assembly & assembly, const ADResidualsPacket & packet)
     661             : {
     662     3372130 :   addJacobian(assembly, packet.residuals, packet.dof_indices, packet.scaling_factor);
     663     3372130 : }

Generated by: LCOV version 1.14