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 "ParallelUniqueId.h"
13 : #include "MooseMesh.h"
14 : #include "MooseTypes.h"
15 : #include "MooseException.h"
16 : #include "MooseVariableFieldBase.h"
17 : #include "FVFluxKernel.h"
18 : #include "FVFluxBC.h"
19 : #include "FVInterfaceKernel.h"
20 : #include "FEProblem.h"
21 : #include "DisplacedProblem.h"
22 : #include "SwapBackSentinel.h"
23 : #include "MaterialBase.h"
24 : #include "libmesh/libmesh_exceptions.h"
25 : #include "libmesh/elem.h"
26 :
27 : // C++
28 : #include <cstring> // for "Jacobian" exception test
29 : #include <set>
30 :
31 : class MooseVariableFVBase;
32 :
33 : /// This works like the Sentinel helper classes in MOOSE, except it is more
34 : /// flexible and concise to use. You just initialize it with a lambda and can
35 : /// return it from within another function scope. So a function that has other
36 : /// cleanup functions associated with it can just wrap those cleanup funcs in
37 : /// an OnScopeExit object and return it. The caller of the function needing
38 : /// cleanup simply catches its return value and voila - the cleanup function
39 : /// are called at the end of the calling functions scope. Like this:
40 : ///
41 : /// OnScopeExit doSomethingSpecial()
42 : /// {
43 : /// ...
44 : /// return OnScopeExit([]{cleanupSomethingSpecial();});
45 : /// }
46 : /// void cleanupSomethingSpecial() {...}
47 : ///
48 : /// void main()
49 : /// {
50 : /// auto cleaner_uper = doSomethingSpecial();
51 : /// }
52 : struct OnScopeExit
53 : {
54 : std::function<void()> _f;
55 : OnScopeExit(std::function<void()> f) noexcept : _f(std::move(f)) {}
56 : OnScopeExit(OnScopeExit && other) : _f(std::move(other._f)) {}
57 : OnScopeExit & operator=(OnScopeExit && other)
58 : {
59 : _f = std::move(other._f);
60 : return *this;
61 : }
62 : ~OnScopeExit()
63 : {
64 : if (_f)
65 : _f();
66 : }
67 : };
68 :
69 : /**
70 : * This loops over a set of mesh faces (i.e. FaceInfo objects). Callback
71 : * routines are provided for visiting each face, for visiting boundary faces,
72 : * for sudomain changes, and pre/post many of these events.
73 : */
74 : template <typename RangeType>
75 : class ThreadedFaceLoop
76 : {
77 : public:
78 : ThreadedFaceLoop(FEProblemBase & fe_problem,
79 : const unsigned int nl_system_num,
80 : const std::set<TagID> & tags,
81 : bool on_displaced);
82 :
83 : ThreadedFaceLoop(ThreadedFaceLoop & x, Threads::split split);
84 :
85 : virtual ~ThreadedFaceLoop();
86 :
87 : virtual void operator()(const RangeType & range, bool bypass_threading = false);
88 :
89 : void join(const ThreadedFaceLoop & y);
90 :
91 : virtual void onFace(const FaceInfo & fi) = 0;
92 : /// This is called once for each face after all face and boundary callbacks have been
93 : /// finished for that face.
94 0 : virtual void postFace(const FaceInfo & /*fi*/) {}
95 : /// This is called once before all face-looping
96 0 : virtual void pre() {}
97 : /// This is called once after all face-looping is finished.
98 0 : virtual void post() {}
99 :
100 : /// This is called once for every face that is on a boundary *after* onFace
101 : /// is called for the face.
102 : virtual void onBoundary(const FaceInfo & fi, BoundaryID boundary) = 0;
103 :
104 : /// Called every time the current subdomain changes (i.e. the subdomain of *this* face's elem element
105 : /// is not the same as the subdomain of the last face's elem element).
106 148941 : virtual void subdomainChanged() { _fe_problem.subdomainSetup(_subdomain, _tid); }
107 :
108 : /// Called every time the neighbor subdomain changes (i.e. the subdomain of *this* face's neighbor element
109 : /// is not the same as the subdomain of the last face's neighbor element).
110 2901180 : virtual void neighborSubdomainChanged()
111 : {
112 2901180 : _fe_problem.neighborSubdomainSetup(_neighbor_subdomain, _tid);
113 2901180 : }
114 :
115 : /// Called if a MooseException is caught anywhere during the computation.
116 0 : void caughtMooseException(MooseException & e)
117 : {
118 0 : Threads::spin_mutex::scoped_lock lock(threaded_element_mutex);
119 0 : std::string what(e.what());
120 0 : _fe_problem.setException(what);
121 0 : }
122 :
123 : protected:
124 : /// Print list of object types executed and in which order
125 0 : virtual void printGeneralExecutionInformation() const {}
126 :
127 : /// Print ordering of objects executed on each block
128 0 : virtual void printBlockExecutionInformation() const {}
129 :
130 : /// Print ordering of objects exected on each boundary
131 0 : virtual void printBoundaryExecutionInformation(const BoundaryID /* bnd_id */) const {}
132 :
133 : /// Reset lists of blocks and boundaries for which execution printing has been done
134 102815 : void resetExecutionPrinting()
135 : {
136 102815 : _blocks_exec_printed.clear();
137 102815 : _boundaries_exec_printed.clear();
138 102815 : }
139 :
140 : FEProblemBase & _fe_problem;
141 : MooseMesh & _mesh;
142 : const std::set<TagID> & _tags;
143 : THREAD_ID _tid;
144 : const unsigned int _nl_system_num;
145 :
146 : /// Whether this loop is operating on the displaced mesh
147 : const bool _on_displaced;
148 :
149 : /// FEProblemBase or DisplacedProblem depending on \p _on_displaced
150 : SubProblem & _subproblem;
151 :
152 : /// The subdomain for the current element
153 : SubdomainID _subdomain;
154 :
155 : /// The subdomain for the last element
156 : SubdomainID _old_subdomain;
157 :
158 : /// The subdomain for the current neighbor
159 : SubdomainID _neighbor_subdomain;
160 :
161 : /// The subdomain for the last neighbor
162 : SubdomainID _old_neighbor_subdomain;
163 :
164 : /// Set to keep track of blocks for which we have printed the execution pattern
165 : mutable std::set<std::pair<const SubdomainID, const SubdomainID>> _blocks_exec_printed;
166 :
167 : /// Set to keep track of boundaries for which we have printed the execution pattern
168 : mutable std::set<BoundaryID> _boundaries_exec_printed;
169 :
170 : /// Holds caught runtime error messages
171 : std::string _error_message;
172 :
173 : private:
174 : /// Whether this is the zeroth threaded copy of this body
175 : const bool _zeroth_copy;
176 :
177 : /// The value of \p Moose::_throw_on_error at the time of construction. This data member only has
178 : /// meaning and will only be read if this is the thread 0 copy of the class
179 : const bool _incoming_throw_on_error;
180 : };
181 :
182 : template <typename RangeType>
183 95610 : ThreadedFaceLoop<RangeType>::ThreadedFaceLoop(FEProblemBase & fe_problem,
184 : const unsigned int nl_system_num,
185 : const std::set<TagID> & tags,
186 : const bool on_displaced)
187 95610 : : _fe_problem(fe_problem),
188 191220 : _mesh(fe_problem.mesh()),
189 95610 : _tags(tags),
190 95610 : _nl_system_num(nl_system_num),
191 95610 : _on_displaced(on_displaced),
192 191220 : _subproblem(_on_displaced ? static_cast<SubProblem &>(*_fe_problem.getDisplacedProblem())
193 95610 : : static_cast<SubProblem &>(_fe_problem)),
194 95610 : _zeroth_copy(true),
195 191220 : _incoming_throw_on_error(Moose::_throw_on_error)
196 : {
197 95610 : Moose::_throw_on_error = true;
198 95610 : }
199 :
200 : template <typename RangeType>
201 8658 : ThreadedFaceLoop<RangeType>::ThreadedFaceLoop(ThreadedFaceLoop & x, Threads::split /*split*/)
202 8658 : : _fe_problem(x._fe_problem),
203 8658 : _mesh(x._mesh),
204 8658 : _tags(x._tags),
205 8658 : _nl_system_num(x._nl_system_num),
206 8658 : _on_displaced(x._on_displaced),
207 8658 : _subproblem(x._subproblem),
208 8658 : _zeroth_copy(false),
209 8658 : _incoming_throw_on_error(false)
210 : {
211 8658 : }
212 :
213 : template <typename RangeType>
214 104268 : ThreadedFaceLoop<RangeType>::~ThreadedFaceLoop()
215 : {
216 104268 : if (_zeroth_copy)
217 : {
218 95610 : Moose::_throw_on_error = _incoming_throw_on_error;
219 :
220 95610 : if (!_error_message.empty())
221 6 : mooseError(_error_message);
222 : }
223 199872 : }
224 :
225 : template <typename RangeType>
226 : void
227 8658 : ThreadedFaceLoop<RangeType>::join(const ThreadedFaceLoop<RangeType> & y)
228 : {
229 8658 : if (_error_message.empty() && !y._error_message.empty())
230 0 : _error_message = y._error_message;
231 8658 : }
232 :
233 : // TODO: ensure the vector<faceinfo> data structure needs to be built such
234 : // that for all sides on an interface between two subdomains, the elements of
235 : // the same subdomain are used consistently for all the "elem" (i.e. not
236 : // "neighbor") parameters in order to avoid jumping back and forth along the
237 : // boundary between using one or the other subdomains' FV kernels
238 : // unpredictably.
239 : template <typename RangeType>
240 : void
241 104268 : ThreadedFaceLoop<RangeType>::operator()(const RangeType & range, bool bypass_threading)
242 : {
243 : // TODO: make this query fv flux kernel specific or somehow integrate the
244 : // fv source kernels into this loop. Also this check will need to increase
245 : // in generality if/when other systems and objects besides FV stuff get
246 : // added to this loop.
247 104268 : std::vector<FVKernel *> kernels;
248 104268 : _fe_problem.theWarehouse()
249 : .query()
250 208536 : .template condition<AttribSysNum>(_nl_system_num)
251 104268 : .template condition<AttribSystem>("FVFluxKernel")
252 104268 : .template condition<AttribDisplaced>(_on_displaced)
253 104268 : .queryInto(kernels);
254 104268 : if (kernels.size() == 0)
255 1445 : return;
256 :
257 : try
258 : {
259 : try
260 : {
261 102823 : ParallelUniqueId puid;
262 102823 : _tid = bypass_threading ? 0 : puid.id;
263 :
264 102823 : pre();
265 102823 : printGeneralExecutionInformation();
266 :
267 102823 : _subdomain = Moose::INVALID_BLOCK_ID;
268 102823 : _neighbor_subdomain = Moose::INVALID_BLOCK_ID;
269 :
270 102823 : typename RangeType::const_iterator faceinfo = range.begin();
271 16808819 : for (faceinfo = range.begin(); faceinfo != range.end(); ++faceinfo)
272 : {
273 16705996 : const Elem & elem = (*faceinfo)->elem();
274 :
275 16705996 : _fe_problem.setCurrentSubdomainID(&elem, _tid);
276 :
277 16705996 : _old_subdomain = _subdomain;
278 16705996 : _subdomain = elem.subdomain_id();
279 16705996 : if (_subdomain != _old_subdomain)
280 : {
281 148941 : subdomainChanged();
282 148941 : printBlockExecutionInformation();
283 : }
284 :
285 16705996 : _old_neighbor_subdomain = _neighbor_subdomain;
286 16705996 : if (const Elem * const neighbor = (*faceinfo)->neighborPtr())
287 : {
288 15144178 : _fe_problem.setNeighborSubdomainID(neighbor, _tid);
289 15144178 : _neighbor_subdomain = neighbor->subdomain_id();
290 : }
291 : else
292 1561818 : _neighbor_subdomain = Moose::INVALID_BLOCK_ID;
293 :
294 16705996 : if (_neighbor_subdomain != _old_neighbor_subdomain)
295 : {
296 2901180 : neighborSubdomainChanged();
297 : // This is going to cause a lot more printing
298 2901180 : printBlockExecutionInformation();
299 : }
300 :
301 16705996 : onFace(**faceinfo);
302 : // Cache data now because onBoundary may clear it. E.g. there was a nasty bug for two
303 : // variable FV systems where if one variable was executing an FVFluxKernel on a boundary
304 : // while the other was executing an FVFluxBC, the FVFluxKernel data would get lost because
305 : // onBoundary would clear the residual/Jacobian data before it was cached
306 16705996 : postFace(**faceinfo);
307 :
308 16705996 : const std::set<BoundaryID> boundary_ids = (*faceinfo)->boundaryIDs();
309 18501299 : for (auto & it : boundary_ids)
310 : {
311 1795311 : printBoundaryExecutionInformation(it);
312 1795311 : onBoundary(**faceinfo, it);
313 : }
314 :
315 16705988 : postFace(**faceinfo);
316 :
317 : } // range
318 102815 : post();
319 :
320 : // Clear execution printing sets to start printing on every block and boundary again
321 102815 : resetExecutionPrinting();
322 102823 : }
323 8 : catch (MetaPhysicL::LogicError & e)
324 : {
325 0 : moose::translateMetaPhysicLError(e);
326 : }
327 16 : catch (std::exception & e)
328 : {
329 : // Continue if we find a libMesh degenerate map exception, but
330 : // just throw for any real error
331 16 : if (!strstr(e.what(), "Jacobian") && !strstr(e.what(), "singular") &&
332 8 : !strstr(e.what(), "det != 0"))
333 8 : throw;
334 :
335 0 : mooseException("We caught a libMesh degeneracy exception in ComputeFVFluxThread:\n",
336 : e.what());
337 : }
338 : }
339 8 : catch (MooseException & e)
340 : {
341 0 : caughtMooseException(e);
342 : }
343 16 : catch (std::runtime_error & e)
344 : {
345 8 : _error_message = e.what();
346 : }
347 104268 : }
348 :
349 : /**
350 : * Base class for assembly-like calculations.
351 : */
352 : template <typename RangeType, typename AttributeTagType>
353 : class ComputeFVFluxThread : public ThreadedFaceLoop<RangeType>
354 : {
355 : public:
356 : ComputeFVFluxThread(FEProblemBase & fe_problem,
357 : const unsigned int nl_system_num,
358 : const std::set<TagID> & tags,
359 : bool on_displaced);
360 :
361 : ComputeFVFluxThread(ComputeFVFluxThread & x, Threads::split split);
362 :
363 : virtual ~ComputeFVFluxThread();
364 :
365 : virtual void onFace(const FaceInfo & fi) override;
366 : virtual void pre() override;
367 : virtual void post() override;
368 :
369 : virtual void onBoundary(const FaceInfo & fi, BoundaryID boundary) override;
370 :
371 : virtual void subdomainChanged() override;
372 : virtual void neighborSubdomainChanged() override;
373 :
374 : protected:
375 : /**
376 : * call either computeResidual, computeJacobian, or computeResidualAndJacobian on the provided
377 : * residual object depending on what derived class of this class is instantiated
378 : */
379 : virtual void compute(FVFaceResidualObject & ro, const FaceInfo & fi) = 0;
380 :
381 : /**
382 : * call either residualSetup or jacobianSetup depending on what derived class of this class is
383 : * instantiated
384 : */
385 : virtual void setup(SetupInterface & obj) = 0;
386 :
387 : /**
388 : * call either addCachedJacobian or addCachedResidual or both depending on what derived class of
389 : * this class is instantiated
390 : */
391 : virtual void addCached() = 0;
392 :
393 : unsigned int _num_cached = 0;
394 :
395 : using ThreadedFaceLoop<RangeType>::_fe_problem;
396 : using ThreadedFaceLoop<RangeType>::_mesh;
397 : using ThreadedFaceLoop<RangeType>::_tid;
398 : using ThreadedFaceLoop<RangeType>::_tags;
399 : using ThreadedFaceLoop<RangeType>::_nl_system_num;
400 : using ThreadedFaceLoop<RangeType>::_on_displaced;
401 : using ThreadedFaceLoop<RangeType>::_subdomain;
402 : using ThreadedFaceLoop<RangeType>::_neighbor_subdomain;
403 : using ThreadedFaceLoop<RangeType>::_blocks_exec_printed;
404 : using ThreadedFaceLoop<RangeType>::_boundaries_exec_printed;
405 :
406 : private:
407 : void reinitVariables(const FaceInfo & fi);
408 : void finalizeContainers();
409 :
410 : /// Print list of object types executed and in which order
411 : virtual void printGeneralExecutionInformation() const override;
412 :
413 : /// Print ordering of objects executed on each block
414 : virtual void printBlockExecutionInformation() const override;
415 :
416 : /// Print ordering of objects exected on each boundary
417 : virtual void printBoundaryExecutionInformation(const BoundaryID bnd_id) const override;
418 :
419 : /// Utility to get the subdomain names from the ids
420 : std::pair<SubdomainName, SubdomainName> getBlockNames() const;
421 :
422 : /// Variables
423 : std::set<MooseVariableFieldBase *> _fv_vars;
424 : std::set<MooseVariableFieldBase *> _elem_sub_fv_vars;
425 : std::set<MooseVariableFieldBase *> _neigh_sub_fv_vars;
426 :
427 : /// FVFluxKernels
428 : std::set<FVFluxKernel *> _fv_flux_kernels;
429 : std::set<FVFluxKernel *> _elem_sub_fv_flux_kernels;
430 : std::set<FVFluxKernel *> _neigh_sub_fv_flux_kernels;
431 :
432 : /// Element face materials
433 : std::vector<std::shared_ptr<MaterialBase>> _elem_face_mats;
434 : std::vector<std::shared_ptr<MaterialBase>> _elem_sub_elem_face_mats;
435 : std::vector<std::shared_ptr<MaterialBase>> _neigh_sub_elem_face_mats;
436 :
437 : // Neighbor face materials
438 : std::vector<std::shared_ptr<MaterialBase>> _neigh_face_mats;
439 : std::vector<std::shared_ptr<MaterialBase>> _elem_sub_neigh_face_mats;
440 : std::vector<std::shared_ptr<MaterialBase>> _neigh_sub_neigh_face_mats;
441 :
442 : const bool _scaling_jacobian;
443 : const bool _scaling_residual;
444 : };
445 :
446 : template <typename RangeType, typename AttributeTagType>
447 95610 : ComputeFVFluxThread<RangeType, AttributeTagType>::ComputeFVFluxThread(FEProblemBase & fe_problem,
448 : unsigned int nl_system_num,
449 : const std::set<TagID> & tags,
450 : bool on_displaced)
451 : : ThreadedFaceLoop<RangeType>(fe_problem, nl_system_num, tags, on_displaced),
452 95610 : _scaling_jacobian(fe_problem.computingScalingJacobian()),
453 191220 : _scaling_residual(fe_problem.computingScalingResidual())
454 : {
455 95610 : }
456 :
457 : template <typename RangeType, typename AttributeTagType>
458 8658 : ComputeFVFluxThread<RangeType, AttributeTagType>::ComputeFVFluxThread(ComputeFVFluxThread & x,
459 : Threads::split split)
460 : : ThreadedFaceLoop<RangeType>(x, split),
461 8658 : _fv_vars(x._fv_vars),
462 8658 : _scaling_jacobian(x._scaling_jacobian),
463 8658 : _scaling_residual(x._scaling_residual)
464 : {
465 8658 : }
466 :
467 : template <typename RangeType, typename AttributeTagType>
468 104268 : ComputeFVFluxThread<RangeType, AttributeTagType>::~ComputeFVFluxThread()
469 : {
470 104268 : }
471 :
472 : template <typename RangeType, typename AttributeTagType>
473 : void
474 16705996 : ComputeFVFluxThread<RangeType, AttributeTagType>::reinitVariables(const FaceInfo & fi)
475 : {
476 : // TODO: this skips necessary FE reinit. In addition to this call, we need
477 : // to conditionally do some FE-specific reinit here if we have any active FE
478 : // variables. However, we still want to keep/do FV-style quadrature.
479 : // Figure out how to do all this some day.
480 16705996 : this->_subproblem.reinitFVFace(_tid, fi);
481 :
482 : // TODO: for FE variables, this is handled via setting needed vars through
483 : // fe problem API which passes the value on to the system class. Then
484 : // reinit is called on fe problem which forwards its calls to the system
485 : // function and then fe problem also calls displaced problem reinit. All
486 : // the call forwarding seems silly, but it does allow the displaced problem
487 : // to be easily kept in sync. However, the displaced problem has different
488 : // pointers for its own face info objects, etc, we can't just pass in fe
489 : // problem's face info down to the sub-problem -- centroids are different,
490 : // volumes are different, etc. To support displaced meshes correctly, we
491 : // need to be able to reinit the subproblems variables using its equivalent
492 : // face info object. How? See https://github.com/idaholab/moose/issues/15064
493 :
494 32857032 : for (auto var : _fv_vars)
495 16151036 : var->computeFaceValues(fi);
496 :
497 16705996 : _fe_problem.resizeMaterialData(Moose::MaterialDataType::FACE_MATERIAL_DATA, /*nqp=*/1, _tid);
498 :
499 16717564 : for (std::shared_ptr<MaterialBase> mat : _elem_face_mats)
500 : {
501 11568 : mat->setFaceInfo(fi);
502 11568 : mat->computeProperties();
503 : }
504 :
505 16705996 : if (fi.neighborPtr())
506 : {
507 15144178 : _fe_problem.resizeMaterialData(
508 : Moose::MaterialDataType::NEIGHBOR_MATERIAL_DATA, /*nqp=*/1, _tid);
509 :
510 15152002 : for (std::shared_ptr<MaterialBase> mat : _neigh_face_mats)
511 : {
512 7824 : mat->setFaceInfo(fi);
513 7824 : mat->computeProperties();
514 : }
515 : }
516 16705996 : }
517 :
518 : template <typename RangeType, typename AttributeTagType>
519 : void
520 16705996 : ComputeFVFluxThread<RangeType, AttributeTagType>::onFace(const FaceInfo & fi)
521 : {
522 16705996 : reinitVariables(fi);
523 :
524 35339851 : for (auto * const k : _fv_flux_kernels)
525 18633855 : compute(*k, fi);
526 16705996 : }
527 :
528 : template <typename RangeType, typename AttributeTagType>
529 : void
530 1795311 : ComputeFVFluxThread<RangeType, AttributeTagType>::onBoundary(const FaceInfo & fi, BoundaryID bnd_id)
531 : {
532 1795311 : if (_scaling_jacobian || _scaling_residual)
533 208 : return;
534 :
535 1795103 : std::vector<FVFluxBC *> bcs;
536 1795103 : _fe_problem.theWarehouse()
537 : .query()
538 3590206 : .template condition<AttribSysNum>(_nl_system_num)
539 1795103 : .template condition<AttribSystem>("FVFluxBC")
540 1795103 : .template condition<AttribDisplaced>(_on_displaced)
541 1795103 : .template condition<AttribThread>(_tid)
542 1795103 : .template condition<AttributeTagType>(_tags)
543 1795103 : .template condition<AttribBoundaries>(bnd_id)
544 1795103 : .queryInto(bcs);
545 :
546 2317860 : for (auto * const bc : bcs)
547 522765 : compute(*bc, fi);
548 :
549 1795095 : std::vector<FVInterfaceKernel *> iks;
550 1795095 : _fe_problem.theWarehouse()
551 : .query()
552 3590190 : .template condition<AttribSysNum>(_nl_system_num)
553 1795095 : .template condition<AttribSystem>("FVInterfaceKernel")
554 1795095 : .template condition<AttribDisplaced>(_on_displaced)
555 1795095 : .template condition<AttribThread>(_tid)
556 1795095 : .template condition<AttributeTagType>(_tags)
557 1795095 : .template condition<AttribBoundaries>(bnd_id)
558 1795095 : .queryInto(iks);
559 :
560 1808738 : for (auto * const ik : iks)
561 13643 : compute(*ik, fi);
562 1795103 : }
563 :
564 : template <typename RangeType, typename AttributeTagType>
565 : void
566 102815 : ComputeFVFluxThread<RangeType, AttributeTagType>::post()
567 : {
568 : // make sure we add any remaining cached residuals/jacobians to add/record
569 102815 : Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
570 102815 : addCached();
571 :
572 102815 : _fe_problem.clearActiveElementalMooseVariables(_tid);
573 102815 : _fe_problem.clearActiveMaterialProperties(_tid);
574 102815 : }
575 :
576 : template <typename RangeType, typename AttributeTagType>
577 : void
578 3050121 : ComputeFVFluxThread<RangeType, AttributeTagType>::finalizeContainers()
579 : {
580 : //
581 : // Finalize our variables
582 : //
583 6100242 : std::set_union(_elem_sub_fv_vars.begin(),
584 : _elem_sub_fv_vars.end(),
585 : _neigh_sub_fv_vars.begin(),
586 : _neigh_sub_fv_vars.end(),
587 3050121 : std::inserter(_fv_vars, _fv_vars.begin()));
588 :
589 : //
590 : // Finalize our kernels
591 : //
592 3050121 : const bool same_kernels = _elem_sub_fv_flux_kernels == _neigh_sub_fv_flux_kernels;
593 3050121 : if (same_kernels)
594 1523368 : _fv_flux_kernels = _elem_sub_fv_flux_kernels;
595 : else
596 3053506 : std::set_union(_elem_sub_fv_flux_kernels.begin(),
597 : _elem_sub_fv_flux_kernels.end(),
598 : _neigh_sub_fv_flux_kernels.begin(),
599 : _neigh_sub_fv_flux_kernels.end(),
600 1526753 : std::inserter(_fv_flux_kernels, _fv_flux_kernels.begin()));
601 3050121 : const bool need_ghosting = !same_kernels;
602 :
603 : //
604 : // Finalize our element face materials
605 : //
606 3050121 : _elem_face_mats = _elem_sub_elem_face_mats;
607 :
608 3050121 : if (need_ghosting)
609 : // Add any element face materials from the neighboring subdomain that do not exist on the
610 : // element subdomain
611 1526753 : for (std::shared_ptr<MaterialBase> neigh_sub_elem_face_mat : _neigh_sub_elem_face_mats)
612 0 : if (std::find(_elem_sub_elem_face_mats.begin(),
613 : _elem_sub_elem_face_mats.end(),
614 0 : neigh_sub_elem_face_mat) == _elem_sub_elem_face_mats.end())
615 0 : _elem_face_mats.push_back(neigh_sub_elem_face_mat);
616 :
617 : //
618 : // Finalize our neighbor face materials
619 : //
620 3050121 : _neigh_face_mats = _neigh_sub_neigh_face_mats;
621 :
622 3050121 : if (need_ghosting)
623 : // Add any neighbor face materials from the element subdomain that do not exist on the
624 : // neighbor subdomain
625 1529287 : for (std::shared_ptr<MaterialBase> elem_sub_neigh_face_mat : _elem_sub_neigh_face_mats)
626 2534 : if (std::find(_neigh_sub_neigh_face_mats.begin(),
627 : _neigh_sub_neigh_face_mats.end(),
628 5068 : elem_sub_neigh_face_mat) == _neigh_sub_neigh_face_mats.end())
629 2534 : _neigh_face_mats.push_back(elem_sub_neigh_face_mat);
630 3050121 : }
631 :
632 : template <typename RangeType, typename AttributeTagType>
633 : void
634 148941 : ComputeFVFluxThread<RangeType, AttributeTagType>::subdomainChanged()
635 : {
636 148941 : ThreadedFaceLoop<RangeType>::subdomainChanged();
637 :
638 : // Clear variables
639 148941 : _fv_vars.clear();
640 148941 : _elem_sub_fv_vars.clear();
641 :
642 : // Clear kernels
643 148941 : _fv_flux_kernels.clear();
644 148941 : _elem_sub_fv_flux_kernels.clear();
645 :
646 : // Clear element face materials
647 148941 : _elem_face_mats.clear();
648 148941 : _elem_sub_elem_face_mats.clear();
649 :
650 : // Clear neighbor face materials
651 148941 : _neigh_face_mats.clear();
652 148941 : _elem_sub_neigh_face_mats.clear();
653 :
654 : // TODO: do this for other relevant objects - like FV BCs, FV source term
655 : // kernels, etc. - but we don't need to add them for other types of objects
656 : // like FE or DG kernels because those kernels don't run in this loop. Do we
657 : // really want to integrate fv source kernels into this loop?
658 148941 : std::vector<FVFluxKernel *> kernels;
659 148941 : _fe_problem.theWarehouse()
660 : .query()
661 297882 : .template condition<AttribSysNum>(_nl_system_num)
662 148941 : .template condition<AttribSystem>("FVFluxKernel")
663 148941 : .template condition<AttribDisplaced>(_on_displaced)
664 148941 : .template condition<AttribSubdomains>(_subdomain)
665 148941 : .template condition<AttribThread>(_tid)
666 148941 : .template condition<AttributeTagType>(_tags)
667 148941 : .queryInto(kernels);
668 :
669 148941 : _elem_sub_fv_flux_kernels = std::set<FVFluxKernel *>(kernels.begin(), kernels.end());
670 :
671 293872 : for (auto * k : _elem_sub_fv_flux_kernels)
672 : {
673 : // TODO: we need a better way to do this - especially when FE objects begin to
674 : // couple to FV vars. This code shoud be refactored out into one place
675 : // where it is easy for users to say initialize all materials and
676 : // variables needed by these objects for me.
677 144931 : const auto & deps = k->getMooseVariableDependencies();
678 289862 : for (auto var : deps)
679 : {
680 : mooseAssert(var->isFV(),
681 : "We do not currently support coupling of FE variables into FV objects");
682 144931 : _elem_sub_fv_vars.insert(var);
683 : }
684 : }
685 :
686 148941 : _fe_problem.getFVMatsAndDependencies(
687 148941 : _subdomain, _elem_sub_elem_face_mats, _elem_sub_neigh_face_mats, _elem_sub_fv_vars, _tid);
688 :
689 148941 : finalizeContainers();
690 148941 : }
691 :
692 : template <typename RangeType, typename AttributeTagType>
693 : void
694 2901180 : ComputeFVFluxThread<RangeType, AttributeTagType>::neighborSubdomainChanged()
695 : {
696 2901180 : ThreadedFaceLoop<RangeType>::neighborSubdomainChanged();
697 :
698 : // Clear variables
699 2901180 : _fv_vars.clear();
700 2901180 : _neigh_sub_fv_vars.clear();
701 :
702 : // Clear kernels
703 2901180 : _fv_flux_kernels.clear();
704 2901180 : _neigh_sub_fv_flux_kernels.clear();
705 :
706 : // Clear element face materials
707 2901180 : _elem_face_mats.clear();
708 2901180 : _neigh_sub_elem_face_mats.clear();
709 :
710 : // Clear neighbor face materials
711 2901180 : _neigh_face_mats.clear();
712 2901180 : _neigh_sub_neigh_face_mats.clear();
713 :
714 : // TODO: do this for other relevant objects - like FV BCs, FV source term
715 : // kernels, etc. - but we don't need to add them for other types of objects
716 : // like FE or DG kernels because those kernels don't run in this loop. Do we
717 : // really want to integrate fv source kernels into this loop?
718 2901180 : std::vector<FVFluxKernel *> kernels;
719 2901180 : _fe_problem.theWarehouse()
720 : .query()
721 5802360 : .template condition<AttribSysNum>(_nl_system_num)
722 2901180 : .template condition<AttribSystem>("FVFluxKernel")
723 2901180 : .template condition<AttribDisplaced>(_on_displaced)
724 2901180 : .template condition<AttribSubdomains>(_neighbor_subdomain)
725 2901180 : .template condition<AttribThread>(_tid)
726 2901180 : .template condition<AttributeTagType>(_tags)
727 2901180 : .queryInto(kernels);
728 :
729 2901180 : _neigh_sub_fv_flux_kernels = std::set<FVFluxKernel *>(kernels.begin(), kernels.end());
730 :
731 4435154 : for (auto * k : _neigh_sub_fv_flux_kernels)
732 : {
733 : // TODO: we need a better way to do this - especially when FE objects begin to
734 : // couple to FV vars. This code shoud be refactored out into one place
735 : // where it is easy for users to say initialize all materials and
736 : // variables needed by these objects for me.
737 1533974 : const auto & deps = k->getMooseVariableDependencies();
738 3067948 : for (auto var : deps)
739 : {
740 : mooseAssert(var->isFV(),
741 : "We do not currently support coupling of FE variables into FV objects");
742 1533974 : _neigh_sub_fv_vars.insert(var);
743 : }
744 : }
745 :
746 2901180 : _fe_problem.getFVMatsAndDependencies(_neighbor_subdomain,
747 2901180 : _neigh_sub_elem_face_mats,
748 2901180 : _neigh_sub_neigh_face_mats,
749 2901180 : _neigh_sub_fv_vars,
750 : _tid);
751 :
752 2901180 : finalizeContainers();
753 2901180 : }
754 :
755 : template <typename RangeType, typename AttributeTagType>
756 : void
757 102823 : ComputeFVFluxThread<RangeType, AttributeTagType>::pre()
758 : {
759 102823 : std::vector<FVFluxBC *> bcs;
760 102823 : _fe_problem.theWarehouse()
761 : .query()
762 205646 : .template condition<AttribSysNum>(_nl_system_num)
763 102823 : .template condition<AttribSystem>("FVFluxBC")
764 102823 : .template condition<AttribDisplaced>(_on_displaced)
765 102823 : .template condition<AttribThread>(_tid)
766 102823 : .template condition<AttributeTagType>(_tags)
767 102823 : .queryInto(bcs);
768 :
769 102823 : std::vector<FVInterfaceKernel *> iks;
770 102823 : _fe_problem.theWarehouse()
771 : .query()
772 205646 : .template condition<AttribSysNum>(_nl_system_num)
773 102823 : .template condition<AttribSystem>("FVInterfaceKernel")
774 102823 : .template condition<AttribDisplaced>(_on_displaced)
775 102823 : .template condition<AttribThread>(_tid)
776 102823 : .template condition<AttributeTagType>(_tags)
777 102823 : .queryInto(iks);
778 :
779 102823 : std::vector<FVFluxKernel *> kernels;
780 102823 : _fe_problem.theWarehouse()
781 : .query()
782 205646 : .template condition<AttribSysNum>(_nl_system_num)
783 102823 : .template condition<AttribSystem>("FVFluxKernel")
784 102823 : .template condition<AttribDisplaced>(_on_displaced)
785 102823 : .template condition<AttribThread>(_tid)
786 102823 : .template condition<AttributeTagType>(_tags)
787 102823 : .queryInto(kernels);
788 :
789 163172 : for (auto * const bc : bcs)
790 60349 : setup(*bc);
791 118306 : for (auto * const ik : iks)
792 15483 : setup(*ik);
793 230871 : for (auto * const kernel : kernels)
794 128048 : setup(*kernel);
795 :
796 : // Clear variables
797 102823 : _fv_vars.clear();
798 102823 : _elem_sub_fv_vars.clear();
799 102823 : _neigh_sub_fv_vars.clear();
800 :
801 : // Clear kernels
802 102823 : _fv_flux_kernels.clear();
803 102823 : _elem_sub_fv_flux_kernels.clear();
804 102823 : _neigh_sub_fv_flux_kernels.clear();
805 :
806 : // Clear element face materials
807 102823 : _elem_face_mats.clear();
808 102823 : _elem_sub_elem_face_mats.clear();
809 102823 : _neigh_sub_elem_face_mats.clear();
810 :
811 : // Clear neighbor face materials
812 102823 : _neigh_face_mats.clear();
813 102823 : _elem_sub_neigh_face_mats.clear();
814 102823 : _neigh_sub_neigh_face_mats.clear();
815 102823 : }
816 :
817 : template <typename RangeType>
818 : class ComputeFVFluxResidualThread : public ComputeFVFluxThread<RangeType, AttribVectorTags>
819 : {
820 : public:
821 : ComputeFVFluxResidualThread(FEProblemBase & fe_problem,
822 : const unsigned int nl_system_num,
823 : const std::set<TagID> & tags,
824 : bool on_displaced);
825 :
826 : ComputeFVFluxResidualThread(ComputeFVFluxResidualThread & x, Threads::split split);
827 :
828 : protected:
829 : using ComputeFVFluxThread<RangeType, AttribVectorTags>::_fe_problem;
830 : using ComputeFVFluxThread<RangeType, AttribVectorTags>::_tid;
831 : using ComputeFVFluxThread<RangeType, AttribVectorTags>::_nl_system_num;
832 : using ComputeFVFluxThread<RangeType, AttribVectorTags>::_on_displaced;
833 : using ComputeFVFluxThread<RangeType, AttribVectorTags>::_num_cached;
834 :
835 : void postFace(const FaceInfo & fi) override;
836 14125693 : void compute(FVFaceResidualObject & ro, const FaceInfo & fi) override { ro.computeResidual(fi); }
837 148519 : void setup(SetupInterface & obj) override { obj.residualSetup(); }
838 75799 : void addCached() override { this->_subproblem.SubProblem::addCachedResidual(_tid); }
839 : };
840 :
841 : template <typename RangeType>
842 70281 : ComputeFVFluxResidualThread<RangeType>::ComputeFVFluxResidualThread(
843 : FEProblemBase & fe_problem,
844 : const unsigned int nl_system_num,
845 : const std::set<TagID> & tags,
846 : bool on_displaced)
847 70281 : : ComputeFVFluxThread<RangeType, AttribVectorTags>(fe_problem, nl_system_num, tags, on_displaced)
848 : {
849 70281 : }
850 :
851 : template <typename RangeType>
852 6395 : ComputeFVFluxResidualThread<RangeType>::ComputeFVFluxResidualThread(ComputeFVFluxResidualThread & x,
853 : Threads::split split)
854 6395 : : ComputeFVFluxThread<RangeType, AttribVectorTags>(x, split)
855 : {
856 6395 : }
857 :
858 : template <typename RangeType>
859 : void
860 24717650 : ComputeFVFluxResidualThread<RangeType>::postFace(const FaceInfo & /*fi*/)
861 : {
862 24717650 : _num_cached++;
863 : // TODO: do we need both calls - or just the neighbor one? - confirm this
864 24717650 : this->_subproblem.SubProblem::cacheResidual(_tid);
865 24717650 : this->_subproblem.SubProblem::cacheResidualNeighbor(_tid);
866 :
867 24717650 : this->_subproblem.SubProblem::addCachedResidual(_tid);
868 24717650 : }
869 :
870 : template <typename RangeType>
871 : class ComputeFVFluxJacobianThread : public ComputeFVFluxThread<RangeType, AttribMatrixTags>
872 : {
873 : public:
874 : ComputeFVFluxJacobianThread(FEProblemBase & fe_problem,
875 : const unsigned int nl_system_num,
876 : const std::set<TagID> & tags,
877 : bool on_displaced);
878 :
879 : ComputeFVFluxJacobianThread(ComputeFVFluxJacobianThread & x, Threads::split split);
880 :
881 : protected:
882 : using ComputeFVFluxThread<RangeType, AttribMatrixTags>::_fe_problem;
883 : using ComputeFVFluxThread<RangeType, AttribMatrixTags>::_tid;
884 : using ComputeFVFluxThread<RangeType, AttribMatrixTags>::_num_cached;
885 :
886 : void postFace(const FaceInfo & fi) override;
887 4886886 : void compute(FVFaceResidualObject & ro, const FaceInfo & fi) override { ro.computeJacobian(fi); }
888 54120 : void setup(SetupInterface & obj) override { obj.jacobianSetup(); }
889 25959 : void addCached() override { this->_subproblem.SubProblem::addCachedJacobian(_tid); }
890 : };
891 :
892 : template <typename RangeType>
893 24025 : ComputeFVFluxJacobianThread<RangeType>::ComputeFVFluxJacobianThread(
894 : FEProblemBase & fe_problem,
895 : const unsigned int nl_system_num,
896 : const std::set<TagID> & tags,
897 : bool on_displaced)
898 24025 : : ComputeFVFluxThread<RangeType, AttribMatrixTags>(fe_problem, nl_system_num, tags, on_displaced)
899 : {
900 24025 : }
901 :
902 : template <typename RangeType>
903 2120 : ComputeFVFluxJacobianThread<RangeType>::ComputeFVFluxJacobianThread(ComputeFVFluxJacobianThread & x,
904 : Threads::split split)
905 2120 : : ComputeFVFluxThread<RangeType, AttribMatrixTags>(x, split)
906 : {
907 2120 : }
908 :
909 : template <typename RangeType>
910 : void
911 8378966 : ComputeFVFluxJacobianThread<RangeType>::postFace(const FaceInfo & /*fi*/)
912 : {
913 8378966 : _num_cached++;
914 : // FV objects do not store their Jacobian data in TaggingInterface data structures; instead they
915 : // add into the cache directly. So we do not need the same cache calls for the Jacobian that we
916 : // need for the residual
917 8378966 : if (_num_cached % 20 == 0)
918 : {
919 409724 : Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
920 409724 : this->_subproblem.SubProblem::addCachedJacobian(_tid);
921 409724 : }
922 8378966 : }
923 :
924 : template <typename RangeType>
925 : class ComputeFVFluxRJThread : public ComputeFVFluxThread<RangeType, AttribVectorTags>
926 : {
927 : public:
928 : ComputeFVFluxRJThread(FEProblemBase & fe_problem,
929 : const unsigned int nl_system_num,
930 : const std::set<TagID> & vector_tags,
931 : const std::set<TagID> & /*matrix_tags*/,
932 : bool on_displaced);
933 :
934 : ComputeFVFluxRJThread(ComputeFVFluxRJThread & x, Threads::split split);
935 :
936 : protected:
937 : using ComputeFVFluxThread<RangeType, AttribVectorTags>::_fe_problem;
938 : using ComputeFVFluxThread<RangeType, AttribVectorTags>::_tid;
939 : using ComputeFVFluxThread<RangeType, AttribVectorTags>::_num_cached;
940 :
941 : void postFace(const FaceInfo & fi) override;
942 157684 : void compute(FVFaceResidualObject & ro, const FaceInfo & fi) override
943 : {
944 157684 : ro.computeResidualAndJacobian(fi);
945 157684 : }
946 1241 : void setup(SetupInterface & obj) override { obj.residualSetup(); }
947 1057 : void addCached() override
948 : {
949 1057 : this->_subproblem.SubProblem::addCachedResidual(_tid);
950 1057 : this->_subproblem.SubProblem::addCachedJacobian(_tid);
951 1057 : }
952 : };
953 :
954 : template <typename RangeType>
955 1304 : ComputeFVFluxRJThread<RangeType>::ComputeFVFluxRJThread(FEProblemBase & fe_problem,
956 : const unsigned int nl_system_num,
957 : const std::set<TagID> & vector_tags,
958 : const std::set<TagID> & /*matrix_tags*/,
959 : bool on_displaced)
960 : : ComputeFVFluxThread<RangeType, AttribVectorTags>(
961 1304 : fe_problem, nl_system_num, vector_tags, on_displaced)
962 : {
963 1304 : }
964 :
965 : template <typename RangeType>
966 143 : ComputeFVFluxRJThread<RangeType>::ComputeFVFluxRJThread(ComputeFVFluxRJThread & x,
967 : Threads::split split)
968 143 : : ComputeFVFluxThread<RangeType, AttribVectorTags>(x, split)
969 : {
970 143 : }
971 :
972 : template <typename RangeType>
973 : void
974 315368 : ComputeFVFluxRJThread<RangeType>::postFace(const FaceInfo & /*fi*/)
975 : {
976 315368 : _num_cached++;
977 : // TODO: do we need both calls - or just the neighbor one? - confirm this
978 315368 : this->_subproblem.SubProblem::cacheResidual(_tid);
979 315368 : this->_subproblem.SubProblem::cacheResidualNeighbor(_tid);
980 : // FV objects do not store their Jacobian data in TaggingInterface data structures; instead they
981 : // add into the cache directly. So we do not need the same cache calls for the Jacobian that we
982 : // need for the residual
983 315368 : if (_num_cached % 20 == 0)
984 : {
985 15511 : Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
986 15511 : this->_subproblem.SubProblem::addCachedResidual(_tid);
987 15511 : this->_subproblem.SubProblem::addCachedJacobian(_tid);
988 15511 : }
989 315368 : }
990 :
991 : template <typename RangeType, typename AttributeTagType>
992 : void
993 102823 : ComputeFVFluxThread<RangeType, AttributeTagType>::printGeneralExecutionInformation() const
994 : {
995 102823 : if (!_fe_problem.shouldPrintExecution(_tid))
996 102791 : return;
997 32 : auto & console = _fe_problem.console();
998 32 : auto execute_on = _fe_problem.getCurrentExecuteOnFlag();
999 32 : console << "[DBG] Beginning finite volume flux objects loop on " << execute_on << std::endl;
1000 32 : mooseDoOnce(console << "[DBG] Loop on faces (FaceInfo), objects ordered on each face: "
1001 : << std::endl;
1002 : console << "[DBG] - (finite volume) flux kernels" << std::endl;
1003 : console << "[DBG] - (finite volume) flux boundary conditions" << std::endl;
1004 : console << "[DBG] - (finite volume) interface kernels" << std::endl;);
1005 32 : }
1006 :
1007 : template <typename RangeType, typename AttributeTagType>
1008 : void
1009 3050121 : ComputeFVFluxThread<RangeType, AttributeTagType>::printBlockExecutionInformation() const
1010 : {
1011 3050121 : if (!_fe_problem.shouldPrintExecution(_tid) || !_fv_flux_kernels.size())
1012 3050025 : return;
1013 :
1014 : // Print the location of the execution
1015 96 : const auto block_pair = std::make_pair(_subdomain, _neighbor_subdomain);
1016 96 : const auto block_pair_names = this->getBlockNames();
1017 96 : if (_blocks_exec_printed.count(block_pair))
1018 0 : return;
1019 96 : auto & console = _fe_problem.console();
1020 96 : console << "[DBG] Flux kernels on block " << block_pair_names.first;
1021 96 : if (_neighbor_subdomain != Moose::INVALID_BLOCK_ID)
1022 64 : console << " and neighbor " << block_pair_names.second << std::endl;
1023 : else
1024 32 : console << " with no neighbor block" << std::endl;
1025 :
1026 : // Print the list of objects
1027 96 : std::vector<MooseObject *> fv_flux_kernels;
1028 288 : for (const auto & fv_kernel : _fv_flux_kernels)
1029 192 : fv_flux_kernels.push_back(dynamic_cast<MooseObject *>(fv_kernel));
1030 384 : console << ConsoleUtils::formatString(ConsoleUtils::mooseObjectVectorToString(fv_flux_kernels),
1031 : "[DBG]")
1032 96 : << std::endl;
1033 96 : _blocks_exec_printed.insert(block_pair);
1034 96 : }
1035 :
1036 : template <typename RangeType, typename AttributeTagType>
1037 : void
1038 1795311 : ComputeFVFluxThread<RangeType, AttributeTagType>::printBoundaryExecutionInformation(
1039 : const BoundaryID bnd_id) const
1040 : {
1041 1795311 : if (!_fe_problem.shouldPrintExecution(_tid))
1042 1795223 : return;
1043 88 : if (_boundaries_exec_printed.count(bnd_id))
1044 0 : return;
1045 88 : std::vector<MooseObject *> bcs;
1046 88 : _fe_problem.theWarehouse()
1047 : .query()
1048 176 : .template condition<AttribSystem>("FVFluxBC")
1049 88 : .template condition<AttribDisplaced>(_on_displaced)
1050 88 : .template condition<AttribThread>(_tid)
1051 88 : .template condition<AttributeTagType>(_tags)
1052 88 : .template condition<AttribBoundaries>(bnd_id)
1053 88 : .queryInto(bcs);
1054 :
1055 88 : std::vector<MooseObject *> iks;
1056 88 : _fe_problem.theWarehouse()
1057 : .query()
1058 176 : .template condition<AttribSystem>("FVInterfaceKernel")
1059 88 : .template condition<AttribDisplaced>(_on_displaced)
1060 88 : .template condition<AttribThread>(_tid)
1061 88 : .template condition<AttributeTagType>(_tags)
1062 88 : .template condition<AttribBoundaries>(bnd_id)
1063 88 : .queryInto(iks);
1064 :
1065 88 : const auto block_pair_names = this->getBlockNames();
1066 88 : if (bcs.size())
1067 : {
1068 64 : auto & console = _fe_problem.console();
1069 64 : console << "[DBG] FVBCs on boundary " << bnd_id << " between subdomain "
1070 64 : << block_pair_names.first;
1071 64 : if (_neighbor_subdomain != Moose::INVALID_BLOCK_ID)
1072 32 : console << " and neighbor " << block_pair_names.second << std::endl;
1073 : else
1074 32 : console << " and the exterior of the mesh " << std::endl;
1075 128 : const std::string fv_bcs = ConsoleUtils::mooseObjectVectorToString(bcs);
1076 64 : console << ConsoleUtils::formatString(fv_bcs, "[DBG]") << std::endl;
1077 64 : }
1078 88 : if (iks.size())
1079 : {
1080 32 : auto & console = _fe_problem.console();
1081 32 : console << "[DBG] FVIKs on boundary " << bnd_id << " between subdomain "
1082 32 : << block_pair_names.first;
1083 32 : if (_neighbor_subdomain != Moose::INVALID_BLOCK_ID)
1084 32 : console << " and neighbor " << block_pair_names.second << std::endl;
1085 : else
1086 0 : console << " and the exterior of the mesh " << std::endl;
1087 64 : const std::string fv_iks = ConsoleUtils::mooseObjectVectorToString(iks);
1088 32 : console << ConsoleUtils::formatString(fv_iks, "[DBG]") << std::endl;
1089 32 : }
1090 88 : _boundaries_exec_printed.insert(bnd_id);
1091 88 : }
1092 :
1093 : template <typename RangeType, typename AttributeTagType>
1094 : std::pair<SubdomainName, SubdomainName>
1095 184 : ComputeFVFluxThread<RangeType, AttributeTagType>::getBlockNames() const
1096 : {
1097 184 : auto block_names = std::make_pair(_mesh.getSubdomainName(_subdomain),
1098 184 : _mesh.getSubdomainName(_neighbor_subdomain));
1099 184 : if (block_names.first == "")
1100 184 : block_names.first = Moose::stringify(_subdomain);
1101 184 : if (block_names.second == "")
1102 184 : block_names.second = Moose::stringify(_neighbor_subdomain);
1103 368 : return block_names;
1104 184 : }
|