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 286584 : 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 248391 : 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 489559326 : bool hasVectorTags() const { return !_vector_tags.empty(); }
131 :
132 283251 : const std::set<TagID> & getVectorTags(VectorTagsKey) const { return _vector_tags; }
133 :
134 242998 : 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 : /**
270 : * Add the provided incoming residuals corresponding to the provided dof indices
271 : */
272 : template <typename T, typename Indices>
273 : void addResiduals(Assembly & assembly,
274 : const DenseVector<T> & residuals,
275 : const Indices & dof_indices,
276 : Real scaling_factor);
277 :
278 : /**
279 : * Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the
280 : * provided dof indices
281 : */
282 : template <typename Residuals, typename Indices>
283 : void addResidualsAndJacobian(Assembly & assembly,
284 : const Residuals & residuals,
285 : const Indices & dof_indices,
286 : Real scaling_factor);
287 :
288 : /**
289 : * Add the provided residual derivatives into the Jacobian for the provided dof indices
290 : */
291 : template <typename Residuals, typename Indices>
292 : void addJacobian(Assembly & assembly,
293 : const Residuals & residuals,
294 : const Indices & dof_indices,
295 : Real scaling_factor);
296 :
297 : /**
298 : * Add the provided incoming residuals corresponding to the provided dof indices
299 : */
300 : void addResiduals(Assembly & assembly, const ADResidualsPacket & packet);
301 :
302 : /**
303 : * Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the
304 : * provided dof indices
305 : */
306 : void addResidualsAndJacobian(Assembly & assembly, const ADResidualsPacket & packet);
307 :
308 : /**
309 : * Add the provided residual derivatives into the Jacobian for the provided dof indices
310 : */
311 : void addJacobian(Assembly & assembly, const ADResidualsPacket & packet);
312 :
313 : /**
314 : * Add the provided incoming residuals corresponding to the provided dof indices. This API should
315 : * only be used if the caller knows that no libMesh-level constraints (hanging nodes or periodic
316 : * boundary conditions) apply to the provided dof indices
317 : */
318 : template <typename Residuals, typename Indices>
319 : void addResidualsWithoutConstraints(Assembly & assembly,
320 : const Residuals & residuals,
321 : const Indices & dof_indices,
322 : Real scaling_factor);
323 :
324 : /**
325 : * Add the provided incoming residuals and derivatives for the Jacobian, corresponding to the
326 : * provided dof indices. This API should only be used if the caller knows that no libMesh-level
327 : * constraints (hanging nodes or periodic boundary conditions) apply to the provided dof indices
328 : */
329 : template <typename Residuals, typename Indices>
330 : void addResidualsAndJacobianWithoutConstraints(Assembly & assembly,
331 : const Residuals & residuals,
332 : const Indices & dof_indices,
333 : Real scaling_factor);
334 :
335 : /**
336 : * Add the provided residual derivatives into the Jacobian for the provided dof indices. This API
337 : * should only be used if the caller knows that no libMesh-level constraints (hanging nodes or
338 : * periodic boundary conditions) apply to the provided dof indices
339 : */
340 : template <typename Residuals, typename Indices>
341 : void addJacobianWithoutConstraints(Assembly & assembly,
342 : const Residuals & residuals,
343 : const Indices & dof_indices,
344 : Real scaling_factor);
345 :
346 : /**
347 : * Add into a single Jacobian element
348 : */
349 : void addJacobianElement(Assembly & assembly,
350 : Real value,
351 : dof_id_type row_index,
352 : dof_id_type column_index,
353 : Real scaling_factor);
354 :
355 : /**
356 : * Add a local Jacobian matrix
357 : */
358 : void addJacobian(Assembly & assembly,
359 : DenseMatrix<Real> & local_k,
360 : const std::vector<dof_id_type> & row_indices,
361 : const std::vector<dof_id_type> & column_indices,
362 : Real scaling_factor);
363 :
364 : /**
365 : * Set residual using the variables' insertion API
366 : */
367 : template <typename T>
368 : void setResidual(SystemBase & sys, const T & residual, MooseVariableFE<T> & var);
369 :
370 : /**
371 : * Set residual at a specified degree of freedom index
372 : */
373 : void setResidual(SystemBase & sys, Real residual, dof_id_type dof_index);
374 :
375 : /**
376 : * Set residuals using the provided functor
377 : */
378 : template <typename SetResidualFunctor>
379 : void setResidual(SystemBase & sys, SetResidualFunctor set_residual_functor);
380 :
381 : /// SubProblem that contains tag info
382 : SubProblem & _subproblem;
383 :
384 : /// Holds local residual entries as they are accumulated by this Kernel
385 : DenseVector<Number> _local_re;
386 :
387 : /// Holds local Jacobian entries as they are accumulated by this Kernel
388 : DenseMatrix<Number> _local_ke;
389 :
390 : /// Holds nonlocal Jacobian entries as they are accumulated by this Kernel
391 : DenseMatrix<Number> _nonlocal_ke;
392 :
393 : private:
394 : /**
395 : * Prepare data for computing element residual according to the specified tags
396 : * Residual blocks for different tags will be extracted from Assembly.
397 : * A local residual will be zeroed. It should be called
398 : * right before the local element vector is computed.
399 : */
400 : void prepareVectorTagInternal(Assembly & assembly,
401 : unsigned int ivar,
402 : const std::set<TagID> & vector_tags,
403 : const std::set<TagID> & absolute_value_vector_tags);
404 :
405 : /// The vector tag ids this Kernel will contribute to
406 : std::set<TagID> _vector_tags;
407 :
408 : /// The absolute value residual tag ids
409 : std::set<TagID> _abs_vector_tags;
410 :
411 : /// The matrices this Kernel will contribute to
412 : std::set<TagID> _matrix_tags;
413 :
414 : /// A set to hold vector tags excluding the reference residual tag. If there is no reference
415 : /// residual problem, this container is the same as \p _vector_tags;
416 : std::set<TagID> _non_ref_vector_tags;
417 :
418 : /// A set to hold absolute value vector tags excluding the reference residual tag. If there is no
419 : /// reference residual problem, this container is the same as \p _abs_vector_tags;
420 : std::set<TagID> _non_ref_abs_vector_tags;
421 :
422 : /// A set of either size 1 or 0. If we have a reference residual problem and \p _vector_tags holds
423 : /// the reference vector tag, then this set holds the reference vector tags, otherwise it holds
424 : /// nothing
425 : std::set<TagID> _ref_vector_tags;
426 :
427 : /// A set of either size 1 or 0. If we have a reference residual problem and \p _abs_vector_tags
428 : /// holds the reference vector tag, then this set holds the reference vector tags, otherwise it
429 : /// holds nothing
430 : std::set<TagID> _ref_abs_vector_tags;
431 :
432 : /// Moose objct this tag works on
433 : const MooseObject & _moose_object;
434 :
435 : /// Parameters from moose object
436 : const InputParameters & _tag_params;
437 :
438 : /// Residual blocks Vectors For each Tag
439 : std::vector<DenseVector<Number> *> _re_blocks;
440 :
441 : /// Residual blocks for absolute value residual tags
442 : std::vector<DenseVector<Number> *> _absre_blocks;
443 :
444 : /// Kernel blocks Vectors For each Tag
445 : std::vector<DenseMatrix<Number> *> _ke_blocks;
446 :
447 : /// A container to hold absolute values of residuals passed into \p addResiduals. We maintain
448 : /// this data member to avoid constant dynamic heap allocations
449 : std::vector<Real> _absolute_residuals;
450 :
451 : friend class NonlinearSystemBase;
452 : };
453 :
454 : #define usingTaggingInterfaceMembers \
455 : using TaggingInterface::_subproblem; \
456 : using TaggingInterface::accumulateTaggedLocalResidual; \
457 : using TaggingInterface::accumulateTaggedLocalMatrix; \
458 : using TaggingInterface::prepareVectorTag; \
459 : using TaggingInterface::prepareMatrixTag; \
460 : using TaggingInterface::prepareVectorTagNeighbor; \
461 : using TaggingInterface::_local_re; \
462 : using TaggingInterface::prepareVectorTagLower; \
463 : using TaggingInterface::prepareMatrixTagNeighbor; \
464 : using TaggingInterface::prepareMatrixTagLower; \
465 : using TaggingInterface::_local_ke
466 :
467 : template <typename Residuals, typename Indices>
468 : void
469 43495012 : TaggingInterface::addResiduals(Assembly & assembly,
470 : const Residuals & residuals,
471 : const Indices & dof_indices,
472 : const Real scaling_factor)
473 : {
474 43495012 : assembly.cacheResiduals(
475 43495012 : residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _vector_tags);
476 43495012 : if (!_abs_vector_tags.empty())
477 : {
478 16380 : _absolute_residuals.resize(residuals.size());
479 49140 : for (const auto i : index_range(residuals))
480 32760 : _absolute_residuals[i] = std::abs(MetaPhysicL::raw_value(residuals[i]));
481 :
482 16380 : assembly.cacheResiduals(_absolute_residuals,
483 : dof_indices,
484 : scaling_factor,
485 16380 : Assembly::LocalDataKey{},
486 16380 : _abs_vector_tags);
487 : }
488 43495012 : }
489 :
490 : template <typename T, typename Indices>
491 : void
492 11631784 : TaggingInterface::addResiduals(Assembly & assembly,
493 : const DenseVector<T> & residuals,
494 : const Indices & dof_indices,
495 : const Real scaling_factor)
496 : {
497 11631784 : addResiduals(assembly, residuals.get_values(), dof_indices, scaling_factor);
498 11631784 : }
499 :
500 : template <typename Residuals, typename Indices>
501 : void
502 175103 : TaggingInterface::addResidualsWithoutConstraints(Assembly & assembly,
503 : const Residuals & residuals,
504 : const Indices & dof_indices,
505 : const Real scaling_factor)
506 : {
507 175103 : assembly.cacheResidualsWithoutConstraints(
508 175103 : residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _vector_tags);
509 175103 : if (!_abs_vector_tags.empty())
510 : {
511 216 : _absolute_residuals.resize(residuals.size());
512 936 : for (const auto i : index_range(residuals))
513 720 : _absolute_residuals[i] = std::abs(MetaPhysicL::raw_value(residuals[i]));
514 :
515 216 : assembly.cacheResidualsWithoutConstraints(_absolute_residuals,
516 : dof_indices,
517 : scaling_factor,
518 216 : Assembly::LocalDataKey{},
519 216 : _abs_vector_tags);
520 : }
521 175103 : }
522 :
523 : template <typename Residuals, typename Indices>
524 : void
525 15034901 : TaggingInterface::addResidualsAndJacobian(Assembly & assembly,
526 : const Residuals & residuals,
527 : const Indices & dof_indices,
528 : Real scaling_factor)
529 : {
530 15034901 : addResiduals(assembly, residuals, dof_indices, scaling_factor);
531 15034901 : addJacobian(assembly, residuals, dof_indices, scaling_factor);
532 15034901 : }
533 :
534 : template <typename Residuals, typename Indices>
535 : void
536 25776217 : TaggingInterface::addJacobian(Assembly & assembly,
537 : const Residuals & residuals,
538 : const Indices & dof_indices,
539 : Real scaling_factor)
540 : {
541 25776217 : assembly.cacheJacobian(
542 25776217 : residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _matrix_tags);
543 25776217 : }
544 :
545 : template <typename Residuals, typename Indices>
546 : void
547 175103 : TaggingInterface::addResidualsAndJacobianWithoutConstraints(Assembly & assembly,
548 : const Residuals & residuals,
549 : const Indices & dof_indices,
550 : Real scaling_factor)
551 : {
552 175103 : addResidualsWithoutConstraints(assembly, residuals, dof_indices, scaling_factor);
553 175103 : addJacobianWithoutConstraints(assembly, residuals, dof_indices, scaling_factor);
554 175103 : }
555 :
556 : template <typename Residuals, typename Indices>
557 : void
558 175103 : TaggingInterface::addJacobianWithoutConstraints(Assembly & assembly,
559 : const Residuals & residuals,
560 : const Indices & dof_indices,
561 : Real scaling_factor)
562 : {
563 175103 : assembly.cacheJacobianWithoutConstraints(
564 175103 : residuals, dof_indices, scaling_factor, Assembly::LocalDataKey{}, _matrix_tags);
565 175103 : }
566 :
567 : inline void
568 10034546 : TaggingInterface::addJacobianElement(Assembly & assembly,
569 : const Real value,
570 : const dof_id_type row_index,
571 : const dof_id_type column_index,
572 : const Real scaling_factor)
573 : {
574 10034546 : assembly.cacheJacobian(
575 10034546 : row_index, column_index, value * scaling_factor, Assembly::LocalDataKey{}, _matrix_tags);
576 10034546 : }
577 :
578 : inline void
579 6719819 : TaggingInterface::addJacobian(Assembly & assembly,
580 : DenseMatrix<Real> & local_k,
581 : const std::vector<dof_id_type> & row_indices,
582 : const std::vector<dof_id_type> & column_indices,
583 : const Real scaling_factor)
584 : {
585 13439638 : for (const auto matrix_tag : _matrix_tags)
586 6719819 : assembly.cacheJacobianBlock(
587 13439638 : local_k, row_indices, column_indices, scaling_factor, Assembly::LocalDataKey{}, matrix_tag);
588 6719819 : }
589 :
590 : template <typename T>
591 : void
592 61502232 : TaggingInterface::setResidual(SystemBase & sys, const T & residual, MooseVariableFE<T> & var)
593 : {
594 123845715 : for (const auto tag_id : _vector_tags)
595 62343483 : if (sys.hasVector(tag_id))
596 61550090 : var.insertNodalValue(sys.getVector(tag_id), residual);
597 61502232 : }
598 :
599 : inline void
600 918313 : TaggingInterface::setResidual(SystemBase & sys, const Real residual, const dof_id_type dof_index)
601 : {
602 1836626 : for (const auto tag_id : _vector_tags)
603 918313 : if (sys.hasVector(tag_id))
604 918313 : sys.getVector(tag_id).set(dof_index, residual);
605 918313 : }
606 :
607 : template <typename SetResidualFunctor>
608 : void
609 : TaggingInterface::setResidual(SystemBase & sys, const SetResidualFunctor set_residual_functor)
610 : {
611 : for (const auto tag_id : _vector_tags)
612 : if (sys.hasVector(tag_id))
613 : set_residual_functor(sys.getVector(tag_id));
614 : }
615 :
616 : inline void
617 6861552 : TaggingInterface::addResiduals(Assembly & assembly, const ADResidualsPacket & packet)
618 : {
619 6861552 : addResiduals(assembly, packet.residuals, packet.dof_indices, packet.scaling_factor);
620 6861552 : }
621 :
622 : inline void
623 0 : TaggingInterface::addResidualsAndJacobian(Assembly & assembly, const ADResidualsPacket & packet)
624 : {
625 0 : addResidualsAndJacobian(assembly, packet.residuals, packet.dof_indices, packet.scaling_factor);
626 0 : }
627 :
628 : inline void
629 3372130 : TaggingInterface::addJacobian(Assembly & assembly, const ADResidualsPacket & packet)
630 : {
631 3372130 : addJacobian(assembly, packet.residuals, packet.dof_indices, packet.scaling_factor);
632 3372130 : }
|