https://mooseframework.inl.gov
NonlinearThread.C
Go to the documentation of this file.
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 "NonlinearThread.h"
11 #include "NonlinearSystem.h"
12 #include "Problem.h"
13 #include "FEProblem.h"
14 #include "KernelBase.h"
15 #include "IntegratedBCBase.h"
16 #include "DGKernelBase.h"
17 #include "InterfaceKernelBase.h"
18 #include "Material.h"
19 #include "TimeKernel.h"
20 #include "SwapBackSentinel.h"
21 #include "FVTimeKernel.h"
22 #include "ComputeJacobianThread.h"
23 
24 #include "libmesh/threads.h"
25 
27  : ThreadedElementLoop<ConstElemRange>(fe_problem),
28  _nl(fe_problem.currentNonlinearSystem()),
29  _num_cached(0),
30  _integrated_bcs(_nl.getIntegratedBCWarehouse()),
31  _dg_kernels(_nl.getDGKernelWarehouse()),
32  _interface_kernels(_nl.getInterfaceKernelWarehouse()),
33  _kernels(_nl.getKernelWarehouse()),
34  _hdg_kernels(_nl.getHDGKernelWarehouse()),
35  _has_active_objects(_integrated_bcs.hasActiveObjects() || _dg_kernels.hasActiveObjects() ||
36  _interface_kernels.hasActiveObjects() || _kernels.hasActiveObjects() ||
37  _fe_problem.haveFV()),
38  _should_execute_dg(false)
39 {
40 }
41 
42 // Splitting Constructor
45  _nl(x._nl),
46  _num_cached(x._num_cached),
47  _integrated_bcs(x._integrated_bcs),
48  _dg_kernels(x._dg_kernels),
49  _interface_kernels(x._interface_kernels),
50  _kernels(x._kernels),
51  _tag_kernels(x._tag_kernels),
52  _hdg_kernels(x._hdg_kernels),
53  _has_active_objects(x._has_active_objects),
54  _should_execute_dg(x._should_execute_dg)
55 {
56 }
57 
59 
60 void
61 NonlinearThread::operator()(const ConstElemRange & range, bool bypass_threading)
62 {
64  ThreadedElementLoop<ConstElemRange>::operator()(range, bypass_threading);
65 }
66 
67 void
69 {
70  // This should come first to setup the residual objects before we do dependency determination of
71  // material properties and variables
73 
75 
76  // Update variable Dependencies
77  std::set<MooseVariableFEBase *> needed_moose_vars;
82 
83  // Update FE variable coupleable vector tags
84  std::set<TagID> needed_fe_var_vector_tags;
86  _subdomain, needed_fe_var_vector_tags, _tid);
88  _subdomain, needed_fe_var_vector_tags, _tid);
90  _subdomain, needed_fe_var_vector_tags, _tid);
91 
92  // Update material dependencies
93  std::unordered_set<unsigned int> needed_mat_props;
98 
99  if (_fe_problem.haveFV())
100  {
101  // Re-query the finite volume elemental kernels
102  _fv_kernels.clear();
104  .query()
105  .template condition<AttribSysNum>(_nl.number())
106  .template condition<AttribSystem>("FVElementalKernel")
107  .template condition<AttribSubdomains>(_subdomain)
108  .template condition<AttribThread>(_tid)
110  for (const auto fv_kernel : _fv_kernels)
111  {
112  const auto & fv_mv_deps = fv_kernel->getMooseVariableDependencies();
113  needed_moose_vars.insert(fv_mv_deps.begin(), fv_mv_deps.end());
114  const auto & fv_mp_deps = fv_kernel->getMatPropDependencies();
115  needed_mat_props.insert(fv_mp_deps.begin(), fv_mp_deps.end());
116  }
117  }
118 
120  _fe_problem.setActiveFEVariableCoupleableVectorTags(needed_fe_var_vector_tags, _tid);
121  _fe_problem.prepareMaterials(needed_mat_props, _subdomain, _tid);
122 }
123 
124 void
125 NonlinearThread::onElement(const Elem * const elem)
126 {
127  // Set up Sentinel class so that, even if reinitMaterials() throws in prepareElement, we
128  // still remember to swap back during stack unwinding.
130 
131  prepareElement(elem);
132 
133  if (dynamic_cast<ComputeJacobianThread *>(this))
134  if (_nl.getScalarVariables(_tid).size() > 0)
136 
138 }
139 
140 void
142 {
144  {
145  const auto & kernels = _tag_kernels->getActiveBlockObjects(_subdomain, _tid);
146  for (const auto & kernel : kernels)
147  compute(*kernel);
148  }
149 
150  if (_fe_problem.haveFV())
151  for (auto kernel : _fv_kernels)
152  compute(*kernel);
153 }
154 
155 void
156 NonlinearThread::onBoundary(const Elem * const elem,
157  const unsigned int side,
158  const BoundaryID bnd_id,
159  const Elem * const lower_d_elem /*=nullptr*/)
160 {
162  {
163  // Set up Sentinel class so that, after we swap in reinitMaterialsFace in prepareFace, even if
164  // one of our callees throws we remember to swap back during stack unwinding. We put our
165  // sentinel here as opposed to in prepareFace because we certainly don't want our materials
166  // swapped back before we proceed to residual/Jacobian computation
168 
169  prepareFace(elem, side, bnd_id, lower_d_elem);
170  computeOnBoundary(bnd_id, lower_d_elem);
171 
172  if (lower_d_elem)
173  {
174  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
175  accumulateLower();
176  }
177  }
178 }
179 
180 void
181 NonlinearThread::computeOnBoundary(BoundaryID bnd_id, const Elem * /*lower_d_elem*/)
182 {
183  const auto & bcs = _ibc_warehouse->getActiveBoundaryObjects(bnd_id, _tid);
184  for (const auto & bc : bcs)
185  if (bc->shouldApply())
186  compute(*bc);
187 }
188 
189 void
190 NonlinearThread::onInterface(const Elem * elem, unsigned int side, BoundaryID bnd_id)
191 {
193  {
194 
195  // Pointer to the neighbor we are currently working on.
196  const Elem * neighbor = elem->neighbor_ptr(side);
197 
198  if (neighbor->active())
199  {
200  _fe_problem.reinitNeighbor(elem, side, _tid);
201 
202  // Set up Sentinels so that, even if one of the reinitMaterialsXXX() calls throws, we
203  // still remember to swap back during stack unwinding. Note that face, boundary, and interface
204  // all operate with the same MaterialData object
206  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
208 
210  _fe_problem.reinitMaterialsNeighbor(neighbor->subdomain_id(), _tid);
211 
212  // Has to happen after face and neighbor properties have been computed. Note that we don't use
213  // a sentinel here because FEProblem::swapBackMaterialsFace is going to handle face materials,
214  // boundary materials, and interface materials (e.g. it queries the boundary material data
215  // with the current element and side
217 
218  computeOnInterface(bnd_id);
219 
220  {
221  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
223  }
224  }
225  }
226 }
227 
228 void
230 {
231  const auto & int_ks = _ik_warehouse->getActiveBoundaryObjects(bnd_id, _tid);
232  for (const auto & interface_kernel : int_ks)
233  compute(*interface_kernel);
234 }
235 
236 void
237 NonlinearThread::onInternalSide(const Elem * elem, unsigned int side)
238 {
240  {
241  // Pointer to the neighbor we are currently working on.
242  const Elem * neighbor = elem->neighbor_ptr(side);
243 
245 
246  // Set up Sentinels so that, even if one of the reinitMaterialsXXX() calls throws, we
247  // still remember to swap back during stack unwinding.
249  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
250 
252  _fe_problem.reinitMaterialsNeighbor(neighbor->subdomain_id(), _tid);
253 
254  computeOnInternalFace(neighbor);
255 
256  {
257  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
259  }
260  }
262  {
263  // Set up Sentinel class so that, after we swap in reinitMaterialsFace in prepareFace, even if
264  // one of our callees throws we remember to swap back during stack unwinding. We put our
265  // sentinel here as opposed to in prepareFace because we certainly don't want our materials
266  // swapped back before we proceed to residual/Jacobian computation
268 
269  prepareFace(elem, side, Moose::INVALID_BOUNDARY_ID, nullptr);
271  }
272 }
273 
274 void
276 {
277  const auto & dgks = _dg_warehouse->getActiveBlockObjects(_subdomain, _tid);
278  for (const auto & dg_kernel : dgks)
279  if (dg_kernel->hasBlocks(neighbor->subdomain_id()))
280  compute(*dg_kernel, neighbor);
281 }
282 
283 void
285 {
286  compute(static_cast<ResidualObject &>(kernel));
287 }
288 
289 void
291 {
292  compute(static_cast<ResidualObject &>(kernel));
293 }
294 
295 void
297 {
298  compute(static_cast<ResidualObject &>(bc));
299 }
300 
301 void
302 NonlinearThread::compute(DGKernelBase & dg, const Elem * /*neighbor*/)
303 {
304  compute(static_cast<ResidualObject &>(dg));
305 }
306 
307 void
309 {
310  compute(static_cast<ResidualObject &>(ik));
311 }
312 
313 void
314 NonlinearThread::postElement(const Elem * /*elem*/)
315 {
316  accumulate();
317 }
318 
319 void
321 {
323 }
324 
325 void
327 {
329  {
330  const auto & console = _fe_problem.console();
331  const auto execute_on = _fe_problem.getCurrentExecuteOnFlag();
332  console << "[DBG] Beginning elemental loop to compute " + objectType() + " on " << execute_on
333  << std::endl;
334  mooseDoOnce(
335  console << "[DBG] Execution order on each element:" << std::endl;
336  console << "[DBG] - kernels on element quadrature points" << std::endl;
337  console << "[DBG] - finite volume elemental kernels on element" << std::endl;
338  console << "[DBG] - integrated boundary conditions on element side quadrature points"
339  << std::endl;
340  console << "[DBG] - DG kernels on element side quadrature points" << std::endl;
341  console << "[DBG] - interface kernels on element side quadrature points" << std::endl;);
342  }
343 }
344 
345 void
347 {
348  // Number of objects executing is approximated by size of warehouses
349  const int num_objects = _kernels.size() + _fv_kernels.size() + _integrated_bcs.size() +
351  const auto & console = _fe_problem.console();
352  const auto block_name = _mesh.getSubdomainName(_subdomain);
353 
354  if (_fe_problem.shouldPrintExecution(_tid) && num_objects > 0)
355  {
356  if (_blocks_exec_printed.count(_subdomain))
357  return;
358  console << "[DBG] Ordering of " + objectType() + " Objects on block " << block_name << " ("
359  << _subdomain << ")" << std::endl;
361  {
362  console << "[DBG] Ordering of kernels:" << std::endl;
363  console << _kernels.activeObjectsToFormattedString() << std::endl;
364  }
365  if (_fv_kernels.size())
366  {
367  console << "[DBG] Ordering of FV elemental kernels:" << std::endl;
368  std::string fvkernels =
369  std::accumulate(_fv_kernels.begin() + 1,
370  _fv_kernels.end(),
371  _fv_kernels[0]->name(),
372  [](const std::string & str_out, FVElementalKernel * kernel)
373  { return str_out + " " + kernel->name(); });
374  console << ConsoleUtils::formatString(fvkernels, "[DBG]") << std::endl;
375  }
377  {
378  console << "[DBG] Ordering of DG kernels:" << std::endl;
379  console << _dg_kernels.activeObjectsToFormattedString() << std::endl;
380  }
381  }
382  else if (_fe_problem.shouldPrintExecution(_tid) && num_objects == 0 &&
384  console << "[DBG] No Active " + objectType() + " Objects on block " << block_name << " ("
385  << _subdomain << ")" << std::endl;
386 
388 }
389 
390 void
392 {
396  return;
397 
398  const auto & console = _fe_problem.console();
399  const auto b_name = _mesh.getBoundaryName(bid);
400  console << "[DBG] Ordering of " + objectType() + " Objects on boundary " << b_name << " (" << bid
401  << ")" << std::endl;
402 
404  {
405  console << "[DBG] Ordering of integrated boundary conditions:" << std::endl;
406  console << _integrated_bcs.activeObjectsToFormattedString() << std::endl;
407  }
408 
409  // We have not checked if we have a neighbor. This could be premature for saying we are executing
410  // interface kernels. However, we should assume the execution will happen on another side of the
411  // same boundary
413  {
414  console << "[DBG] Ordering of interface kernels:" << std::endl;
415  console << _interface_kernels.activeObjectsToFormattedString() << std::endl;
416  }
417 
418  _boundaries_exec_printed.insert(bid);
419 }
420 
421 void
422 NonlinearThread::prepareFace(const Elem * const elem,
423  const unsigned int side,
424  const BoundaryID bnd_id,
425  const Elem * const lower_d_elem)
426 {
427  _fe_problem.reinitElemFace(elem, side, _tid);
428 
429  // Needed to use lower-dimensional variables on Materials
430  if (lower_d_elem)
431  _fe_problem.reinitLowerDElem(lower_d_elem, _tid);
432 
433  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
434  if (bnd_id != Moose::INVALID_BOUNDARY_ID)
436 }
437 
438 bool
439 NonlinearThread::shouldComputeInternalSide(const Elem & elem, const Elem & neighbor) const
440 {
444 }
std::string activeObjectsToFormattedString(THREAD_ID tid=0, const std::string &prefix="[DBG]") const
Output the active content of the warehouse to a string, meant to be output to the console...
unsigned int size(THREAD_ID tid=0) const
Return how many kernels we store in the current warehouse.
virtual bool shouldComputeInternalSide(const Elem &elem, const Elem &neighbor) const
Whether to compute the internal side for the provided element-neighbor pair.
Base class for assembly-like calculations.
virtual void accumulate()=0
Add element residual/Jacobian into assembly global data.
virtual void accumulateLower()=0
Add lower-d residual/Jacobian into assembly global data.
const std::vector< MooseVariableScalar * > & getScalarVariables(THREAD_ID tid)
Definition: SystemBase.h:757
bool hasActiveBlockObjects(THREAD_ID tid=0) const
virtual void compute(ResidualObject &ro)=0
Will dispatch to computeResidual/computeJacobian/computeResidualAndJacobian based on the derived clas...
const std::map< SubdomainID, std::vector< std::shared_ptr< T > > > & getActiveBlockObjects(THREAD_ID tid=0) const
MooseObjectWarehouse< InterfaceKernelBase > * _ik_warehouse
virtual void computeOnBoundary(BoundaryID bnd_id, const Elem *lower_d_elem)
virtual bool haveFV() const override
returns true if this problem includes/needs finite volume functionality.
const BoundaryID INVALID_BOUNDARY_ID
Definition: MooseTypes.C:22
MooseObjectTagWarehouse< IntegratedBCBase > & _integrated_bcs
Reference to BC storage structures.
const ExecFlagType & getCurrentExecuteOnFlag() const
Return/set the current execution flag.
std::vector< T * > & queryInto(std::vector< T *> &results, Args &&... args)
queryInto executes the query and stores the results in the given vector.
Definition: TheWarehouse.h:311
const std::string & getBoundaryName(BoundaryID boundary_id)
Return the name of the boundary given the id.
Definition: MooseMesh.C:1787
void reinitMaterialsBoundary(BoundaryID boundary_id, const THREAD_ID tid, bool swap_stateful=true, const std::deque< MaterialBase *> *reinit_mats=nullptr)
reinit materials on a boundary
const MaterialWarehouse & getMaterialWarehouse() const
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
void prepareFace(const Elem *elem, unsigned int side, BoundaryID bnd_id=Moose::INVALID_BOUNDARY_ID, const Elem *lower_d_elem=nullptr)
Reinitialize variables and materials on a face.
void updateBlockVariableDependency(SubdomainID id, std::set< MooseVariableFieldBase *> &needed_moose_vars, THREAD_ID tid=0) const
const std::string & getSubdomainName(SubdomainID subdomain_id) const
Return the name of a block given an id.
Definition: MooseMesh.C:1758
void reinitMaterialsFace(SubdomainID blk_id, const THREAD_ID tid, bool swap_stateful=true, const std::deque< MaterialBase *> *reinit_mats=nullptr)
reinit materials on element faces
bool shouldComputeInternalSide(const Elem &elem, const Elem &neighbor) const override
Whether to compute the internal side for the provided element-neighbor pair.
Serves as a base class for DGKernel and ADDGKernel.
Definition: DGKernelBase.h:32
virtual void reinitElemNeighborAndLowerD(const Elem *elem, unsigned int side, const THREAD_ID tid) override
virtual void onBoundary(const Elem *elem, unsigned int side, BoundaryID bnd_id, const Elem *lower_d_elem=nullptr) override
Called when doing boundary assembling.
virtual void accumulateNeighbor()=0
Add neighbor residual/Jacobian into assembly global data.
bool hasActiveBoundaryObjects(THREAD_ID tid=0) const
virtual void swapBackMaterialsFace(const THREAD_ID tid)
FVElemental is used for calculating residual contributions from volume integral terms of a PDE where ...
virtual void setActiveElementalMooseVariables(const std::set< MooseVariableFEBase *> &moose_vars, const THREAD_ID tid) override
Set the MOOSE variables to be reinited on each element.
const bool _has_active_objects
Whether there are any active residual objects; otherwise we will do an early return.
bool shouldPrintExecution(const THREAD_ID tid) const
Check whether the problem should output execution orders at this time.
MooseObjectWarehouse< HDGKernel > * _hdg_warehouse
This is the common base class for the three main kernel types implemented in MOOSE, Kernel, VectorKernel and ArrayKernel.
Definition: KernelBase.h:23
TheWarehouse & theWarehouse() const
void printBlockExecutionInformation() const override
Print list of specific objects executed on each block and in which order.
void reinitMaterialsNeighbor(SubdomainID blk_id, const THREAD_ID tid, bool swap_stateful=true, const std::deque< MaterialBase *> *reinit_mats=nullptr)
reinit materials on the neighboring element face
MooseObjectTagWarehouse< KernelBase > & _kernels
std::string formatString(std::string message, const std::string &prefix)
Add new lines and prefixes to a string for pretty display in output NOTE: This makes a copy of the st...
Definition: ConsoleUtils.C:582
boundary_id_type BoundaryID
virtual void computeOnElement()
virtual std::string objectType() const
Return what the loops is meant to compute.
NonlinearThread(FEProblemBase &fe_problem)
void printGeneralExecutionInformation() const override
Print information about the loop, mostly order of execution of objects.
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1159
virtual void computeOnInterface(BoundaryID bnd_id)
NonlinearSystemBase & _nl
virtual void postElement(const Elem *) override
Called after the element assembly is done (including surface assembling)
virtual void operator()(const RangeType &range, bool bypass_threading=false)
virtual void operator()(const ConstElemRange &range, bool bypass_threading=false) override
void updateBoundaryMatPropDependency(std::unordered_set< unsigned int > &needed_mat_props, THREAD_ID tid=0) const
const std::map< BoundaryID, std::vector< std::shared_ptr< T > > > & getActiveBoundaryObjects(THREAD_ID tid=0) const
virtual void swapBackMaterialsNeighbor(const THREAD_ID tid)
virtual void onInternalSide(const Elem *elem, unsigned int side) override
Called when doing internal edge assembling.
virtual void computeOnInternalFace()=0
std::set< SubdomainID > _blocks_exec_printed
Keep track of which blocks were visited.
tbb::split split
MooseObjectWarehouse< IntegratedBCBase > * _ibc_warehouse
virtual void reinitLowerDElem(const Elem *lower_d_elem, const THREAD_ID tid, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
const ConsoleStream & console() const
Return console handle.
Definition: Problem.h:48
virtual void swapBackMaterials(const THREAD_ID tid)
void updateBlockMatPropDependency(SubdomainID id, std::unordered_set< unsigned int > &needed_mat_props, THREAD_ID tid=0) const
void updateBlockFEVariableCoupledVectorTagDependency(SubdomainID id, std::set< TagID > &needed_fe_var_vector_tags, THREAD_ID tid=0) const
Update FE variable coupleable vector tag vector.
std::set< BoundaryID > _boundaries_exec_printed
Keep track of which boundaries were visited.
virtual void subdomainSetup(SubdomainID subdomain, const THREAD_ID tid)
bool _should_execute_dg
Whether DG kernels should be executed for a given elem-neighbor pairing.
Query query()
query creates and returns an initialized a query object for querying objects from the warehouse...
Definition: TheWarehouse.h:466
virtual void accumulateNeighborLower()=0
Add neighbor and lower residual/Jacobian into assembly global data.
MooseObjectTagWarehouse< InterfaceKernelBase > & _interface_kernels
Reference to interface kernel storage structure.
Base class for deriving any boundary condition of a integrated type.
void printBoundaryExecutionInformation(const unsigned int bid) const override
Print list of specific objects executed on each boundary and in which order.
virtual void onInterface(const Elem *elem, unsigned int side, BoundaryID bnd_id) override
Called when doing interface assembling.
void reinitElemFace(const Elem *elem, unsigned int side, BoundaryID, const THREAD_ID tid)
InterfaceKernelBase is the base class for all InterfaceKernel type classes.
virtual void determineObjectWarehouses()=0
Determine the objects we will actually compute based on vector/matrix tag information.
SubdomainID _subdomain
The subdomain for the current element.
MooseObjectWarehouse< KernelBase > * _tag_kernels
MooseObjectWarehouse< DGKernelBase > * _dg_warehouse
virtual void reinitNeighbor(const Elem *elem, unsigned int side, const THREAD_ID tid) override
virtual void setActiveFEVariableCoupleableVectorTags(std::set< TagID > &vtags, const THREAD_ID tid) override
virtual void onElement(const Elem *elem) override
Assembly of the element (not including surface assembly)
void reinitMaterialsInterface(BoundaryID boundary_id, const THREAD_ID tid, bool swap_stateful=true)
void updateBoundaryVariableDependency(std::set< MooseVariableFieldBase *> &needed_moose_vars, THREAD_ID tid=0) const
std::vector< FVElementalKernel * > _fv_kernels
Current subdomain FVElementalKernels.
void prepareMaterials(const std::unordered_set< unsigned int > &consumer_needed_mat_props, const SubdomainID blk_id, const THREAD_ID tid)
Add the MooseVariables and the material properties that the current materials depend on to the depend...
virtual void reinitOffDiagScalars(const THREAD_ID tid) override
virtual void post() override
Called after the element range loop.
MooseObjectTagWarehouse< DGKernelBase > & _dg_kernels
Reference to DGKernel storage structure.
virtual ~NonlinearThread()
virtual void subdomainChanged() override
Called every time the current subdomain changes (i.e.
The "SwapBackSentinel" class&#39;s destructor guarantees that FEProblemBase::swapBackMaterials{Face,Neighbor}() is called even when an exception is thrown from FEProblemBase::reinitMaterials{Face,Neighbor}.