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 "MooseLinearVariableFV.h"
11 : #include "TimeIntegrator.h"
12 : #include "NonlinearSystemBase.h"
13 : #include "DisplacedSystem.h"
14 : #include "SystemBase.h"
15 : #include "LinearSystem.h"
16 : #include "AuxiliarySystem.h"
17 : #include "SubProblem.h"
18 : #include "Assembly.h"
19 : #include "MathFVUtils.h"
20 : #include "FVUtils.h"
21 : #include "FVFluxBC.h"
22 : #include "FVDirichletBCBase.h"
23 : #include "GreenGaussGradient.h"
24 : #include "LinearFVBoundaryCondition.h"
25 : #include "LinearFVAdvectionDiffusionFunctorDirichletBC.h"
26 : #include "GradientLimiterType.h"
27 :
28 : #include "libmesh/numeric_vector.h"
29 :
30 : #include <climits>
31 : #include <typeinfo>
32 :
33 : using namespace Moose;
34 :
35 : registerMooseObject("MooseApp", MooseLinearVariableFVReal);
36 :
37 : namespace
38 : {
39 : const std::vector<std::unique_ptr<libMesh::NumericVector<libMesh::Number>>> &
40 1228 : linearFVGradientContainer(SystemBase & sys)
41 : {
42 1228 : if (auto * const linear_system = dynamic_cast<LinearSystem *>(&sys))
43 1032 : return linear_system->linearFVGradientContainer();
44 :
45 196 : if (auto * const auxiliary_system = dynamic_cast<AuxiliarySystem *>(&sys))
46 196 : return auxiliary_system->linearFVGradientContainer();
47 :
48 0 : mooseError("The assigned system is not a linear or an auxiliary system. Linear variables can "
49 : "only be assigned to linear or auxiliary systems.");
50 : }
51 : }
52 :
53 : template <typename OutputType>
54 : InputParameters
55 5517 : MooseLinearVariableFV<OutputType>::validParams()
56 : {
57 5517 : InputParameters params = MooseVariableField<OutputType>::validParams();
58 11034 : params.set<bool>("fv") = true;
59 22068 : params.set<MooseEnum>("family") = "MONOMIAL";
60 16551 : params.set<MooseEnum>("order") = "CONSTANT";
61 5517 : return params;
62 0 : }
63 :
64 : template <typename OutputType>
65 1228 : MooseLinearVariableFV<OutputType>::MooseLinearVariableFV(const InputParameters & parameters)
66 : : MooseVariableField<OutputType>(parameters),
67 1228 : _needs_cell_gradients(false),
68 1228 : _linear_system(dynamic_cast<LinearSystem *>(&this->_sys)),
69 1228 : _auxiliary_system(dynamic_cast<AuxiliarySystem *>(&this->_sys)),
70 1228 : _grad_container(linearFVGradientContainer(this->_sys)),
71 1228 : _sys_num(this->_sys.number()),
72 1228 : _solution(this->_sys.currentSolution()),
73 : // The following members are needed to be able to interface with the postprocessor and
74 : // auxiliary systems
75 1228 : _phi(this->_assembly.template fePhi<OutputShape>(FEType(CONSTANT, MONOMIAL))),
76 1228 : _grad_phi(this->_assembly.template feGradPhi<OutputShape>(FEType(CONSTANT, MONOMIAL))),
77 1228 : _phi_face(this->_assembly.template fePhiFace<OutputShape>(FEType(CONSTANT, MONOMIAL))),
78 1228 : _grad_phi_face(this->_assembly.template feGradPhiFace<OutputShape>(FEType(CONSTANT, MONOMIAL))),
79 1228 : _phi_face_neighbor(
80 1228 : this->_assembly.template fePhiFaceNeighbor<OutputShape>(FEType(CONSTANT, MONOMIAL))),
81 1228 : _grad_phi_face_neighbor(
82 1228 : this->_assembly.template feGradPhiFaceNeighbor<OutputShape>(FEType(CONSTANT, MONOMIAL))),
83 1228 : _phi_neighbor(this->_assembly.template fePhiNeighbor<OutputShape>(FEType(CONSTANT, MONOMIAL))),
84 1228 : _grad_phi_neighbor(
85 3684 : this->_assembly.template feGradPhiNeighbor<OutputShape>(FEType(CONSTANT, MONOMIAL)))
86 : {
87 1228 : if (!_linear_system && !_auxiliary_system)
88 0 : this->paramError("solver_sys",
89 : "The assigned system is not a linear or an auxiliary system! Linear variables "
90 : "can only be assigned to linear or auxiliary systems!");
91 1228 : _element_data = std::make_unique<MooseVariableDataLinearFV<OutputType>>(
92 1228 : *this, _sys, _tid, Moose::ElementType::Element, this->_assembly.elem());
93 1228 : _neighbor_data = std::make_unique<MooseVariableDataLinearFV<OutputType>>(
94 1228 : *this, _sys, _tid, Moose::ElementType::Neighbor, this->_assembly.neighbor());
95 :
96 1228 : if (libMesh::n_threads() > 1)
97 0 : mooseError("MooseLinearVariableFV does not support threading at the moment!");
98 1228 : }
99 :
100 : template <typename OutputType>
101 : void
102 189 : MooseLinearVariableFV<OutputType>::computeCellGradients(
103 : const Moose::FV::GradientLimiterType limiter_type)
104 : {
105 189 : if (limiter_type == Moose::FV::GradientLimiterType::None)
106 126 : computeCellGradients();
107 : else
108 63 : computeCellLimitedGradients(limiter_type);
109 189 : }
110 :
111 : template <typename OutputType>
112 : void
113 63 : MooseLinearVariableFV<OutputType>::computeCellLimitedGradients(
114 : const Moose::FV::GradientLimiterType limiter_type)
115 : {
116 63 : computeCellGradients();
117 :
118 63 : if (_linear_system)
119 63 : _linear_system->requestLinearFVLimitedGradients(limiter_type, this->_var_num);
120 : else
121 0 : _auxiliary_system->requestLinearFVLimitedGradients(limiter_type, this->_var_num);
122 63 : }
123 :
124 : template <typename OutputType>
125 : bool
126 0 : MooseLinearVariableFV<OutputType>::isExtrapolatedBoundaryFace(
127 : const FaceInfo & /*fi*/, const Elem * const /*elem*/, const Moose::StateArg & /*state*/) const
128 : {
129 : /// This is not used by this variable at this point.
130 0 : return false;
131 : }
132 :
133 : template <typename OutputType>
134 : Real
135 114270963 : MooseLinearVariableFV<OutputType>::getElemValue(const ElemInfo & elem_info,
136 : const StateArg & state) const
137 : {
138 : mooseAssert(
139 : this->hasBlocks(elem_info.subdomain_id()),
140 : "The variable should be defined on the element's subdomain! This typically occurs when the "
141 : "user wants to evaluate the elements right next to the boundary of two variables (block "
142 : "boundary). The subdomain which is queried: " +
143 : Moose::stringify(this->activeSubdomains()) + " the subdomain of the element " +
144 : std::to_string(elem_info.subdomain_id()));
145 :
146 : // It's not safe to use solutionState(0) because it returns the libMesh System solution member
147 : // which is wrong during things like finite difference Jacobian evaluation, e.g. when PETSc
148 : // perturbs the solution vector we feed these perturbations into the current_local_solution
149 : // while the libMesh solution is frozen in the non-perturbed state
150 228541926 : const auto & global_soln = (state.state == 0)
151 114270963 : ? *this->_sys.currentSolution()
152 0 : : this->_sys.solutionState(state.state, state.iteration_type);
153 :
154 114270963 : return global_soln(elem_info.dofIndices()[this->_sys_num][this->_var_num]);
155 : }
156 :
157 : template <typename OutputType>
158 : VectorValue<Real>
159 44963694 : MooseLinearVariableFV<OutputType>::gradSln(const ElemInfo & elem_info, const StateArg & state) const
160 : {
161 44963694 : if (state.state != 0)
162 0 : gradientStateError(state);
163 :
164 44963694 : if (_needs_cell_gradients)
165 : {
166 44963694 : _cell_gradient.zero();
167 134819185 : for (const auto i : make_range(this->_mesh.dimension()))
168 89855491 : _cell_gradient(i) =
169 89855491 : (*_grad_container[i])(elem_info.dofIndices()[this->_sys_num][this->_var_num]);
170 : }
171 :
172 44963694 : return _cell_gradient;
173 : }
174 :
175 : template <typename OutputType>
176 : Real
177 0 : MooseLinearVariableFV<OutputType>::gradSlnComponent(const ElemInfo & elem_info,
178 : const unsigned int component) const
179 : {
180 : mooseAssert(_needs_cell_gradients,
181 : "Gradient component requested without calling computeCellGradients().");
182 : mooseAssert(component < _grad_container.size(), "Gradient component index out of range.");
183 :
184 0 : return (*_grad_container[component])(elem_info.dofIndices()[this->_sys_num][this->_var_num]);
185 : }
186 :
187 : template <typename OutputType>
188 : VectorValue<Real>
189 70514970 : MooseLinearVariableFV<OutputType>::gradSln(const ElemInfo & elem_info,
190 : const StateArg & state,
191 : const Moose::FV::GradientLimiterType limiter_type) const
192 : {
193 : return (limiter_type == Moose::FV::GradientLimiterType::None)
194 70514970 : ? gradSln(elem_info, state)
195 70514970 : : limitedGradSln(elem_info, state, limiter_type);
196 : }
197 :
198 : template <typename OutputType>
199 : VectorValue<Real>
200 28016490 : MooseLinearVariableFV<OutputType>::limitedGradSln(
201 : const ElemInfo & elem_info,
202 : const StateArg & state,
203 : const Moose::FV::GradientLimiterType limiter_type) const
204 : {
205 28016490 : if (state.state != 0)
206 0 : gradientStateError(state);
207 :
208 28016490 : _cell_gradient.zero();
209 28016490 : const auto & limited_grad_container =
210 28016490 : _linear_system ? _linear_system->linearFVLimitedGradientContainer(limiter_type)
211 0 : : _auxiliary_system->linearFVLimitedGradientContainer(limiter_type);
212 84014928 : for (const auto i : make_range(this->_mesh.dimension()))
213 55998438 : _cell_gradient(i) =
214 55998438 : (*limited_grad_container[i])(elem_info.dofIndices()[this->_sys_num][this->_var_num]);
215 :
216 28016490 : return _cell_gradient;
217 : }
218 :
219 : template <typename OutputType>
220 : VectorValue<Real>
221 140 : MooseLinearVariableFV<OutputType>::gradSln(const FaceInfo & fi, const StateArg & state) const
222 : {
223 140 : const auto face_type = fi.faceType(std::make_pair(this->_var_num, this->_sys_num));
224 : mooseAssert(face_type != FaceInfo::VarFaceNeighbors::NEITHER,
225 : "Gradient requested on a face where the variable is defined on neither side.");
226 :
227 140 : const bool var_defined_on_elem = (face_type == FaceInfo::VarFaceNeighbors::BOTH) ||
228 : (face_type == FaceInfo::VarFaceNeighbors::ELEM);
229 140 : const auto * const elem_one = var_defined_on_elem ? fi.elemInfo() : fi.neighborInfo();
230 140 : const auto * const elem_two = var_defined_on_elem ? fi.neighborInfo() : fi.elemInfo();
231 :
232 140 : const auto elem_one_grad = gradSln(*elem_one, state);
233 :
234 : // If we have a neighbor then we interpolate between the two to the face.
235 140 : if (face_type == FaceInfo::VarFaceNeighbors::BOTH)
236 : {
237 : mooseAssert(elem_two, "Face type indicates BOTH but neighbor information is missing.");
238 0 : const auto elem_two_grad = gradSln(*elem_two, state);
239 0 : return Moose::FV::linearInterpolation(elem_one_grad, elem_two_grad, fi, var_defined_on_elem);
240 : }
241 : else
242 140 : return elem_one_grad;
243 : }
244 :
245 : template <typename OutputType>
246 : VectorValue<Real>
247 0 : MooseLinearVariableFV<OutputType>::gradSln(const FaceInfo & fi,
248 : const StateArg & state,
249 : const Moose::FV::GradientLimiterType limiter_type) const
250 : {
251 : return (limiter_type == Moose::FV::GradientLimiterType::None)
252 0 : ? gradSln(fi, state)
253 0 : : limitedGradSln(fi, state, limiter_type);
254 : }
255 :
256 : template <typename OutputType>
257 : VectorValue<Real>
258 0 : MooseLinearVariableFV<OutputType>::limitedGradSln(
259 : const FaceInfo & fi,
260 : const StateArg & state,
261 : const Moose::FV::GradientLimiterType limiter_type) const
262 : {
263 0 : const auto face_type = fi.faceType(std::make_pair(this->_var_num, this->_sys_num));
264 : mooseAssert(face_type != FaceInfo::VarFaceNeighbors::NEITHER,
265 : "Limited gradient requested on a face where the variable is defined on neither "
266 : "side.");
267 :
268 0 : const bool var_defined_on_elem = (face_type == FaceInfo::VarFaceNeighbors::BOTH) ||
269 : (face_type == FaceInfo::VarFaceNeighbors::ELEM);
270 0 : const auto * const elem_one = var_defined_on_elem ? fi.elemInfo() : fi.neighborInfo();
271 0 : const auto * const elem_two = var_defined_on_elem ? fi.neighborInfo() : fi.elemInfo();
272 :
273 0 : const auto elem_one_grad = limitedGradSln(*elem_one, state, limiter_type);
274 :
275 0 : if (face_type == FaceInfo::VarFaceNeighbors::BOTH)
276 : {
277 : mooseAssert(elem_two, "Face type indicates BOTH but neighbor information is missing.");
278 0 : const auto elem_two_grad = limitedGradSln(*elem_two, state, limiter_type);
279 0 : return Moose::FV::linearInterpolation(elem_one_grad, elem_two_grad, fi, var_defined_on_elem);
280 : }
281 : else
282 0 : return elem_one_grad;
283 : }
284 :
285 : template <typename OutputType>
286 : void
287 1226 : MooseLinearVariableFV<OutputType>::initialSetup()
288 : {
289 1226 : MooseVariableField<OutputType>::initialSetup();
290 1226 : cacheBoundaryBCMap();
291 1226 : }
292 :
293 : template <typename OutputType>
294 : void
295 1388 : MooseLinearVariableFV<OutputType>::timestepSetup()
296 : {
297 1388 : MooseVariableField<OutputType>::timestepSetup();
298 1388 : cacheBoundaryBCMap();
299 1388 : }
300 :
301 : template <typename OutputType>
302 : typename MooseLinearVariableFV<OutputType>::ValueType
303 1578 : MooseLinearVariableFV<OutputType>::evaluate(const FaceArg & face, const StateArg & state) const
304 : {
305 1578 : const FaceInfo * const fi = face.fi;
306 :
307 : mooseAssert(fi, "The face information must be non-null");
308 :
309 1578 : const auto face_type = fi->faceType(std::make_pair(this->_var_num, this->_sys_num));
310 :
311 1578 : if (face_type == FaceInfo::VarFaceNeighbors::BOTH)
312 14 : return Moose::FV::interpolate(*this, face, state);
313 1564 : else if (auto * bc_pointer = this->getBoundaryCondition(*fi->boundaryIDs().begin()))
314 : {
315 : mooseAssert(fi->boundaryIDs().size() == 1, "We should only have one boundary on every face.");
316 1557 : bc_pointer->setupFaceData(fi, face_type);
317 1557 : return bc_pointer->computeBoundaryValue();
318 : }
319 : // If no boundary condition is defined but we are evaluating on a boundary, just return the
320 : // element value
321 7 : else if (face_type == FaceInfo::VarFaceNeighbors::ELEM)
322 : {
323 7 : const auto & elem_info = this->_mesh.elemInfo(fi->elemPtr()->id());
324 7 : return getElemValue(elem_info, state);
325 : }
326 0 : else if (face_type == FaceInfo::VarFaceNeighbors::NEIGHBOR)
327 : {
328 0 : const auto & elem_info = this->_mesh.elemInfo(fi->neighborPtr()->id());
329 0 : return getElemValue(elem_info, state);
330 : }
331 : else
332 0 : mooseError("We should never get here!");
333 : }
334 :
335 : template <typename OutputType>
336 : typename MooseLinearVariableFV<OutputType>::ValueType
337 0 : MooseLinearVariableFV<OutputType>::evaluate(const NodeArg & /*node_arg*/,
338 : const StateArg & /*state*/) const
339 : {
340 0 : mooseError("Not implemented yet");
341 : }
342 :
343 : template <typename OutputType>
344 : typename MooseLinearVariableFV<OutputType>::DotType
345 : MooseLinearVariableFV<OutputType>::evaluateDot(const ElemArg &, const StateArg &) const
346 : {
347 : timeIntegratorError();
348 : }
349 :
350 : template <>
351 : ADReal
352 0 : MooseLinearVariableFV<Real>::evaluateDot(const ElemArg & /*elem_arg*/,
353 : const StateArg & /*state*/) const
354 : {
355 0 : timeIntegratorError();
356 : }
357 :
358 : template <typename OutputType>
359 : void
360 2614 : MooseLinearVariableFV<OutputType>::cacheBoundaryBCMap()
361 : {
362 2614 : _boundary_id_to_bc.clear();
363 2614 : std::vector<LinearFVBoundaryCondition *> bcs;
364 :
365 : // I believe because query() returns by value but condition returns by reference that binding to a
366 : // const lvalue reference results in the query() getting destructed and us holding onto a dangling
367 : // reference. I think that condition returned by value we would be able to bind to a const lvalue
368 : // reference here. But as it is we'll bind to a regular lvalue
369 2614 : auto base_query = this->_subproblem.getMooseApp()
370 2614 : .theWarehouse()
371 : .query()
372 2614 : .template condition<AttribSystem>("LinearFVBoundaryCondition")
373 2614 : .template condition<AttribThread>(_tid)
374 2614 : .template condition<AttribVar>(_var_num)
375 2614 : .template condition<AttribSysNum>(this->_sys.number());
376 :
377 10765 : for (const auto bnd_id : this->_mesh.getBoundaryIDs())
378 : {
379 8151 : auto base_query_copy = base_query;
380 16302 : base_query_copy.template condition<AttribBoundaries>(std::set<BoundaryID>({bnd_id}))
381 8151 : .queryInto(bcs);
382 : mooseAssert(bcs.size() <= 1, "cannot have multiple BCs on the same boundary");
383 8151 : if (!bcs.empty())
384 6767 : _boundary_id_to_bc.emplace(bnd_id, bcs[0]);
385 : }
386 2614 : }
387 :
388 : template <typename OutputType>
389 : LinearFVBoundaryCondition *
390 5311658 : MooseLinearVariableFV<OutputType>::getBoundaryCondition(const BoundaryID bd_id) const
391 : {
392 5311658 : const auto iter = _boundary_id_to_bc.find(bd_id);
393 5311658 : if (iter == _boundary_id_to_bc.end())
394 23367 : return nullptr;
395 : else
396 5288291 : return iter->second;
397 : }
398 :
399 : template <typename OutputType>
400 : const Elem * const &
401 0 : MooseLinearVariableFV<OutputType>::currentElem() const
402 : {
403 0 : return this->_assembly.elem();
404 : }
405 :
406 : template <typename OutputType>
407 : bool
408 0 : MooseLinearVariableFV<OutputType>::isDirichletBoundaryFace(const FaceInfo & fi) const
409 : {
410 0 : for (const auto bnd_id : fi.boundaryIDs())
411 0 : if (auto it = _boundary_id_to_bc.find(bnd_id); it != _boundary_id_to_bc.end())
412 0 : if (dynamic_cast<LinearFVAdvectionDiffusionFunctorDirichletBC *>(it->second))
413 0 : return true;
414 :
415 0 : return false;
416 : }
417 :
418 : // ****************************************************************************
419 : // The functions below are used for interfacing the auxiliary and
420 : // postprocessor/userobject systems. Most of the postprocessors/
421 : // auxkernels require quadrature-based evaluations and we provide that
422 : // interface with the functions below.
423 : // ****************************************************************************
424 :
425 : template <typename OutputType>
426 : void
427 0 : MooseLinearVariableFV<OutputType>::getDofIndices(const Elem * elem,
428 : std::vector<dof_id_type> & dof_indices) const
429 : {
430 0 : dof_indices.clear();
431 0 : const auto & elem_info = this->_mesh.elemInfo(elem->id());
432 0 : dof_indices.push_back(elem_info.dofIndices()[this->_sys_num][this->number()]);
433 0 : }
434 :
435 : template <typename OutputType>
436 : void
437 378 : MooseLinearVariableFV<OutputType>::setDofValue(const DofValue & value, unsigned int index)
438 : {
439 378 : _element_data->setDofValue(value, index);
440 378 : }
441 :
442 : template <typename OutputType>
443 : void
444 126 : MooseLinearVariableFV<OutputType>::setDofValues(const DenseVector<DofValue> & values)
445 : {
446 126 : _element_data->setDofValues(values);
447 126 : }
448 :
449 : template <typename OutputType>
450 : void
451 11922 : MooseLinearVariableFV<OutputType>::clearDofIndices()
452 : {
453 11922 : _element_data->clearDofIndices();
454 11922 : }
455 :
456 : template <typename OutputType>
457 : const typename MooseLinearVariableFV<OutputType>::DofValues &
458 0 : MooseLinearVariableFV<OutputType>::dofValues() const
459 : {
460 0 : return _element_data->dofValues();
461 : }
462 :
463 : template <typename OutputType>
464 : const typename MooseLinearVariableFV<OutputType>::DofValues &
465 0 : MooseLinearVariableFV<OutputType>::dofValuesNeighbor() const
466 : {
467 0 : return _neighbor_data->dofValues();
468 : }
469 :
470 : template <typename OutputType>
471 : const typename MooseLinearVariableFV<OutputType>::DofValues &
472 0 : MooseLinearVariableFV<OutputType>::dofValuesOld() const
473 : {
474 0 : return _element_data->dofValuesOld();
475 : }
476 :
477 : template <typename OutputType>
478 : const typename MooseLinearVariableFV<OutputType>::DofValues &
479 0 : MooseLinearVariableFV<OutputType>::dofValuesOlder() const
480 : {
481 0 : return _element_data->dofValuesOlder();
482 : }
483 :
484 : template <typename OutputType>
485 : const typename MooseLinearVariableFV<OutputType>::DofValues &
486 0 : MooseLinearVariableFV<OutputType>::dofValuesPreviousNL() const
487 : {
488 0 : return _element_data->dofValuesPreviousNL();
489 : }
490 :
491 : template <typename OutputType>
492 : const typename MooseLinearVariableFV<OutputType>::DofValues &
493 0 : MooseLinearVariableFV<OutputType>::dofValuesOldNeighbor() const
494 : {
495 0 : return _neighbor_data->dofValuesOld();
496 : }
497 :
498 : template <typename OutputType>
499 : const typename MooseLinearVariableFV<OutputType>::DofValues &
500 0 : MooseLinearVariableFV<OutputType>::dofValuesOlderNeighbor() const
501 : {
502 0 : return _neighbor_data->dofValuesOlder();
503 : }
504 :
505 : template <typename OutputType>
506 : const typename MooseLinearVariableFV<OutputType>::DofValues &
507 0 : MooseLinearVariableFV<OutputType>::dofValuesPreviousNLNeighbor() const
508 : {
509 0 : return _neighbor_data->dofValuesPreviousNL();
510 : }
511 :
512 : template <typename OutputType>
513 : const typename MooseLinearVariableFV<OutputType>::DofValues &
514 0 : MooseLinearVariableFV<OutputType>::dofValuesDot() const
515 : {
516 0 : timeIntegratorError();
517 : }
518 :
519 : template <typename OutputType>
520 : const typename MooseLinearVariableFV<OutputType>::DofValues &
521 0 : MooseLinearVariableFV<OutputType>::dofValuesDotDot() const
522 : {
523 0 : timeIntegratorError();
524 : }
525 :
526 : template <typename OutputType>
527 : const typename MooseLinearVariableFV<OutputType>::DofValues &
528 0 : MooseLinearVariableFV<OutputType>::dofValuesDotOld() const
529 : {
530 0 : timeIntegratorError();
531 : }
532 :
533 : template <typename OutputType>
534 : const typename MooseLinearVariableFV<OutputType>::DofValues &
535 0 : MooseLinearVariableFV<OutputType>::dofValuesDotDotOld() const
536 : {
537 0 : timeIntegratorError();
538 : }
539 :
540 : template <typename OutputType>
541 : const typename MooseLinearVariableFV<OutputType>::DofValues &
542 0 : MooseLinearVariableFV<OutputType>::dofValuesDotNeighbor() const
543 : {
544 0 : timeIntegratorError();
545 : }
546 :
547 : template <typename OutputType>
548 : const typename MooseLinearVariableFV<OutputType>::DofValues &
549 0 : MooseLinearVariableFV<OutputType>::dofValuesDotDotNeighbor() const
550 : {
551 0 : timeIntegratorError();
552 : }
553 :
554 : template <typename OutputType>
555 : const typename MooseLinearVariableFV<OutputType>::DofValues &
556 0 : MooseLinearVariableFV<OutputType>::dofValuesDotOldNeighbor() const
557 : {
558 0 : timeIntegratorError();
559 : }
560 :
561 : template <typename OutputType>
562 : const typename MooseLinearVariableFV<OutputType>::DofValues &
563 0 : MooseLinearVariableFV<OutputType>::dofValuesDotDotOldNeighbor() const
564 : {
565 0 : timeIntegratorError();
566 : }
567 :
568 : template <typename OutputType>
569 : const MooseArray<Number> &
570 0 : MooseLinearVariableFV<OutputType>::dofValuesDuDotDu() const
571 : {
572 0 : timeIntegratorError();
573 : }
574 :
575 : template <typename OutputType>
576 : const MooseArray<Number> &
577 0 : MooseLinearVariableFV<OutputType>::dofValuesDuDotDotDu() const
578 : {
579 0 : timeIntegratorError();
580 : }
581 :
582 : template <typename OutputType>
583 : const MooseArray<Number> &
584 0 : MooseLinearVariableFV<OutputType>::dofValuesDuDotDuNeighbor() const
585 : {
586 0 : timeIntegratorError();
587 : }
588 :
589 : template <typename OutputType>
590 : const MooseArray<Number> &
591 0 : MooseLinearVariableFV<OutputType>::dofValuesDuDotDotDuNeighbor() const
592 : {
593 0 : timeIntegratorError();
594 : }
595 :
596 : template <typename OutputType>
597 : void
598 575485 : MooseLinearVariableFV<OutputType>::computeElemValues()
599 : {
600 575485 : _element_data->setGeometry(Moose::Volume);
601 575485 : _element_data->computeValues();
602 575485 : }
603 :
604 : template <typename OutputType>
605 : void
606 405 : MooseLinearVariableFV<OutputType>::computeElemValuesFace()
607 : {
608 405 : _element_data->setGeometry(Moose::Face);
609 405 : _element_data->computeValues();
610 405 : }
611 :
612 : template <typename OutputType>
613 : void
614 0 : MooseLinearVariableFV<OutputType>::computeNeighborValuesFace()
615 : {
616 0 : _neighbor_data->setGeometry(Moose::Face);
617 0 : _neighbor_data->computeValues();
618 0 : }
619 :
620 : template <typename OutputType>
621 : void
622 0 : MooseLinearVariableFV<OutputType>::computeNeighborValues()
623 : {
624 0 : _neighbor_data->setGeometry(Moose::Volume);
625 0 : _neighbor_data->computeValues();
626 0 : }
627 :
628 : template <typename OutputType>
629 : void
630 0 : MooseLinearVariableFV<OutputType>::computeLowerDValues()
631 : {
632 0 : lowerDError();
633 : }
634 :
635 : template <typename OutputType>
636 : void
637 0 : MooseLinearVariableFV<OutputType>::computeNodalNeighborValues()
638 : {
639 0 : nodalError();
640 : }
641 :
642 : template <typename OutputType>
643 : void
644 0 : MooseLinearVariableFV<OutputType>::computeNodalValues()
645 : {
646 0 : nodalError();
647 : }
648 :
649 : template <typename OutputType>
650 : const std::vector<dof_id_type> &
651 0 : MooseLinearVariableFV<OutputType>::dofIndices() const
652 : {
653 0 : return _element_data->dofIndices();
654 : }
655 :
656 : template <typename OutputType>
657 : const std::vector<dof_id_type> &
658 0 : MooseLinearVariableFV<OutputType>::dofIndicesNeighbor() const
659 : {
660 0 : return _neighbor_data->dofIndices();
661 : }
662 :
663 : template <typename OutputType>
664 : void
665 0 : MooseLinearVariableFV<OutputType>::setLowerDofValues(const DenseVector<DofValue> &)
666 : {
667 0 : lowerDError();
668 : }
669 :
670 : template <typename OutputType>
671 : unsigned int
672 1228 : MooseLinearVariableFV<OutputType>::oldestSolutionStateRequested() const
673 : {
674 1228 : unsigned int state = 0;
675 1228 : state = std::max(state, _element_data->oldestSolutionStateRequested());
676 1228 : state = std::max(state, _neighbor_data->oldestSolutionStateRequested());
677 1228 : return state;
678 : }
679 :
680 : template <typename OutputType>
681 : void
682 18 : MooseLinearVariableFV<OutputType>::clearAllDofIndices()
683 : {
684 18 : _element_data->clearDofIndices();
685 18 : _neighbor_data->clearDofIndices();
686 18 : }
687 :
688 : template <typename OutputType>
689 : void
690 0 : MooseLinearVariableFV<OutputType>::setNodalValue(const OutputType & /*value*/, unsigned int /*idx*/)
691 : {
692 0 : nodalError();
693 : }
694 :
695 : template <typename OutputType>
696 : const typename MooseLinearVariableFV<OutputType>::DofValues &
697 0 : MooseLinearVariableFV<OutputType>::nodalVectorTagValue(TagID) const
698 : {
699 0 : nodalError();
700 : }
701 :
702 : template <typename OutputType>
703 : const typename MooseLinearVariableFV<OutputType>::DofValues &
704 0 : MooseLinearVariableFV<OutputType>::nodalMatrixTagValue(TagID) const
705 : {
706 0 : nodalError();
707 : }
708 :
709 : template <typename OutputType>
710 : const std::vector<dof_id_type> &
711 0 : MooseLinearVariableFV<OutputType>::dofIndicesLower() const
712 : {
713 0 : lowerDError();
714 : }
715 :
716 : template <typename OutputType>
717 : const typename MooseLinearVariableFV<OutputType>::FieldVariablePhiValue &
718 0 : MooseLinearVariableFV<OutputType>::phiLower() const
719 : {
720 0 : lowerDError();
721 : }
722 :
723 : template <typename OutputType>
724 : void
725 504 : MooseLinearVariableFV<OutputType>::insert(NumericVector<Number> & vector)
726 : {
727 504 : _element_data->insert(vector);
728 504 : }
729 :
730 : template <typename OutputType>
731 : void
732 0 : MooseLinearVariableFV<OutputType>::insertLower(NumericVector<Number> & /*residual*/)
733 : {
734 0 : mooseError("We don't support value insertion to residuals in MooseLinearVariableFV!");
735 : }
736 :
737 : template <typename OutputType>
738 : void
739 0 : MooseLinearVariableFV<OutputType>::add(NumericVector<Number> & /*residual*/)
740 : {
741 0 : mooseError("We don't support value addition to residuals in MooseLinearVariableFV!");
742 : }
743 :
744 : template <typename OutputType>
745 : void
746 2332 : MooseLinearVariableFV<OutputType>::setActiveTags(const std::set<TagID> & vtags)
747 : {
748 2332 : _element_data->setActiveTags(vtags);
749 2332 : _neighbor_data->setActiveTags(vtags);
750 2332 : }
751 :
752 : template <typename OutputType>
753 : const MooseArray<OutputType> &
754 0 : MooseLinearVariableFV<OutputType>::nodalValueArray() const
755 : {
756 0 : nodalError();
757 : }
758 :
759 : template <typename OutputType>
760 : const MooseArray<OutputType> &
761 0 : MooseLinearVariableFV<OutputType>::nodalValueOldArray() const
762 : {
763 0 : nodalError();
764 : }
765 :
766 : template <typename OutputType>
767 : const MooseArray<OutputType> &
768 0 : MooseLinearVariableFV<OutputType>::nodalValueOlderArray() const
769 : {
770 0 : nodalError();
771 : }
772 :
773 : template <typename OutputType>
774 : std::size_t
775 0 : MooseLinearVariableFV<OutputType>::phiLowerSize() const
776 : {
777 0 : lowerDError();
778 : }
779 :
780 : template <typename OutputType>
781 : const typename MooseLinearVariableFV<OutputType>::FieldVariableValue &
782 66 : MooseLinearVariableFV<OutputType>::vectorTagValue(TagID tag) const
783 : {
784 66 : return _element_data->vectorTagValue(tag);
785 : }
786 :
787 : template <typename OutputType>
788 : const typename MooseLinearVariableFV<OutputType>::DofValues &
789 66 : MooseLinearVariableFV<OutputType>::vectorTagDofValue(TagID tag) const
790 : {
791 66 : return _element_data->vectorTagDofValue(tag);
792 : }
793 :
794 : template <typename OutputType>
795 : const typename MooseLinearVariableFV<OutputType>::FieldVariableValue &
796 22 : MooseLinearVariableFV<OutputType>::matrixTagValue(TagID tag) const
797 : {
798 22 : return _element_data->matrixTagValue(tag);
799 : }
800 :
801 : template <typename OutputType>
802 : const typename MooseLinearVariableFV<OutputType>::FieldVariablePhiSecond &
803 0 : MooseLinearVariableFV<OutputType>::secondPhi() const
804 : {
805 0 : mooseError("We don't currently implement second derivatives for FV");
806 : }
807 :
808 : template <typename OutputType>
809 : const typename MooseLinearVariableFV<OutputType>::FieldVariablePhiValue &
810 0 : MooseLinearVariableFV<OutputType>::curlPhi() const
811 : {
812 0 : mooseError("We don't currently implement curl for FV");
813 : }
814 :
815 : template <typename OutputType>
816 : const typename MooseLinearVariableFV<OutputType>::FieldVariablePhiDivergence &
817 0 : MooseLinearVariableFV<OutputType>::divPhi() const
818 : {
819 0 : mooseError("We don't currently implement divergence for FV");
820 : }
821 :
822 : template <typename OutputType>
823 : const typename MooseLinearVariableFV<OutputType>::FieldVariablePhiSecond &
824 0 : MooseLinearVariableFV<OutputType>::secondPhiFace() const
825 : {
826 0 : mooseError("We don't currently implement second derivatives for FV");
827 : }
828 :
829 : template <typename OutputType>
830 : const typename MooseLinearVariableFV<OutputType>::FieldVariablePhiSecond &
831 0 : MooseLinearVariableFV<OutputType>::secondPhiFaceNeighbor() const
832 : {
833 0 : mooseError("We don't currently implement second derivatives for FV");
834 : }
835 :
836 : template <typename OutputType>
837 : const typename MooseLinearVariableFV<OutputType>::FieldVariablePhiSecond &
838 0 : MooseLinearVariableFV<OutputType>::secondPhiNeighbor() const
839 : {
840 0 : mooseError("We don't currently implement second derivatives for FV");
841 : }
842 :
843 : template <typename OutputType>
844 : const typename MooseLinearVariableFV<OutputType>::FieldVariableValue &
845 11011 : MooseLinearVariableFV<OutputType>::sln() const
846 : {
847 11011 : return _element_data->sln(Moose::Current);
848 : }
849 :
850 : template <typename OutputType>
851 : const typename MooseLinearVariableFV<OutputType>::FieldVariableValue &
852 0 : MooseLinearVariableFV<OutputType>::slnOld() const
853 : {
854 0 : return _element_data->sln(Moose::Old);
855 : }
856 :
857 : template <typename OutputType>
858 : const typename MooseLinearVariableFV<OutputType>::FieldVariableValue &
859 0 : MooseLinearVariableFV<OutputType>::slnOlder() const
860 : {
861 0 : return _element_data->sln(Moose::Older);
862 : }
863 :
864 : template <typename OutputType>
865 : const typename MooseLinearVariableFV<OutputType>::FieldVariableGradient &
866 102 : MooseLinearVariableFV<OutputType>::gradSln() const
867 : {
868 102 : return _element_data->gradSln(Moose::Current);
869 : }
870 :
871 : template <typename OutputType>
872 : const typename MooseLinearVariableFV<OutputType>::FieldVariableGradient &
873 0 : MooseLinearVariableFV<OutputType>::gradSlnOld() const
874 : {
875 0 : return _element_data->gradSln(Moose::Old);
876 : }
877 :
878 : template <typename OutputType>
879 : const typename MooseLinearVariableFV<OutputType>::FieldVariableValue &
880 0 : MooseLinearVariableFV<OutputType>::slnNeighbor() const
881 : {
882 0 : return _neighbor_data->sln(Moose::Current);
883 : }
884 :
885 : template <typename OutputType>
886 : const typename MooseLinearVariableFV<OutputType>::FieldVariableValue &
887 0 : MooseLinearVariableFV<OutputType>::slnOldNeighbor() const
888 : {
889 0 : return _neighbor_data->sln(Moose::Old);
890 : }
891 :
892 : template <typename OutputType>
893 : const typename MooseLinearVariableFV<OutputType>::FieldVariableGradient &
894 0 : MooseLinearVariableFV<OutputType>::gradSlnNeighbor() const
895 : {
896 0 : return _neighbor_data->gradSln(Moose::Current);
897 : }
898 :
899 : template <typename OutputType>
900 : const typename MooseLinearVariableFV<OutputType>::FieldVariableGradient &
901 0 : MooseLinearVariableFV<OutputType>::gradSlnOldNeighbor() const
902 : {
903 0 : return _neighbor_data->gradSln(Moose::Old);
904 : }
905 :
906 : template <typename OutputType>
907 : const ADTemplateVariableSecond<OutputType> &
908 0 : MooseLinearVariableFV<OutputType>::adSecondSln() const
909 : {
910 0 : adError();
911 : }
912 :
913 : template <typename OutputType>
914 : const ADTemplateVariableValue<OutputType> &
915 0 : MooseLinearVariableFV<OutputType>::adUDot() const
916 : {
917 0 : adError();
918 : }
919 :
920 : template <typename OutputType>
921 : const ADTemplateVariableValue<OutputType> &
922 0 : MooseLinearVariableFV<OutputType>::adUDotDot() const
923 : {
924 0 : adError();
925 : }
926 :
927 : template <typename OutputType>
928 : const ADTemplateVariableGradient<OutputType> &
929 0 : MooseLinearVariableFV<OutputType>::adGradSlnDot() const
930 : {
931 0 : adError();
932 : }
933 :
934 : template <typename OutputType>
935 : const ADTemplateVariableValue<OutputType> &
936 0 : MooseLinearVariableFV<OutputType>::adSlnNeighbor() const
937 : {
938 0 : adError();
939 : }
940 :
941 : template <typename OutputType>
942 : const ADTemplateVariableGradient<OutputType> &
943 0 : MooseLinearVariableFV<OutputType>::adGradSlnNeighbor() const
944 : {
945 0 : adError();
946 : }
947 :
948 : template <typename OutputType>
949 : const ADTemplateVariableSecond<OutputType> &
950 0 : MooseLinearVariableFV<OutputType>::adSecondSlnNeighbor() const
951 : {
952 0 : adError();
953 : }
954 :
955 : template <typename OutputType>
956 : const ADTemplateVariableValue<OutputType> &
957 0 : MooseLinearVariableFV<OutputType>::adUDotNeighbor() const
958 : {
959 0 : adError();
960 : }
961 :
962 : template <typename OutputType>
963 : const ADTemplateVariableValue<OutputType> &
964 0 : MooseLinearVariableFV<OutputType>::adUDotDotNeighbor() const
965 : {
966 0 : adError();
967 : }
968 :
969 : template <typename OutputType>
970 : const ADTemplateVariableGradient<OutputType> &
971 0 : MooseLinearVariableFV<OutputType>::adGradSlnNeighborDot() const
972 : {
973 0 : adError();
974 : }
975 :
976 : template <typename OutputType>
977 : const ADTemplateVariableValue<OutputType> &
978 0 : MooseLinearVariableFV<OutputType>::adSln() const
979 : {
980 0 : adError();
981 : }
982 :
983 : template <typename OutputType>
984 : const ADTemplateVariableGradient<OutputType> &
985 0 : MooseLinearVariableFV<OutputType>::adGradSln() const
986 : {
987 0 : adError();
988 : }
989 :
990 : template <typename OutputType>
991 : const ADTemplateVariableCurl<OutputType> &
992 0 : MooseLinearVariableFV<OutputType>::adCurlSln() const
993 : {
994 0 : adError();
995 : }
996 :
997 : template <typename OutputType>
998 : const ADTemplateVariableCurl<OutputType> &
999 0 : MooseLinearVariableFV<OutputType>::adCurlSlnNeighbor() const
1000 : {
1001 0 : adError();
1002 : }
1003 :
1004 : template <typename OutputType>
1005 : const typename MooseLinearVariableFV<OutputType>::ADDofValues &
1006 0 : MooseLinearVariableFV<OutputType>::adDofValues() const
1007 : {
1008 0 : adError();
1009 : }
1010 :
1011 : template <typename OutputType>
1012 : const typename MooseLinearVariableFV<OutputType>::ADDofValues &
1013 0 : MooseLinearVariableFV<OutputType>::adDofValuesNeighbor() const
1014 : {
1015 0 : adError();
1016 : }
1017 :
1018 : template <typename OutputType>
1019 : const typename MooseLinearVariableFV<OutputType>::ADDofValues &
1020 0 : MooseLinearVariableFV<OutputType>::adDofValuesDot() const
1021 : {
1022 0 : adError();
1023 : }
1024 :
1025 : template <typename OutputType>
1026 : const dof_id_type &
1027 0 : MooseLinearVariableFV<OutputType>::nodalDofIndex() const
1028 : {
1029 0 : nodalError();
1030 : }
1031 :
1032 : template <typename OutputType>
1033 : const dof_id_type &
1034 0 : MooseLinearVariableFV<OutputType>::nodalDofIndexNeighbor() const
1035 : {
1036 0 : nodalError();
1037 : }
1038 :
1039 : template <typename OutputType>
1040 : void
1041 1228 : MooseLinearVariableFV<OutputType>::sizeMatrixTagData()
1042 : {
1043 1228 : _element_data->sizeMatrixTagData();
1044 1228 : }
1045 :
1046 : template class MooseLinearVariableFV<Real>;
|