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