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 "MooseArray.h"
13 : #include "MooseTypes.h"
14 : #include "MooseVariableDataBase.h"
15 : #include "Conversion.h"
16 :
17 : #include "libmesh/tensor_tools.h"
18 : #include "libmesh/vector_value.h"
19 : #include "libmesh/tensor_value.h"
20 : #include "libmesh/type_n_tensor.h"
21 : #include "libmesh/enum_fe_family.h"
22 : #include "libmesh/fe_type.h"
23 : #include "ADUtils.h"
24 :
25 : #include <functional>
26 : #include <vector>
27 :
28 : namespace libMesh
29 : {
30 : class QBase;
31 : class DofMap;
32 : template <typename>
33 : class DenseVector;
34 : }
35 :
36 : using libMesh::DenseVector;
37 : using libMesh::DofMap;
38 : using libMesh::QBase;
39 : using libMesh::VectorValue;
40 :
41 : class TimeIntegrator;
42 : class Assembly;
43 : class SubProblem;
44 : template <typename>
45 : class MooseVariableFE;
46 : class SystemBase;
47 :
48 : template <typename OutputType>
49 : class MooseVariableData : public MooseVariableDataBase<OutputType>
50 : {
51 : public:
52 : // type for gradient, second and divergence of template class OutputType
53 : using typename MooseVariableDataBase<OutputType>::OutputGradient;
54 : using typename MooseVariableDataBase<OutputType>::OutputSecond;
55 : using typename MooseVariableDataBase<OutputType>::OutputDivergence;
56 :
57 : // shortcut for types storing values on quadrature points
58 : using typename MooseVariableDataBase<OutputType>::FieldVariableValue;
59 : using typename MooseVariableDataBase<OutputType>::FieldVariableGradient;
60 : typedef MooseArray<OutputSecond> FieldVariableSecond;
61 : typedef MooseArray<OutputType> FieldVariableCurl;
62 : typedef MooseArray<OutputDivergence> FieldVariableDivergence;
63 :
64 : // shape function type for the template class OutputType
65 : typedef typename Moose::ShapeType<OutputType>::type OutputShape;
66 :
67 : // type for gradient, second and divergence of shape functions of template class OutputType
68 : typedef typename libMesh::TensorTools::IncrementRank<OutputShape>::type OutputShapeGradient;
69 : typedef typename libMesh::TensorTools::IncrementRank<OutputShapeGradient>::type OutputShapeSecond;
70 : typedef typename libMesh::TensorTools::DecrementRank<OutputShape>::type OutputShapeDivergence;
71 :
72 : // shortcut for types storing shape function values on quadrature points
73 : typedef MooseArray<std::vector<OutputShape>> FieldVariablePhiValue;
74 : typedef MooseArray<std::vector<OutputShapeGradient>> FieldVariablePhiGradient;
75 : typedef MooseArray<std::vector<OutputShapeSecond>> FieldVariablePhiSecond;
76 : typedef MooseArray<std::vector<OutputShape>> FieldVariablePhiCurl;
77 : typedef MooseArray<std::vector<OutputShapeDivergence>> FieldVariablePhiDivergence;
78 :
79 : // shortcut for types storing test function values on quadrature points
80 : // Note: here we assume the types are the same as of shape functions.
81 : typedef MooseArray<std::vector<OutputShape>> FieldVariableTestValue;
82 : typedef MooseArray<std::vector<OutputShapeGradient>> FieldVariableTestGradient;
83 : typedef MooseArray<std::vector<OutputShapeSecond>> FieldVariableTestSecond;
84 : typedef MooseArray<std::vector<OutputShape>> FieldVariableTestCurl;
85 : typedef MooseArray<std::vector<OutputShapeDivergence>> FieldVariableTestDivergence;
86 :
87 : // DoF value type for the template class OutputType
88 : using typename MooseVariableDataBase<OutputType>::DofValue;
89 : using typename MooseVariableDataBase<OutputType>::DofValues;
90 : using typename MooseVariableDataBase<OutputType>::ADDofValue;
91 : using typename MooseVariableDataBase<OutputType>::ADDofValues;
92 :
93 : MooseVariableData(const MooseVariableFE<OutputType> & var,
94 : SystemBase & sys,
95 : THREAD_ID tid,
96 : Moose::ElementType element_type,
97 : const QBase * const & qrule_in,
98 : const QBase * const & qrule_face_in,
99 : const Node * const & node,
100 : const Elem * const & elem);
101 :
102 : /**
103 : * Returns whether this data structure needs automatic differentiation calculations
104 : */
105 0 : bool needsAD() const { return _need_ad; }
106 :
107 : /**
108 : * Set the geometry type before calculating variables values
109 : * @param gm_type The type type of geometry; either Volume or Face
110 : */
111 : void setGeometry(Moose::GeometryType gm_type);
112 :
113 : //////////////// Heavy lifting computational routines //////////////////////////////
114 :
115 : /**
116 : * compute the variable values
117 : */
118 : void computeValues();
119 :
120 : /**
121 : * compute the values for const monomial variables
122 : */
123 : void computeConstantMonomialValues();
124 :
125 : /**
126 : * compute AD things
127 : */
128 : template <bool constant_monomial>
129 : void computeAD(const unsigned int num_dofs, const unsigned int nqp);
130 :
131 : /**
132 : * compute nodal things
133 : */
134 : void computeNodalValues();
135 :
136 : ///////////////////////////// Shape functions /////////////////////////////////////
137 :
138 : /**
139 : * phi getter
140 : */
141 126769272 : const FieldVariablePhiValue & phi() const { return *_phi; }
142 :
143 : /**
144 : * phi_face getter
145 : */
146 1000098 : const FieldVariablePhiValue & phiFace() const { return *_phi_face; }
147 :
148 : /**
149 : * grad_phi getter
150 : */
151 126705330 : const FieldVariablePhiGradient & gradPhi() const { return *_grad_phi; }
152 :
153 : /**
154 : * mapped_grad_phi getter
155 : */
156 1910 : const MappedArrayVariablePhiGradient & arrayGradPhi() const
157 : {
158 : mooseAssert(var().fieldType() == Moose::VarFieldType::VAR_FIELD_ARRAY, "Not an array variable");
159 1910 : return _mapped_grad_phi;
160 : }
161 :
162 : /**
163 : * grad_phi_face getter
164 : */
165 970014 : const FieldVariablePhiGradient & gradPhiFace() const { return *_grad_phi_face; }
166 :
167 : /**
168 : * mapped_grad_phi_face getter
169 : */
170 712 : const MappedArrayVariablePhiGradient & arrayGradPhiFace() const
171 : {
172 : mooseAssert(var().fieldType() == Moose::VarFieldType::VAR_FIELD_ARRAY, "Not an array variable");
173 712 : return _mapped_grad_phi_face;
174 : }
175 :
176 : /**
177 : * second_phi getter
178 : */
179 : const FieldVariablePhiSecond & secondPhi() const;
180 :
181 : /**
182 : * second_phi_face getter
183 : */
184 : const FieldVariablePhiSecond & secondPhiFace() const;
185 :
186 : /**
187 : * curl_phi getter
188 : */
189 : const FieldVariablePhiCurl & curlPhi() const;
190 :
191 : /**
192 : * curl_phi_face getter
193 : */
194 : const FieldVariablePhiCurl & curlPhiFace() const;
195 :
196 : /**
197 : * divergence_phi getter
198 : */
199 : const FieldVariablePhiDivergence & divPhi() const;
200 :
201 : /**
202 : * divergence_phi_face getter
203 : */
204 : const FieldVariablePhiDivergence & divPhiFace() const;
205 :
206 : /**
207 : * ad_grad_phi getter
208 : */
209 : const ADTemplateVariablePhiGradient<OutputShape> & adGradPhi() const;
210 :
211 : /**
212 : * ad_grad_phi_face getter
213 : */
214 : const ADTemplateVariablePhiGradient<OutputShape> & adGradPhiFace() const;
215 :
216 : /**
217 : * Return phi size
218 : */
219 15920842 : std::size_t phiSize() const { return _phi->size(); }
220 :
221 : /**
222 : * Return phiFace size
223 : */
224 0 : std::size_t phiFaceSize() const { return _phi_face->size(); }
225 :
226 : /**
227 : * Whether or not this variable is computing any second derivatives.
228 : */
229 127107304 : bool usesSecondPhi() const
230 : {
231 127107304 : return _need_second || _need_second_old || _need_second_older || _need_second_previous_nl;
232 : }
233 :
234 : /**
235 : * Whether or not this variable is computing the curl
236 : */
237 1873054 : bool computingCurl() const { return _need_curl || _need_curl_old; }
238 :
239 : /**
240 : * Whether or not this variable is computing the divergence
241 : */
242 1873054 : bool computingDiv() const { return _need_div || _need_div_old; }
243 :
244 : //////////////////////////////// Nodal stuff ///////////////////////////////////////////
245 :
246 9293669 : bool isNodal() const override { return _is_nodal; }
247 121 : bool hasDoFsOnNodes() const override { return _continuity != libMesh::DISCONTINUOUS; }
248 198389 : libMesh::FEContinuity getContinuity() const override { return _continuity; };
249 66750 : const Node * const & node() const { return _node; }
250 14591781 : const dof_id_type & nodalDofIndex() const { return _nodal_dof_index; }
251 418254962 : bool isNodalDefined() const { return _has_dof_indices; }
252 :
253 : /**
254 : * The current element
255 : */
256 35286 : const Elem * const & currentElem() const { return _elem; }
257 :
258 : /**
259 : * The current side
260 : */
261 0 : const unsigned int & currentSide() const { return _current_side; }
262 :
263 : /**
264 : * prepare the initial condition
265 : */
266 : void prepareIC();
267 :
268 : //////////////////////////////////// Solution getters /////////////////////////////////////
269 :
270 : /**
271 : * Local time derivative of solution gradient getter
272 : */
273 : const FieldVariableGradient & gradSlnDot() const;
274 :
275 : /**
276 : * Local second time derivative of solution gradient getter
277 : */
278 : const FieldVariableGradient & gradSlnDotDot() const;
279 :
280 : /**
281 : * Local solution second spatial derivative getter
282 : * @param state The state of the simulation: current, old, older, previous nl
283 : */
284 : const FieldVariableSecond & secondSln(Moose::SolutionState state) const;
285 :
286 : /**
287 : * Local solution curl getter
288 : * @param state The state of the simulation: current, old, older
289 : */
290 : const FieldVariableCurl & curlSln(Moose::SolutionState state) const;
291 :
292 : /**
293 : * Local solution divergence getter
294 : * @param state The state of the simulation: current, old, older
295 : */
296 : const FieldVariableDivergence & divSln(Moose::SolutionState state) const;
297 :
298 14398 : const ADTemplateVariableValue<OutputType> & adSln() const
299 : {
300 14398 : _need_ad = _need_ad_u = true;
301 14398 : return _ad_u;
302 : }
303 :
304 9785 : const ADTemplateVariableGradient<OutputType> & adGradSln() const
305 : {
306 9785 : _need_ad = _need_ad_grad_u = true;
307 9785 : return _ad_grad_u;
308 : }
309 :
310 0 : const ADTemplateVariableGradient<OutputType> & adGradSlnDot() const
311 : {
312 0 : _need_ad = _need_ad_grad_u_dot = true;
313 :
314 0 : if (!_time_integrator)
315 : // If we don't have a time integrator (this will be the case for variables that are a part of
316 : // the AuxiliarySystem) then we have no way to calculate _ad_grad_u_dot and we are just going
317 : // to copy the values from _grad_u_dot. Of course in order to be able to do that we need to
318 : // calculate _grad_u_dot
319 0 : _need_grad_dot = true;
320 :
321 0 : return _ad_grad_u_dot;
322 : }
323 :
324 0 : const ADTemplateVariableSecond<OutputType> & adSecondSln() const
325 : {
326 0 : _need_ad = _need_ad_second_u = true;
327 0 : secondPhi();
328 0 : secondPhiFace();
329 0 : return _ad_second_u;
330 : }
331 :
332 : const ADTemplateVariableValue<OutputType> & adUDot() const;
333 :
334 : const ADTemplateVariableValue<OutputType> & adUDotDot() const;
335 :
336 : const FieldVariableValue & uDot() const;
337 :
338 : const FieldVariableValue & uDotDot() const;
339 :
340 : const FieldVariableValue & uDotOld() const;
341 :
342 : const FieldVariableValue & uDotDotOld() const;
343 :
344 15724 : const VariableValue & duDotDu() const
345 : {
346 15724 : _need_du_dot_du = true;
347 15724 : return _du_dot_du;
348 : }
349 :
350 156 : const VariableValue & duDotDotDu() const
351 : {
352 156 : _need_du_dotdot_du = true;
353 156 : return _du_dotdot_du;
354 : }
355 :
356 49 : const ADTemplateVariableCurl<OutputType> & adCurlSln() const
357 : {
358 49 : _need_ad = _need_ad_curl_u = true;
359 49 : curlPhi();
360 49 : curlPhiFace();
361 49 : return _ad_curl_u;
362 : }
363 :
364 : ///////////////////////// Nodal value getters ///////////////////////////////////////////
365 :
366 : const OutputType & nodalValueDot() const;
367 : const OutputType & nodalValueDotDot() const;
368 : const OutputType & nodalValueDotOld() const;
369 : const OutputType & nodalValueDotDotOld() const;
370 : const OutputType & nodalValueDuDotDu() const;
371 : const OutputType & nodalValueDuDotDotDu() const;
372 :
373 : const typename Moose::ADType<OutputType>::type & adNodalValue() const;
374 :
375 : /**
376 : * Set local DOF values and evaluate the values on quadrature points
377 : */
378 : void setDofValues(const DenseVector<DofValue> & values);
379 :
380 : ///@{
381 : /**
382 : * dof value setters
383 : */
384 : void setDofValue(const DofValue & value, unsigned int index);
385 : ///@}
386 :
387 : /**
388 : * Write a nodal value to the passed-in solution vector
389 : */
390 : void insertNodalValue(libMesh::NumericVector<libMesh::Number> & residual, const DofValue & v);
391 : DofValue getNodalValue(const Node & node, Moose::SolutionState state) const;
392 : DofValue
393 : getElementalValue(const Elem * elem, Moose::SolutionState state, unsigned int idx = 0) const;
394 :
395 : ///////////////////////////// dof indices ///////////////////////////////////////////////
396 :
397 : void getDofIndices(const Elem * elem, std::vector<dof_id_type> & dof_indices) const;
398 2606237727 : const std::vector<dof_id_type> & dofIndices() const { return _dof_indices; }
399 15456557 : unsigned int numberOfDofs() const { return _dof_indices.size(); }
400 694877108 : void clearDofIndices() { _dof_indices.clear(); }
401 :
402 : /**
403 : * Get the dof indices corresponding to the current element
404 : */
405 : void prepare();
406 :
407 : /**
408 : * Prepare degrees of freedom for the current node
409 : */
410 : void reinitNode();
411 :
412 : /**
413 : * Prepare dof indices and solution values for elemental auxiliary variables
414 : */
415 : void reinitAux();
416 :
417 : /**
418 : * Set _dof_indices to the degrees of freedom existing on the passed-in nodes
419 : * @param nodes The nodes to grab the degrees of freedom for
420 : */
421 : void reinitNodes(const std::vector<dof_id_type> & nodes);
422 :
423 : /**
424 : * Add passed in local DOF values to a solution vector
425 : */
426 : void addSolution(libMesh::NumericVector<libMesh::Number> & sol,
427 : const DenseVector<libMesh::Number> & v) const;
428 :
429 : /////////////////////////// DoF value getters /////////////////////////////////////
430 :
431 : const DofValues & dofValuesDot() const;
432 : const DofValues & dofValuesDotOld() const;
433 : const DofValues & dofValuesDotDot() const;
434 : const DofValues & dofValuesDotDotOld() const;
435 : const MooseArray<libMesh::Number> & dofValuesDuDotDu() const;
436 : const MooseArray<libMesh::Number> & dofValuesDuDotDotDu() const;
437 :
438 : /**
439 : * Return the AD dof values
440 : */
441 : const ADDofValues & adDofValues() const;
442 :
443 : /**
444 : * Return the AD time derivative values of degrees of freedom
445 : */
446 : const ADDofValues & adDofValuesDot() const;
447 :
448 : /////////////////////////////// Increment stuff ///////////////////////////////////////
449 :
450 : /**
451 : * Increment getter
452 : * @return The increment
453 : */
454 178 : const FieldVariableValue & increment() const { return _increment; }
455 :
456 : /**
457 : * Compute and store incremental change in solution at QPs based on increment_vec
458 : */
459 : void computeIncrementAtQps(const libMesh::NumericVector<libMesh::Number> & increment_vec);
460 :
461 : /**
462 : * Compute and store incremental change at the current node based on increment_vec
463 : */
464 : void computeIncrementAtNode(const libMesh::NumericVector<libMesh::Number> & increment_vec);
465 :
466 : private:
467 : /**
468 : * Helper method for assigning the _ad_dof_* arrays
469 : */
470 : void fetchADDofValues();
471 :
472 : /**
473 : * Helper method for assigning nodal values from their corresponding solution values (dof
474 : * values as they're referred to here in this class). This method is only truly meaningful
475 : * for nodal basis families
476 : */
477 : void assignADNodalValue();
478 :
479 : /**
480 : * Internal method for computeValues() and computeConstantMonomialValues()
481 : *
482 : * Monomial is a template parameter so that we get compile time optimization
483 : * for constant monomial vs non-constant monomial
484 : */
485 : template <bool constant_monomial>
486 : void computeValuesInternal();
487 :
488 : template <bool constant_monomial,
489 : typename DestinationType,
490 : typename ShapeType,
491 : typename DofValuesType>
492 : void fill(DestinationType & dest,
493 : const ShapeType & phi,
494 : const DofValuesType & dof_values,
495 : unsigned int nqp,
496 : std::size_t num_shapes);
497 :
498 : /// A const reference to the owning MooseVariableFE object
499 : const MooseVariableFE<OutputType> & _var;
500 :
501 : const libMesh::FEType & _fe_type;
502 :
503 : const unsigned int _var_num;
504 :
505 : const Assembly & _assembly;
506 :
507 : /// The element type this object is storing data for. This is either Element, Neighbor, or Lower
508 : Moose::ElementType _element_type;
509 :
510 : /// if variable is nodal
511 : bool _is_nodal;
512 :
513 : /// The dof index for the current node
514 : dof_id_type _nodal_dof_index;
515 :
516 : /// Continuity type of the variable
517 : libMesh::FEContinuity _continuity;
518 :
519 : /// Increment in the variable used in dampers
520 : FieldVariableValue _increment;
521 :
522 : /// AD nodal value
523 : typename Moose::ADType<OutputType>::type _ad_nodal_value;
524 :
525 : /// A zero AD variable
526 : ADReal _ad_zero;
527 :
528 : /// SolutionState second_u flags
529 : mutable bool _need_second = false;
530 : mutable bool _need_second_old = false;
531 : mutable bool _need_second_older = false;
532 : mutable bool _need_second_previous_nl = false;
533 :
534 : /// curl flags
535 : mutable bool _need_curl = false;
536 : mutable bool _need_curl_old = false;
537 : mutable bool _need_curl_older = false;
538 :
539 : /// divergence flags
540 : mutable bool _need_div = false;
541 : mutable bool _need_div_old = false;
542 : mutable bool _need_div_older = false;
543 :
544 : /// AD flags
545 : mutable bool _need_ad = false;
546 : mutable bool _need_ad_u = false;
547 : mutable bool _need_ad_grad_u = false;
548 : mutable bool _need_ad_grad_u_dot = false;
549 : mutable bool _need_ad_second_u = false;
550 : mutable bool _need_ad_curl_u = false;
551 : mutable bool _need_ad_div_u = false;
552 : mutable bool _need_ad_u_dot = false;
553 : mutable bool _need_ad_u_dotdot = false;
554 :
555 : bool _has_dof_indices;
556 :
557 : /// grad_u dots
558 : FieldVariableGradient _grad_u_dot;
559 : FieldVariableGradient _grad_u_dotdot;
560 :
561 : /// second_u
562 : FieldVariableSecond _second_u;
563 : FieldVariableSecond _second_u_old;
564 : FieldVariableSecond _second_u_older;
565 : FieldVariableSecond _second_u_previous_nl;
566 :
567 : /// curl_u
568 : FieldVariableCurl _curl_u;
569 : FieldVariableCurl _curl_u_old;
570 : FieldVariableCurl _curl_u_older;
571 :
572 : /// divergence_u
573 : FieldVariableDivergence _div_u;
574 : FieldVariableDivergence _div_u_old;
575 : FieldVariableDivergence _div_u_older;
576 :
577 : /// AD u
578 : ADTemplateVariableValue<OutputType> _ad_u;
579 : ADTemplateVariableGradient<OutputType> _ad_grad_u;
580 : ADTemplateVariableSecond<OutputType> _ad_second_u;
581 : ADDofValues _ad_dof_values;
582 : ADDofValues _ad_dofs_dot;
583 : MooseArray<ADRealEigenVector> _ad_dofs_dot_eigen;
584 : ADDofValues _ad_dofs_dotdot;
585 : MooseArray<ADRealEigenVector> _ad_dofs_dotdot_eigen;
586 : ADTemplateVariableValue<OutputType> _ad_u_dot;
587 : ADTemplateVariableValue<OutputType> _ad_u_dotdot;
588 : ADTemplateVariableGradient<OutputType> _ad_grad_u_dot;
589 : ADTemplateVariableCurl<OutputType> _ad_curl_u;
590 :
591 : // time derivatives
592 :
593 : /// u_dot (time derivative)
594 : FieldVariableValue _u_dot;
595 :
596 : /// u_dotdot (second time derivative)
597 : FieldVariableValue _u_dotdot;
598 :
599 : /// u_dot_old (time derivative)
600 : FieldVariableValue _u_dot_old;
601 :
602 : /// u_dotdot_old (second time derivative)
603 : FieldVariableValue _u_dotdot_old;
604 :
605 : /// derivative of u_dot wrt u
606 : VariableValue _du_dot_du;
607 :
608 : /// derivative of u_dotdot wrt u
609 : VariableValue _du_dotdot_du;
610 :
611 : /// The current qrule. This has to be a reference because the current qrule will be constantly
612 : /// changing. If we initialized this to point to one qrule, then in the next calculation we would
613 : /// be pointing to the wrong place!
614 : const QBase * const & _qrule;
615 : const QBase * const & _qrule_face;
616 :
617 : // Shape function values, gradients, second derivatives
618 : const FieldVariablePhiValue * _phi;
619 : const FieldVariablePhiGradient * _grad_phi;
620 : mutable const FieldVariablePhiSecond * _second_phi;
621 : mutable const FieldVariablePhiCurl * _curl_phi;
622 : mutable const FieldVariablePhiDivergence * _div_phi;
623 :
624 : // Mapped array phi
625 : MappedArrayVariablePhiGradient _mapped_grad_phi;
626 : MappedArrayVariablePhiGradient _mapped_grad_phi_face;
627 : MappedArrayVariablePhiGradient _mapped_grad_phi_neighbor;
628 : MappedArrayVariablePhiGradient _mapped_grad_phi_face_neighbor;
629 :
630 : // Values, gradients and second derivatives of shape function on faces
631 : const FieldVariablePhiValue * _phi_face;
632 : const FieldVariablePhiGradient * _grad_phi_face;
633 : mutable const FieldVariablePhiSecond * _second_phi_face;
634 : mutable const FieldVariablePhiCurl * _curl_phi_face;
635 : mutable const FieldVariablePhiDivergence * _div_phi_face;
636 :
637 : const ADTemplateVariablePhiGradient<OutputShape> * _ad_grad_phi;
638 : const ADTemplateVariablePhiGradient<OutputShape> * _ad_grad_phi_face;
639 :
640 : const QBase * _current_qrule;
641 : const FieldVariablePhiValue * _current_phi;
642 : const FieldVariablePhiGradient * _current_grad_phi;
643 : const FieldVariablePhiSecond * _current_second_phi;
644 : const FieldVariablePhiCurl * _current_curl_phi;
645 : const FieldVariablePhiDivergence * _current_div_phi;
646 : const ADTemplateVariablePhiGradient<OutputShape> * _current_ad_grad_phi;
647 :
648 : // dual mortar
649 : const bool _use_dual;
650 :
651 : std::function<const typename OutputTools<OutputType>::VariablePhiValue &(const Assembly &,
652 : libMesh::FEType)>
653 : _phi_assembly_method;
654 : std::function<const typename OutputTools<OutputShape>::VariablePhiValue &(const Assembly &,
655 : libMesh::FEType)>
656 : _phi_face_assembly_method;
657 :
658 : std::function<const typename OutputTools<OutputShape>::VariablePhiGradient &(const Assembly &,
659 : libMesh::FEType)>
660 : _grad_phi_assembly_method;
661 : std::function<const typename OutputTools<OutputShape>::VariablePhiGradient &(const Assembly &,
662 : libMesh::FEType)>
663 : _grad_phi_face_assembly_method;
664 :
665 : std::function<const typename OutputTools<OutputShape>::VariablePhiSecond &(const Assembly &,
666 : libMesh::FEType)>
667 : _second_phi_assembly_method;
668 : std::function<const typename OutputTools<OutputShape>::VariablePhiSecond &(const Assembly &,
669 : libMesh::FEType)>
670 : _second_phi_face_assembly_method;
671 :
672 : std::function<const typename OutputTools<OutputShape>::VariablePhiCurl &(const Assembly &,
673 : libMesh::FEType)>
674 : _curl_phi_assembly_method;
675 : std::function<const typename OutputTools<OutputShape>::VariablePhiCurl &(const Assembly &,
676 : libMesh::FEType)>
677 : _curl_phi_face_assembly_method;
678 :
679 : std::function<const typename OutputTools<OutputShape>::VariablePhiDivergence &(const Assembly &,
680 : libMesh::FEType)>
681 : _div_phi_assembly_method;
682 :
683 : std::function<const typename OutputTools<OutputShape>::VariablePhiDivergence &(const Assembly &,
684 : libMesh::FEType)>
685 : _div_phi_face_assembly_method;
686 :
687 : std::function<const ADTemplateVariablePhiGradient<OutputShape> &(const Assembly &,
688 : libMesh::FEType)>
689 : _ad_grad_phi_assembly_method;
690 : std::function<const ADTemplateVariablePhiGradient<OutputShape> &(const Assembly &,
691 : libMesh::FEType)>
692 : _ad_grad_phi_face_assembly_method;
693 :
694 : /// Pointer to time integrator
695 : const TimeIntegrator * _time_integrator;
696 :
697 : /// The current node. This has to be a reference because the current node will be constantly
698 : /// changing. If we initialized this to point to one node, then in the next calculation we would
699 : /// be pointing to the wrong place!
700 : const Node * const & _node;
701 :
702 : /// The current elem. This has to be a reference because the current elem will be constantly
703 : /// changing. If we initialized this to point to one elem, then in the next calculation we would
704 : /// be pointing to the wrong place!
705 : const Elem * const & _elem;
706 :
707 : /// Whether this variable is being calculated on a displaced system
708 : const bool _displaced;
709 :
710 : /// The current element side
711 : const unsigned int & _current_side;
712 :
713 : /// A dummy ADReal variable
714 522213 : ADReal _ad_real_dummy = 0;
715 :
716 : using MooseVariableDataBase<OutputType>::var;
717 : using MooseVariableDataBase<OutputType>::_sys;
718 : using MooseVariableDataBase<OutputType>::_subproblem;
719 : using MooseVariableDataBase<OutputType>::_need_vector_tag_dof_u;
720 : using MooseVariableDataBase<OutputType>::_need_matrix_tag_dof_u;
721 : using MooseVariableDataBase<OutputType>::_vector_tags_dof_u;
722 : using MooseVariableDataBase<OutputType>::_matrix_tags_dof_u;
723 : using MooseVariableDataBase<OutputType>::_vector_tag_u;
724 : using MooseVariableDataBase<OutputType>::_need_vector_tag_u;
725 : using MooseVariableDataBase<OutputType>::_vector_tag_grad;
726 : using MooseVariableDataBase<OutputType>::_need_vector_tag_grad;
727 : using MooseVariableDataBase<OutputType>::_matrix_tag_u;
728 : using MooseVariableDataBase<OutputType>::_need_matrix_tag_u;
729 : using MooseVariableDataBase<OutputType>::_dof_indices;
730 : using MooseVariableDataBase<OutputType>::_has_dof_values;
731 : using MooseVariableDataBase<OutputType>::fetchDofValues;
732 : using MooseVariableDataBase<OutputType>::assignNodalValue;
733 : using MooseVariableDataBase<OutputType>::zeroSizeDofValues;
734 : using MooseVariableDataBase<OutputType>::_solution_tag;
735 : using MooseVariableDataBase<OutputType>::_old_solution_tag;
736 : using MooseVariableDataBase<OutputType>::_older_solution_tag;
737 : using MooseVariableDataBase<OutputType>::_previous_nl_solution_tag;
738 : using MooseVariableDataBase<OutputType>::_dof_map;
739 : using MooseVariableDataBase<OutputType>::_need_u_dot;
740 : using MooseVariableDataBase<OutputType>::_need_u_dotdot;
741 : using MooseVariableDataBase<OutputType>::_need_u_dot_old;
742 : using MooseVariableDataBase<OutputType>::_need_u_dotdot_old;
743 : using MooseVariableDataBase<OutputType>::_need_du_dot_du;
744 : using MooseVariableDataBase<OutputType>::_need_du_dotdot_du;
745 : using MooseVariableDataBase<OutputType>::_need_grad_dot;
746 : using MooseVariableDataBase<OutputType>::_need_grad_dotdot;
747 : using MooseVariableDataBase<OutputType>::_need_dof_values_dot;
748 : using MooseVariableDataBase<OutputType>::_need_dof_values_dotdot;
749 : using MooseVariableDataBase<OutputType>::_need_dof_values_dot_old;
750 : using MooseVariableDataBase<OutputType>::_need_dof_values_dotdot_old;
751 : using MooseVariableDataBase<OutputType>::_need_dof_du_dot_du;
752 : using MooseVariableDataBase<OutputType>::_need_dof_du_dotdot_du;
753 : using MooseVariableDataBase<OutputType>::_dof_values_dot;
754 : using MooseVariableDataBase<OutputType>::_dof_values_dotdot;
755 : using MooseVariableDataBase<OutputType>::_dof_values_dot_old;
756 : using MooseVariableDataBase<OutputType>::_dof_values_dotdot_old;
757 : using MooseVariableDataBase<OutputType>::_dof_du_dot_du;
758 : using MooseVariableDataBase<OutputType>::_dof_du_dotdot_du;
759 : using MooseVariableDataBase<OutputType>::_tid;
760 : using MooseVariableDataBase<OutputType>::_nodal_value_dot;
761 : using MooseVariableDataBase<OutputType>::_nodal_value_dotdot;
762 : using MooseVariableDataBase<OutputType>::_nodal_value_dot_old;
763 : using MooseVariableDataBase<OutputType>::_nodal_value_dotdot_old;
764 : using MooseVariableDataBase<OutputType>::_required_vector_tags;
765 : using MooseVariableDataBase<OutputType>::_count;
766 : };
767 :
768 : /////////////////////// General template definitions //////////////////////////////////////
769 :
770 : template <typename OutputType>
771 : const typename MooseVariableData<OutputType>::ADDofValues &
772 390 : MooseVariableData<OutputType>::adDofValues() const
773 : {
774 390 : _need_ad = true;
775 390 : return _ad_dof_values;
776 : }
777 :
778 : template <typename OutputType>
779 : const typename MooseVariableData<OutputType>::ADDofValues &
780 13 : MooseVariableData<OutputType>::adDofValuesDot() const
781 : {
782 13 : _need_ad = _need_ad_u_dot = true;
783 13 : if (!_time_integrator)
784 : // See explanation in adUDot() body
785 0 : _need_u_dot = true;
786 13 : return _ad_dofs_dot;
787 : }
788 :
789 : template <typename OutputType>
790 : const typename Moose::ADType<OutputType>::type &
791 1935 : MooseVariableData<OutputType>::adNodalValue() const
792 : {
793 1935 : _need_ad = true;
794 1935 : return _ad_nodal_value;
795 : }
796 :
797 : template <typename OutputType>
798 : const ADTemplateVariableValue<OutputType> &
799 639 : MooseVariableData<OutputType>::adUDot() const
800 : {
801 639 : _need_ad = _need_ad_u_dot = true;
802 :
803 639 : if (!_time_integrator)
804 : // If we don't have a time integrator (this will be the case for variables that are a part of
805 : // the AuxiliarySystem) then we have no way to calculate _ad_u_dot and we are just going to
806 : // copy the values from _u_dot. Of course in order to be able to do that we need to calculate
807 : // _u_dot
808 0 : _need_u_dot = true;
809 :
810 639 : return _ad_u_dot;
811 : }
812 :
813 : template <typename OutputType>
814 : const ADTemplateVariableValue<OutputType> &
815 52 : MooseVariableData<OutputType>::adUDotDot() const
816 : {
817 52 : _need_ad = _need_ad_u_dotdot = true;
818 :
819 52 : if (!_time_integrator)
820 : // If we don't have a time integrator (this will be the case for variables that are a part
821 : // of the AuxiliarySystem) then we have no way to calculate _ad_u_dotdot and we are just
822 : // going to copy the values from _u_dotdot. Of course in order to be able to do that we need
823 : // to calculate _u_dotdot
824 0 : _need_u_dotdot = true;
825 :
826 52 : return _ad_u_dotdot;
827 : }
828 :
829 : template <typename OutputType>
830 : const ADTemplateVariablePhiGradient<typename MooseVariableData<OutputType>::OutputShape> &
831 4739 : MooseVariableData<OutputType>::adGradPhi() const
832 : {
833 4739 : if (_element_type == Moose::ElementType::Neighbor || _element_type == Moose::ElementType::Lower)
834 0 : mooseError("Unsupported element type: ", Moose::stringify(_element_type));
835 : mooseAssert(_ad_grad_phi, "this should be non-null");
836 4739 : return *_ad_grad_phi;
837 : }
838 :
839 : template <typename OutputType>
840 : const ADTemplateVariablePhiGradient<typename MooseVariableData<OutputType>::OutputShape> &
841 1704 : MooseVariableData<OutputType>::adGradPhiFace() const
842 : {
843 1704 : if (_element_type == Moose::ElementType::Neighbor || _element_type == Moose::ElementType::Lower)
844 0 : mooseError("Unsupported element type: ", Moose::stringify(_element_type));
845 : mooseAssert(_ad_grad_phi_face, "this should be non-null");
846 1704 : return *_ad_grad_phi_face;
847 : }
|