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 : #include "MooseVariableDataFV.h"
11 : #include "MooseVariableField.h"
12 : #include "Assembly.h"
13 : #include "MooseError.h"
14 : #include "DisplacedSystem.h"
15 : #include "TimeIntegrator.h"
16 : #include "MooseVariableFV.h"
17 : #include "MooseTypes.h"
18 : #include "MooseMesh.h"
19 : #include "Attributes.h"
20 : #include "FVDirichletBCBase.h"
21 : #include "SubProblem.h"
22 : #include "FVKernel.h"
23 : #include "ADUtils.h"
24 :
25 : #include "libmesh/quadrature.h"
26 : #include "libmesh/fe_base.h"
27 : #include "libmesh/system.h"
28 : #include "libmesh/type_n_tensor.h"
29 :
30 : template <typename OutputType>
31 12748 : MooseVariableDataFV<OutputType>::MooseVariableDataFV(const MooseVariableFV<OutputType> & var,
32 : SystemBase & sys,
33 : THREAD_ID tid,
34 : Moose::ElementType element_type,
35 : const Elem * const & elem)
36 :
37 : : MooseVariableDataBase<OutputType>(var, sys, tid),
38 : MeshChangedInterface(var.parameters()),
39 12748 : _var(var),
40 25496 : _fe_type(_var.feType()),
41 12748 : _var_num(_var.number()),
42 12748 : _assembly(_subproblem.assembly(_tid, var.kind() == Moose::VAR_SOLVER ? sys.number() : 0)),
43 12748 : _element_type(element_type),
44 12748 : _ad_zero(0),
45 12748 : _need_second(false),
46 12748 : _need_second_old(false),
47 12748 : _need_second_older(false),
48 12748 : _need_second_previous_nl(false),
49 12748 : _need_curl(false),
50 12748 : _need_curl_old(false),
51 12748 : _need_curl_older(false),
52 : // for FV variables, they use each other's ad_u values to compute ghost cell
53 : // values - we don't have any way to propagate these inter-variable-data
54 : // dependencies. So if something needs an ad_u value, that need isn't
55 : // propagated through to both the element and the neighbor
56 : // data structures. So instead just set _need_ad+_need_ad_u values to true always.
57 12748 : _need_ad(true),
58 12748 : _need_ad_u(true),
59 12748 : _need_ad_u_dot(false),
60 12748 : _need_ad_u_dotdot(false),
61 12748 : _need_ad_grad_u(false),
62 12748 : _need_ad_grad_u_dot(false),
63 12748 : _need_ad_second_u(false),
64 12748 : _time_integrator(_sys.queryTimeIntegrator(_var_num)),
65 12748 : _elem(elem),
66 12748 : _displaced(dynamic_cast<const DisplacedSystem *>(&_sys) ? true : false),
67 63740 : _qrule(nullptr)
68 : {
69 12748 : _fv_elemental_kernel_query_cache =
70 12748 : _subproblem.getMooseApp().theWarehouse().query().template condition<AttribSystem>(
71 : "FVElementalKernel");
72 12748 : _fv_flux_kernel_query_cache =
73 12748 : _subproblem.getMooseApp().theWarehouse().query().template condition<AttribSystem>(
74 : "FVFluxKernel");
75 12748 : }
76 :
77 : template <typename OutputType>
78 : void
79 81476192 : MooseVariableDataFV<OutputType>::setGeometry(Moose::GeometryType gm_type)
80 : {
81 81476192 : switch (gm_type)
82 : {
83 26958530 : case Moose::Volume:
84 : {
85 26958530 : _qrule = _assembly.qRule();
86 : // TODO: set integration multiplier to cell volume
87 26958530 : break;
88 : }
89 54517662 : case Moose::Face:
90 : {
91 54517662 : _qrule = _assembly.qRuleFace();
92 : // TODO: set integration multiplier to face area
93 54517662 : break;
94 : }
95 : }
96 81476192 : }
97 :
98 : template <typename OutputType>
99 : const typename MooseVariableDataFV<OutputType>::FieldVariableValue &
100 0 : MooseVariableDataFV<OutputType>::uDot() const
101 : {
102 0 : if (_sys.solutionUDot())
103 : {
104 0 : _var.requireQpComputations();
105 0 : _need_u_dot = true;
106 0 : return _u_dot;
107 : }
108 : else
109 0 : mooseError("MooseVariableFE: Time derivative of solution (`u_dot`) is not stored. Please set "
110 : "uDotRequested() to true in FEProblemBase before requesting `u_dot`.");
111 : }
112 :
113 : template <typename OutputType>
114 : const typename MooseVariableDataFV<OutputType>::FieldVariableValue &
115 0 : MooseVariableDataFV<OutputType>::uDotDot() const
116 : {
117 0 : if (_sys.solutionUDotDot())
118 : {
119 0 : _var.requireQpComputations();
120 0 : _need_u_dotdot = true;
121 0 : return _u_dotdot;
122 : }
123 : else
124 0 : mooseError("MooseVariableFE: Second time derivative of solution (`u_dotdot`) is not stored. "
125 : "Please set uDotDotRequested() to true in FEProblemBase before requesting "
126 : "`u_dotdot`.");
127 : }
128 :
129 : template <typename OutputType>
130 : const typename MooseVariableDataFV<OutputType>::FieldVariableValue &
131 0 : MooseVariableDataFV<OutputType>::uDotOld() const
132 : {
133 0 : if (_sys.solutionUDotOld())
134 : {
135 0 : _var.requireQpComputations();
136 0 : _need_u_dot_old = true;
137 0 : return _u_dot_old;
138 : }
139 : else
140 0 : mooseError("MooseVariableFE: Old time derivative of solution (`u_dot_old`) is not stored. "
141 : "Please set uDotOldRequested() to true in FEProblemBase before requesting "
142 : "`u_dot_old`.");
143 : }
144 :
145 : template <typename OutputType>
146 : const typename MooseVariableDataFV<OutputType>::FieldVariableValue &
147 0 : MooseVariableDataFV<OutputType>::uDotDotOld() const
148 : {
149 0 : if (_sys.solutionUDotDotOld())
150 : {
151 0 : _var.requireQpComputations();
152 0 : _need_u_dotdot_old = true;
153 0 : return _u_dotdot_old;
154 : }
155 : else
156 0 : mooseError("MooseVariableFE: Old second time derivative of solution (`u_dotdot_old`) is not "
157 : "stored. Please set uDotDotOldRequested() to true in FEProblemBase before "
158 : "requesting `u_dotdot_old`");
159 : }
160 :
161 : template <typename OutputType>
162 : const typename MooseVariableDataFV<OutputType>::FieldVariableValue &
163 7198 : MooseVariableDataFV<OutputType>::sln(Moose::SolutionState state) const
164 : {
165 7198 : _var.requireQpComputations();
166 7198 : return MooseVariableDataBase<OutputType>::sln(state);
167 : }
168 :
169 : template <typename OutputType>
170 : const typename MooseVariableDataFV<OutputType>::FieldVariableGradient &
171 0 : MooseVariableDataFV<OutputType>::gradSlnDot() const
172 : {
173 0 : if (_sys.solutionUDot())
174 : {
175 0 : _var.requireQpComputations();
176 0 : _need_grad_dot = true;
177 0 : return _grad_u_dot;
178 : }
179 : else
180 0 : mooseError("MooseVariableFE: Time derivative of solution (`u_dot`) is not stored. Please set "
181 : "uDotRequested() to true in FEProblemBase before requesting `u_dot`.");
182 : }
183 :
184 : template <typename OutputType>
185 : const typename MooseVariableDataFV<OutputType>::FieldVariableGradient &
186 0 : MooseVariableDataFV<OutputType>::gradSlnDotDot() const
187 : {
188 0 : if (_sys.solutionUDotDot())
189 : {
190 0 : _var.requireQpComputations();
191 0 : _need_grad_dotdot = true;
192 0 : return _grad_u_dotdot;
193 : }
194 : else
195 0 : mooseError("MooseVariableFE: Second time derivative of solution (`u_dotdot`) is not stored. "
196 : "Please set uDotDotRequested() to true in FEProblemBase before requesting "
197 : "`u_dotdot`.");
198 : }
199 :
200 : template <typename OutputType>
201 : const typename MooseVariableDataFV<OutputType>::FieldVariableSecond &
202 0 : MooseVariableDataFV<OutputType>::secondSln(Moose::SolutionState state) const
203 : {
204 0 : _var.requireQpComputations();
205 0 : switch (state)
206 : {
207 0 : case Moose::Current:
208 : {
209 0 : _need_second = true;
210 0 : return _second_u;
211 : }
212 :
213 0 : case Moose::Old:
214 : {
215 0 : _need_second_old = true;
216 0 : return _second_u_old;
217 : }
218 :
219 0 : case Moose::Older:
220 : {
221 0 : _need_second_older = true;
222 0 : return _second_u_older;
223 : }
224 :
225 0 : case Moose::PreviousNL:
226 : {
227 0 : _need_second_previous_nl = true;
228 0 : return _second_u_previous_nl;
229 : }
230 :
231 0 : default:
232 : // We should never get here but gcc requires that we have a default. See
233 : // htpps://stackoverflow.com/questions/18680378/after-defining-case-for-all-enum-values-compiler-still-says-control-reaches-e
234 0 : mooseError("Unknown SolutionState!");
235 : }
236 : }
237 :
238 : template <typename OutputType>
239 : const typename MooseVariableDataFV<OutputType>::FieldVariableCurl &
240 0 : MooseVariableDataFV<OutputType>::curlSln(Moose::SolutionState state) const
241 : {
242 0 : _var.requireQpComputations();
243 0 : switch (state)
244 : {
245 0 : case Moose::Current:
246 : {
247 0 : _need_curl = true;
248 0 : return _curl_u;
249 : }
250 :
251 0 : case Moose::Old:
252 : {
253 0 : _need_curl_old = true;
254 0 : return _curl_u_old;
255 : }
256 :
257 0 : case Moose::Older:
258 : {
259 0 : _need_curl_older = true;
260 0 : return _curl_u_older;
261 : }
262 :
263 0 : default:
264 0 : mooseError("We don't currently support curl from the previous non-linear iteration");
265 : }
266 : }
267 :
268 : template <typename OutputType>
269 : void
270 79363710 : MooseVariableDataFV<OutputType>::initializeSolnVars()
271 : {
272 : auto && active_coupleable_matrix_tags =
273 79363710 : _sys.subproblem().getActiveFEVariableCoupleableMatrixTags(_tid);
274 : mooseAssert(_qrule, "We should have a non-null qrule");
275 79363710 : const auto nqp = _qrule->n_points();
276 :
277 158754420 : for (auto tag : _required_vector_tags)
278 79390710 : if (_need_vector_tag_u[tag])
279 : {
280 79366710 : _vector_tag_u[tag].resize(nqp);
281 79366710 : assignForAllQps(0, _vector_tag_u[tag], nqp);
282 : }
283 :
284 79363774 : for (auto tag : active_coupleable_matrix_tags)
285 64 : if (_need_matrix_tag_u[tag])
286 : {
287 32 : _matrix_tag_u[tag].resize(nqp);
288 32 : assignForAllQps(0, _matrix_tag_u[tag], nqp);
289 : }
290 :
291 79363710 : if (_need_second)
292 : {
293 0 : _second_u.resize(nqp);
294 0 : assignForAllQps(0, _second_u, nqp);
295 : }
296 :
297 79363710 : if (_need_curl)
298 : {
299 0 : _curl_u.resize(nqp);
300 0 : assignForAllQps(0, _curl_u, nqp);
301 : }
302 :
303 79363710 : if (_need_second_previous_nl)
304 : {
305 0 : _second_u_previous_nl.resize(nqp);
306 0 : assignForAllQps(0, _second_u_previous_nl, nqp);
307 : }
308 :
309 79363710 : if (_subproblem.isTransient())
310 : {
311 47777638 : if (_need_u_dot)
312 : {
313 0 : _u_dot.resize(nqp);
314 0 : assignForAllQps(0, _u_dot, nqp);
315 : }
316 :
317 47777638 : if (_need_u_dotdot)
318 : {
319 0 : _u_dotdot.resize(nqp);
320 0 : assignForAllQps(0, _u_dotdot, nqp);
321 : }
322 :
323 47777638 : if (_need_u_dot_old)
324 : {
325 0 : _u_dot_old.resize(nqp);
326 0 : assignForAllQps(0, _u_dot_old, nqp);
327 : }
328 :
329 47777638 : if (_need_u_dotdot_old)
330 : {
331 0 : _u_dotdot_old.resize(nqp);
332 0 : assignForAllQps(0, _u_dotdot_old, nqp);
333 : }
334 :
335 47777638 : if (_need_du_dot_du)
336 : {
337 0 : _du_dot_du.resize(nqp);
338 0 : assignForAllQps(0, _du_dot_du, nqp);
339 : }
340 :
341 47777638 : if (_need_du_dotdot_du)
342 : {
343 0 : _du_dotdot_du.resize(nqp);
344 0 : assignForAllQps(0, _du_dotdot_du, nqp);
345 : }
346 :
347 47777638 : if (_need_grad_dot)
348 : {
349 0 : _grad_u_dot.resize(nqp);
350 0 : assignForAllQps(0, _grad_u_dot, nqp);
351 : }
352 :
353 47777638 : if (_need_grad_dotdot)
354 : {
355 0 : _grad_u_dotdot.resize(nqp);
356 0 : assignForAllQps(0, _grad_u_dotdot, nqp);
357 : }
358 :
359 47777638 : if (_need_second_old)
360 : {
361 0 : _second_u_old.resize(nqp);
362 0 : assignForAllQps(0, _second_u_old, nqp);
363 : }
364 :
365 47777638 : if (_need_curl_old)
366 : {
367 0 : _curl_u_old.resize(nqp);
368 0 : assignForAllQps(0, _curl_u_old, nqp);
369 : }
370 :
371 47777638 : if (_need_second_older)
372 : {
373 0 : _second_u_older.resize(nqp);
374 0 : assignForAllQps(0, _second_u_older, nqp);
375 : }
376 : }
377 79363710 : }
378 :
379 : template <typename OutputType>
380 : void
381 52060924 : MooseVariableDataFV<OutputType>::computeValuesFace(const FaceInfo & /*fi*/)
382 : {
383 52060924 : _has_dirichlet_bc = false;
384 52060924 : _dof_map.dof_indices(_elem, _dof_indices, _var_num);
385 :
386 : // TODO: compute reconstructed values somehow. For now, just do the trivial
387 : // reconstruction where we take the const cell value from the centroid and
388 : // use that value on the face. How will users affect the reconstruction
389 : // method used here? After reconstruction, we should cache the computed
390 : // soln/gradients per element so we don't have to recompute them again for
391 : // other faces that share an element with this face.
392 : //
393 : // TODO: Also we need to be able to track (AD) derivatives through the
394 : // reconstruction process - how will we do that?
395 :
396 52060924 : computeValues();
397 :
398 : // TODO: maybe initialize a separate _grad_u_interface here that is
399 : // only used for diffusion terms that need an interface gradient. Also -
400 : // it would need cross-diffusion corrections for non-orthogonal meshes
401 : // eventually. Or maybe we just leave this alone zero and have users
402 : // calculate whatever grad_interface value they want and provide some
403 : // helper functions for cross-diffusion correction.
404 :
405 : // TODO: figure out how to store old/older values of reconstructed
406 : // solutions/gradient states without having to re-reconstruct them from
407 : // older dof/soln values.
408 52060924 : }
409 :
410 : template <typename OutputType>
411 : const std::vector<dof_id_type> &
412 1304019828 : MooseVariableDataFV<OutputType>::initDofIndices()
413 : {
414 1304019828 : Moose::initDofIndices(*this, *_elem);
415 1304019828 : return _dof_indices;
416 : }
417 :
418 : template <typename OutputType>
419 : void
420 79363710 : MooseVariableDataFV<OutputType>::computeValues()
421 : {
422 79363710 : initDofIndices();
423 79363710 : initializeSolnVars();
424 :
425 79363710 : unsigned int num_dofs = _dof_indices.size();
426 :
427 79363710 : if (num_dofs > 0)
428 76958095 : fetchDoFValues();
429 : else
430 : // We don't have any dofs. There's nothing to do
431 2405615 : return;
432 :
433 : mooseAssert(num_dofs == 1 && _vector_tags_dof_u[_solution_tag].size() == 1,
434 : "There should only be one dof per elem for FV variables");
435 :
436 76958095 : bool is_transient = _subproblem.isTransient();
437 76958095 : const auto nqp = _qrule->n_points();
438 : auto && active_coupleable_matrix_tags =
439 76958095 : _sys.subproblem().getActiveFEVariableCoupleableMatrixTags(_tid);
440 :
441 76958095 : bool second_required =
442 76958095 : _need_second || _need_second_old || _need_second_older || _need_second_previous_nl;
443 76958095 : bool curl_required = _need_curl || _need_curl_old;
444 :
445 155480781 : for (const auto qp : make_range(nqp))
446 : {
447 78522686 : if (is_transient)
448 : {
449 47765891 : if (_need_u_dot)
450 0 : _u_dot[qp] = _dof_values_dot[0];
451 :
452 47765891 : if (_need_u_dotdot)
453 0 : _u_dotdot[qp] = _dof_values_dotdot[0];
454 :
455 47765891 : if (_need_u_dot_old)
456 0 : _u_dot_old[qp] = _dof_values_dot_old[0];
457 :
458 47765891 : if (_need_u_dotdot_old)
459 0 : _u_dotdot_old[qp] = _dof_values_dotdot_old[0];
460 :
461 47765891 : if (_need_du_dot_du)
462 0 : _du_dot_du[qp] = _dof_du_dot_du[0];
463 :
464 47765891 : if (_need_du_dotdot_du)
465 0 : _du_dotdot_du[qp] = _dof_du_dotdot_du[0];
466 : }
467 :
468 78522686 : if (second_required)
469 : {
470 0 : if (_need_second)
471 0 : _second_u[qp] = 0.;
472 :
473 0 : if (_need_second_previous_nl)
474 0 : _second_u_previous_nl[qp] = 0.;
475 :
476 0 : if (is_transient)
477 : {
478 0 : if (_need_second_old)
479 0 : _second_u_old[qp] = 0.;
480 :
481 0 : if (_need_second_older)
482 0 : _second_u_older[qp] = 0.;
483 : }
484 : }
485 :
486 78522686 : if (curl_required)
487 : {
488 0 : if (_need_curl)
489 0 : _curl_u[qp] = 0.;
490 :
491 0 : if (is_transient && _need_curl_old)
492 0 : _curl_u_old[qp] = 0.;
493 : }
494 :
495 157072372 : for (auto tag : _required_vector_tags)
496 78549686 : if (_need_vector_tag_u[tag])
497 78525686 : _vector_tag_u[tag][qp] = _vector_tags_dof_u[tag][0];
498 :
499 78522750 : for (auto tag : active_coupleable_matrix_tags)
500 64 : if (_need_matrix_tag_u[tag])
501 32 : _matrix_tag_u[tag][qp] = _matrix_tags_dof_u[tag][0];
502 : }
503 :
504 : // Automatic differentiation
505 76958095 : if (_need_ad)
506 76958095 : computeAD(num_dofs, nqp);
507 : }
508 :
509 : template <typename OutputType>
510 : void
511 76958095 : MooseVariableDataFV<OutputType>::computeAD(const unsigned int num_dofs, const unsigned int nqp)
512 : {
513 : // This query and if check prevent running this code when we have FV
514 : // variables, but no kernels. When this happens, maxVarNDofsPerElem is
515 : // not computed (because no kernels) and is zero giving nonsense offsets for
516 : // AD stuff. So we just skip all this when that is the case. Maybe there
517 : // is a better way to do this - like just checking if getMaxVarNDofsPerElem
518 : // returns zero?
519 76958095 : std::vector<FVKernel *> ks;
520 76958095 : _fv_elemental_kernel_query_cache.queryInto(ks);
521 76958095 : if (ks.size() == 0)
522 : {
523 11418297 : _fv_flux_kernel_query_cache.queryInto(ks);
524 11418297 : if (ks.size() == 0)
525 640345 : return;
526 : }
527 :
528 76317750 : _ad_dof_values.resize(num_dofs);
529 76317750 : if (_need_ad_u)
530 76317750 : _ad_u.resize(nqp);
531 :
532 76317750 : if (_need_ad_grad_u)
533 9554144 : _ad_grad_u.resize(nqp);
534 :
535 76317750 : if (_need_ad_second_u)
536 0 : _ad_second_u.resize(nqp);
537 :
538 76317750 : if (_need_ad_u_dot)
539 : {
540 14437453 : _ad_dofs_dot.resize(num_dofs);
541 14437453 : _ad_u_dot.resize(nqp);
542 : }
543 :
544 76317750 : if (_need_ad_u_dotdot)
545 : {
546 580900 : _ad_dofs_dotdot.resize(num_dofs);
547 580900 : _ad_u_dotdot.resize(nqp);
548 : }
549 :
550 76317750 : if (_need_ad_second_u)
551 0 : assignForAllQps(0, _ad_second_u, nqp);
552 :
553 76317750 : if (_need_ad_u_dot)
554 14437453 : assignForAllQps(_ad_zero, _ad_u_dot, nqp);
555 :
556 76317750 : if (_need_ad_u_dotdot)
557 580900 : assignForAllQps(_ad_zero, _ad_u_dotdot, nqp);
558 :
559 76317750 : const bool do_derivatives =
560 76317750 : ADReal::do_derivatives && _sys.number() == _subproblem.currentNlSysNum();
561 :
562 152635500 : for (unsigned int i = 0; i < num_dofs; i++)
563 : {
564 76317750 : _ad_dof_values[i] = (*_sys.currentSolution())(_dof_indices[i]);
565 :
566 : // NOTE! You have to do this AFTER setting the value!
567 76317750 : if (do_derivatives)
568 18322822 : Moose::derivInsert(_ad_dof_values[i].derivatives(), _dof_indices[i], 1.);
569 :
570 76317750 : if (_need_ad_u_dot && safeToComputeADUDot() && _time_integrator->dt())
571 : {
572 14392543 : _ad_dofs_dot[i] = _ad_dof_values[i];
573 28205786 : _time_integrator->computeADTimeDerivatives(_ad_dofs_dot[i],
574 14392543 : _dof_indices[i],
575 14392543 : _need_ad_u_dotdot ? _ad_dofs_dotdot[i]
576 : : _ad_real_dummy);
577 : }
578 : }
579 :
580 76317750 : if (_need_ad_u)
581 76317750 : assignForAllQps(_ad_dof_values[0], _ad_u, nqp);
582 :
583 76317750 : if (_need_ad_grad_u)
584 9554144 : assignForAllQps(static_cast<const MooseVariableFV<OutputType> &>(_var).adGradSln(
585 19108288 : _elem, Moose::currentState()),
586 9554144 : _ad_grad_u,
587 : nqp);
588 :
589 76317750 : if (_need_ad_u_dot)
590 : {
591 14437453 : if (safeToComputeADUDot())
592 : {
593 14437453 : assignForAllQps(_ad_dofs_dot[0], _ad_u_dot, nqp);
594 14437453 : if (_need_ad_u_dotdot)
595 580900 : assignForAllQps(_ad_dofs_dotdot[0], _ad_u_dotdot, nqp);
596 : }
597 : else
598 : {
599 0 : assignForAllQps(_u_dot[0], _ad_u_dot, nqp);
600 0 : if (_need_ad_u_dotdot)
601 0 : assignForAllQps(_u_dotdot[0], _ad_u_dotdot, nqp);
602 : }
603 : }
604 76958095 : }
605 :
606 : template <typename OutputType>
607 : void
608 390600 : MooseVariableDataFV<OutputType>::setDofValue(const OutputData & value, unsigned int index)
609 : {
610 : mooseAssert(index == 0, "We only ever have one dof value locally");
611 390600 : _vector_tags_dof_u[_solution_tag][index] = value;
612 390600 : _has_dof_values = true;
613 :
614 390600 : auto & u = _vector_tag_u[_solution_tag];
615 : // Update the qp values as well
616 1075792 : for (const auto qp : index_range(u))
617 685192 : u[qp] = value;
618 :
619 390600 : if (_need_ad_u)
620 627840 : for (const auto qp : index_range(_ad_u))
621 237240 : _ad_u[qp] = value;
622 390600 : }
623 :
624 : template <typename OutputType>
625 : void
626 4000 : MooseVariableDataFV<OutputType>::setDofValues(const DenseVector<OutputData> & values)
627 : {
628 4000 : auto & dof_values = _vector_tags_dof_u[_solution_tag];
629 8000 : for (unsigned int i = 0; i < values.size(); i++)
630 4000 : dof_values[i] = values(i);
631 4000 : _has_dof_values = true;
632 4000 : }
633 :
634 : template <typename OutputType>
635 : typename MooseVariableDataFV<OutputType>::OutputData
636 2326 : MooseVariableDataFV<OutputType>::getElementalValue(const Elem * elem,
637 : Moose::SolutionState state,
638 : unsigned int idx) const
639 : {
640 2326 : std::vector<dof_id_type> dof_indices;
641 2326 : _dof_map.dof_indices(elem, dof_indices, _var_num);
642 :
643 2326 : switch (state)
644 : {
645 2326 : case Moose::Current:
646 2326 : return (*_sys.currentSolution())(dof_indices[idx]);
647 :
648 0 : case Moose::Old:
649 0 : return _sys.solutionOld()(dof_indices[idx]);
650 :
651 0 : case Moose::Older:
652 0 : return _sys.solutionOlder()(dof_indices[idx]);
653 :
654 0 : default:
655 0 : mooseError("PreviousNL not currently supported for getElementalValue");
656 : }
657 2326 : }
658 :
659 : template <typename OutputType>
660 : void
661 0 : MooseVariableDataFV<OutputType>::getDofIndices(const Elem * elem,
662 : std::vector<dof_id_type> & dof_indices) const
663 : {
664 0 : _dof_map.dof_indices(elem, dof_indices, _var_num);
665 0 : }
666 :
667 : template <typename OutputType>
668 : const typename MooseVariableDataFV<OutputType>::DoFValue &
669 0 : MooseVariableDataFV<OutputType>::dofValuesDot() const
670 : {
671 0 : if (_sys.solutionUDot())
672 : {
673 0 : _need_dof_values_dot = true;
674 0 : return _dof_values_dot;
675 : }
676 : else
677 0 : mooseError(
678 : "MooseVariableDataFV: Time derivative of solution (`u_dot`) is not stored. Please set "
679 : "uDotRequested() to true in FEProblemBase before requesting `u_dot`.");
680 : }
681 :
682 : template <typename OutputType>
683 : const typename MooseVariableDataFV<OutputType>::DoFValue &
684 0 : MooseVariableDataFV<OutputType>::dofValuesDotDot() const
685 : {
686 0 : if (_sys.solutionUDotDot())
687 : {
688 0 : _need_dof_values_dotdot = true;
689 0 : return _dof_values_dotdot;
690 : }
691 : else
692 0 : mooseError(
693 : "MooseVariableDataFV: Second time derivative of solution (`u_dotdot`) is not stored. "
694 : "Please set uDotDotRequested() to true in FEProblemBase before requesting "
695 : "`u_dotdot`.");
696 : }
697 :
698 : template <typename OutputType>
699 : const typename MooseVariableDataFV<OutputType>::DoFValue &
700 0 : MooseVariableDataFV<OutputType>::dofValuesDotOld() const
701 : {
702 0 : if (_sys.solutionUDotOld())
703 : {
704 0 : _need_dof_values_dot_old = true;
705 0 : return _dof_values_dot_old;
706 : }
707 : else
708 0 : mooseError("MooseVariableDataFV: Old time derivative of solution (`u_dot_old`) is not stored. "
709 : "Please set uDotOldRequested() to true in FEProblemBase before requesting "
710 : "`u_dot_old`.");
711 : }
712 :
713 : template <typename OutputType>
714 : const typename MooseVariableDataFV<OutputType>::DoFValue &
715 0 : MooseVariableDataFV<OutputType>::dofValuesDotDotOld() const
716 : {
717 0 : if (_sys.solutionUDotDotOld())
718 : {
719 0 : _need_dof_values_dotdot_old = true;
720 0 : return _dof_values_dotdot_old;
721 : }
722 : else
723 0 : mooseError(
724 : "MooseVariableDataFV: Old second time derivative of solution (`u_dotdot_old`) is not "
725 : "stored. Please set uDotDotOldRequested() to true in FEProblemBase before "
726 : "requesting `u_dotdot_old`.");
727 : }
728 :
729 : template <typename OutputType>
730 : const MooseArray<Number> &
731 0 : MooseVariableDataFV<OutputType>::dofValuesDuDotDu() const
732 : {
733 0 : _need_dof_du_dot_du = true;
734 0 : return _dof_du_dot_du;
735 : }
736 :
737 : template <typename OutputType>
738 : const MooseArray<Number> &
739 0 : MooseVariableDataFV<OutputType>::dofValuesDuDotDotDu() const
740 : {
741 0 : _need_dof_du_dotdot_du = true;
742 0 : return _dof_du_dotdot_du;
743 : }
744 :
745 : template <typename OutputType>
746 : void
747 0 : MooseVariableDataFV<OutputType>::fetchADDoFValues()
748 : {
749 0 : auto n = _dof_indices.size();
750 : libmesh_assert(n);
751 0 : _ad_dof_values.resize(n);
752 :
753 0 : const bool do_derivatives =
754 0 : ADReal::do_derivatives && _sys.number() == _subproblem.currentNlSysNum();
755 :
756 0 : for (decltype(n) i = 0; i < n; ++i)
757 : {
758 0 : _ad_dof_values[i] = _vector_tags_dof_u[_solution_tag][i];
759 0 : if (do_derivatives)
760 0 : Moose::derivInsert(_ad_dof_values[i].derivatives(), _dof_indices[i], 1.);
761 : }
762 0 : }
763 :
764 : template <typename OutputType>
765 : void
766 385530 : MooseVariableDataFV<OutputType>::prepareIC()
767 : {
768 : // TODO: implement this function
769 385530 : initDofIndices();
770 385530 : _vector_tags_dof_u[_solution_tag].resize(_dof_indices.size());
771 :
772 : mooseAssert(_qrule, "We should have a non-null qrule");
773 385530 : const auto nqp = _qrule->n_points();
774 385530 : _vector_tag_u[_solution_tag].resize(nqp);
775 385530 : }
776 :
777 : template <typename OutputType>
778 : void
779 1248 : MooseVariableDataFV<OutputType>::meshChanged()
780 : {
781 1248 : _prev_elem = nullptr;
782 1248 : }
783 :
784 : template class MooseVariableDataFV<Real>;
785 : // TODO: implement vector fv variable support. This will require some template
786 : // specializations for various member functions in this and the FV variable
787 : // classes. And then you will need to uncomment out the line below:
788 : // template class MooseVariableDataFV<RealVectorValue>;
|