https://mooseframework.inl.gov
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ComputeElemAuxVarsThread< AuxKernelType > Class Template Reference

#include <ComputeElemAuxVarsThread.h>

Inheritance diagram for ComputeElemAuxVarsThread< AuxKernelType >:
[legend]

Public Member Functions

 ComputeElemAuxVarsThread (FEProblemBase &problem, const MooseObjectWarehouse< AuxKernelType > &storage, bool need_materials)
 
 ComputeElemAuxVarsThread (ComputeElemAuxVarsThread &x, Threads::split split)
 
virtual ~ComputeElemAuxVarsThread ()
 
virtual void subdomainChanged () override
 Called every time the current subdomain changes (i.e. More...
 
virtual void onElement (const Elem *elem) override
 Assembly of the element (not including surface assembly) More...
 
virtual void post () override
 Called after the element range loop. More...
 
void join (const ComputeElemAuxVarsThread &)
 
virtual void caughtMooseException (MooseException &e) override
 Called if a MooseException is caught anywhere during the computation. More...
 
virtual bool keepGoing () override
 Whether or not the loop should continue. More...
 
virtual void preElement (const Elem *elem) override
 Called before the element assembly. More...
 
virtual void preInternalSide (const Elem *elem, unsigned int side) override
 Called before evaluations on an element internal side. More...
 
virtual void preBoundary (const Elem *elem, unsigned int side, BoundaryID bnd_id, const Elem *lower_d_elem=nullptr) override
 Called before the boundary assembly. More...
 
virtual void neighborSubdomainChanged () override
 Called every time the neighbor subdomain changes (i.e. More...
 
virtual void operator() (const ConstElemRange &range, bool bypass_threading=false)
 
virtual void pre ()
 Called before the element range loop. More...
 
virtual void postElement (const Elem *elem)
 Called after the element assembly is done (including surface assembling) More...
 
virtual void onBoundary (const Elem *elem, unsigned int side, BoundaryID bnd_id, const Elem *lower_d_elem=nullptr)
 Called when doing boundary assembling. More...
 
virtual void postInternalSide (const Elem *elem, unsigned int side)
 Called after evaluations on an element internal side. More...
 
virtual void onInternalSide (const Elem *elem, unsigned int side)
 Called when doing internal edge assembling. More...
 
virtual void onExternalSide (const Elem *elem, unsigned int side)
 Called when iterating over external sides (no side neighbor) More...
 
virtual void onInterface (const Elem *elem, unsigned int side, BoundaryID bnd_id)
 Called when doing interface assembling. More...
 

Protected Member Functions

void printGeneralExecutionInformation () const override
 Print list of object types executed and in which order. More...
 
void printBlockExecutionInformation () const override
 Print list of specific objects executed and in which order. More...
 
void prepareElement (const Elem *elem)
 
void clearVarsAndMaterials ()
 
void printExecutionOrdering (const std::vector< T * > &objs, const bool print_header=true, const std::string &line_prefix="[DBG]") const
 Routine to output the ordering of objects within a vector of pointers to these objects. More...
 
void printExecutionOrdering (const std::vector< std::shared_ptr< T >> &objs_ptrs, const bool print_header=true, const std::string &line_prefix="[DBG]") const
 
virtual void printBoundaryExecutionInformation (const unsigned int) const
 Print information about the particular ordering of objects on each boundary. More...
 
void resetExecPrintedSets () const
 Resets the set of blocks and boundaries visited. More...
 
virtual bool shouldComputeInternalSide (const Elem &elem, const Elem &neighbor) const
 Whether to compute the internal side for the provided element-neighbor pair. More...
 

Protected Attributes

AuxiliarySystem_aux_sys
 
const MooseObjectWarehouse< AuxKernelType > & _aux_kernels
 Storage object containing active AuxKernel objects. More...
 
bool _need_materials
 
FEProblemBase_fe_problem
 
MooseMesh_mesh
 
THREAD_ID _tid
 
SubdomainID _subdomain
 The subdomain for the current element. More...
 
SubdomainID _old_subdomain
 The subdomain for the last element. More...
 
SubdomainID _neighbor_subdomain
 The subdomain for the current neighbor. More...
 
SubdomainID _old_neighbor_subdomain
 The subdomain for the last neighbor. More...
 
std::set< SubdomainID_blocks_exec_printed
 Keep track of which blocks were visited. More...
 
std::set< BoundaryID_boundaries_exec_printed
 Keep track of which boundaries were visited. More...
 

Detailed Description

template<typename AuxKernelType>
class ComputeElemAuxVarsThread< AuxKernelType >

Definition at line 23 of file ComputeElemAuxVarsThread.h.

Constructor & Destructor Documentation

◆ ComputeElemAuxVarsThread() [1/2]

template<typename AuxKernelType >
ComputeElemAuxVarsThread< AuxKernelType >::ComputeElemAuxVarsThread ( FEProblemBase problem,
const MooseObjectWarehouse< AuxKernelType > &  storage,
bool  need_materials 
)

Definition at line 23 of file ComputeElemAuxVarsThread.C.

28  _aux_sys(problem.getAuxiliarySystem()),
29  _aux_kernels(storage),
30  _need_materials(need_materials)
31 {
32 }
const MooseObjectWarehouse< AuxKernelType > & _aux_kernels
Storage object containing active AuxKernel objects.
AuxiliarySystem & getAuxiliarySystem()

◆ ComputeElemAuxVarsThread() [2/2]

template<typename AuxKernelType >
ComputeElemAuxVarsThread< AuxKernelType >::ComputeElemAuxVarsThread ( ComputeElemAuxVarsThread< AuxKernelType > &  x,
Threads::split  split 
)

Definition at line 36 of file ComputeElemAuxVarsThread.C.

39  _aux_sys(x._aux_sys),
42 {
43 }
const MooseObjectWarehouse< AuxKernelType > & _aux_kernels
Storage object containing active AuxKernel objects.
FEProblemBase & _fe_problem

◆ ~ComputeElemAuxVarsThread()

template<typename AuxKernelType >
ComputeElemAuxVarsThread< AuxKernelType >::~ComputeElemAuxVarsThread ( )
virtual

Definition at line 46 of file ComputeElemAuxVarsThread.C.

47 {
48 }

Member Function Documentation

◆ caughtMooseException()

void ThreadedElementLoop< ConstElemRange >::caughtMooseException ( MooseException e)
overridevirtualinherited

Called if a MooseException is caught anywhere during the computation.

The single input parameter taken is a MooseException object.

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 105 of file ThreadedElementLoop.h.

Referenced by ComputeJacobianForScalingThread::operator()().

106 {
107  Threads::spin_mutex::scoped_lock lock(threaded_element_mutex);
108 
109  std::string what(e.what());
111 }
virtual const char * what() const
Get out the error message.
virtual void setException(const std::string &message)
Set an exception, which is stored at this point by toggling a member variable in this class...
static Threads::spin_mutex threaded_element_mutex
This mutex is used by all derived classes of the ThreadedElementLoop.

◆ clearVarsAndMaterials()

void ThreadedElementLoop< ConstElemRange >::clearVarsAndMaterials ( )
protectedinherited

Definition at line 195 of file ThreadedElementLoop.h.

Referenced by NonlinearThread::post().

196 {
199 }
void clearActiveMaterialProperties(const THREAD_ID tid)
Clear the active material properties.
virtual void clearActiveElementalMooseVariables(const THREAD_ID tid) override
Clear the active elemental MooseVariableFEBase.

◆ join()

template<typename AuxKernelType >
void ComputeElemAuxVarsThread< AuxKernelType >::join ( const ComputeElemAuxVarsThread< AuxKernelType > &  )

Definition at line 137 of file ComputeElemAuxVarsThread.C.

138 {
139 }

◆ keepGoing()

virtual bool ThreadedElementLoop< ConstElemRange >::keepGoing ( )
inlineoverridevirtualinherited

Whether or not the loop should continue.

Returns
true to keep going, false to stop.

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 45 of file ThreadedElementLoop.h.

Referenced by ComputeJacobianForScalingThread::operator()().

45 { return !_fe_problem.hasException(); }
virtual bool hasException()
Whether or not an exception has occurred.

◆ neighborSubdomainChanged()

void ThreadedElementLoop< ConstElemRange >::neighborSubdomainChanged ( )
overridevirtualinherited

Called every time the neighbor subdomain changes (i.e.

the subdomain of this neighbor is not the same as the subdomain of the last neighbor). Beware of over-using this! You might think that you can do some expensive stuff in here and get away with it... but there are applications that have TONS of subdomains....

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 139 of file ThreadedElementLoop.h.

140 {
143 }
Base class for assembly-like calculations.
virtual void neighborSubdomainSetup(SubdomainID subdomain, const THREAD_ID tid)

◆ onBoundary()

void ThreadedElementLoopBase< ConstElemRange >::onBoundary ( const Elem *  elem,
unsigned int  side,
BoundaryID  bnd_id,
const Elem *  lower_d_elem = nullptr 
)
virtualinherited

Called when doing boundary assembling.

Parameters
elem- The element we are checking is on the boundary.
side- The side of the element in question.
bnd_id- ID of the boundary we are at
lower_d_elem- Lower dimensional element (e.g. Mortar)

Reimplemented in ComputeUserObjectsThread, NonlinearThread, ComputeIndicatorThread, ComputeMaterialsObjectThread, and ComputeMarkerThread.

Definition at line 366 of file ThreadedElementLoopBase.h.

370 {
371 }

◆ onElement()

template<typename AuxKernelType >
void ComputeElemAuxVarsThread< AuxKernelType >::onElement ( const Elem *  elem)
overridevirtual

Assembly of the element (not including surface assembly)

Parameters
elem- active element

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 91 of file ComputeElemAuxVarsThread.C.

92 {
94  {
95  const std::vector<std::shared_ptr<AuxKernelType>> & kernels =
97  _fe_problem.prepare(elem, _tid);
99 
100  // Set up the sentinel so that, even if reinitMaterials() throws, we
101  // still remember to swap back.
103 
104  if (_need_materials)
105  _fe_problem.reinitMaterials(elem->subdomain_id(), _tid);
106 
107  for (const auto & aux : kernels)
108  {
109  aux->compute();
110  aux->variable().insert(_aux_sys.solution());
111 
112  // update the aux solution vector if writable coupled variables are used
113  if (aux->hasWritableCoupledVariables())
114  {
115  for (auto * var : aux->getWritableCoupledVariables())
116  var->insert(_aux_sys.solution());
117 
118  _fe_problem.reinitElem(elem, _tid);
119  }
120  }
121  }
122 }
bool hasActiveBlockObjects(THREAD_ID tid=0) const
const MooseObjectWarehouse< AuxKernelType > & _aux_kernels
Storage object containing active AuxKernel objects.
const std::map< SubdomainID, std::vector< std::shared_ptr< T > > > & getActiveBlockObjects(THREAD_ID tid=0) const
virtual void prepare(const Elem *elem, const THREAD_ID tid) override
NumericVector< Number > & solution()
Definition: SystemBase.h:195
virtual void reinitElem(const Elem *elem, const THREAD_ID tid) override
virtual void swapBackMaterials(const THREAD_ID tid)
void reinitMaterials(SubdomainID blk_id, const THREAD_ID tid, bool swap_stateful=true)
SubdomainID _subdomain
The subdomain for the current element.
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}.

◆ onExternalSide()

void ThreadedElementLoopBase< ConstElemRange >::onExternalSide ( const Elem *  elem,
unsigned int  side 
)
virtualinherited

Called when iterating over external sides (no side neighbor)

Parameters
elem- Element we are on
side- local side number of the element 'elem'

Reimplemented in ComputeUserObjectsThread.

Definition at line 393 of file ThreadedElementLoopBase.h.

394 {
395 }

◆ onInterface()

void ThreadedElementLoopBase< ConstElemRange >::onInterface ( const Elem *  elem,
unsigned int  side,
BoundaryID  bnd_id 
)
virtualinherited

Called when doing interface assembling.

Parameters
elem- Element we are on
side- local side number of the element 'elem'
bnd_id- ID of the interface we are at

Reimplemented in ComputeUserObjectsThread, NonlinearThread, and ComputeMaterialsObjectThread.

Definition at line 399 of file ThreadedElementLoopBase.h.

402 {
403 }

◆ onInternalSide()

void ThreadedElementLoopBase< ConstElemRange >::onInternalSide ( const Elem *  elem,
unsigned int  side 
)
virtualinherited

Called when doing internal edge assembling.

Parameters
elem- Element we are on
side- local side number of the element 'elem'

Reimplemented in ComputeUserObjectsThread, NonlinearThread, ComputeIndicatorThread, ComputeMaterialsObjectThread, and ComputeMarkerThread.

Definition at line 387 of file ThreadedElementLoopBase.h.

388 {
389 }

◆ operator()()

void ThreadedElementLoopBase< ConstElemRange >::operator() ( const ConstElemRange range,
bool  bypass_threading = false 
)
virtualinherited

Reimplemented in NonlinearThread, and ComputeJacobianForScalingThread.

Definition at line 223 of file ThreadedElementLoopBase.h.

224 {
225  try
226  {
227  try
228  {
229  ParallelUniqueId puid;
230  _tid = bypass_threading ? 0 : puid.id;
231 
232  pre();
234 
237  typename RangeType::const_iterator el = range.begin();
238  for (el = range.begin(); el != range.end(); ++el)
239  {
240  if (!keepGoing())
241  break;
242 
243  const Elem * elem = *el;
244 
245  preElement(elem);
246 
248  _subdomain = elem->subdomain_id();
249  if (_subdomain != _old_subdomain)
250  {
253  }
254 
255  onElement(elem);
256 
257  if (_mesh.interiorLowerDBlocks().count(elem->subdomain_id()) > 0 ||
258  _mesh.boundaryLowerDBlocks().count(elem->subdomain_id()) > 0)
259  {
260  postElement(elem);
261  continue;
262  }
263 
264  for (unsigned int side = 0; side < elem->n_sides(); side++)
265  {
266  std::vector<BoundaryID> boundary_ids = _mesh.getBoundaryIDs(elem, side);
267  const Elem * lower_d_elem = _mesh.getLowerDElem(elem, side);
268 
269  if (boundary_ids.size() > 0)
270  for (std::vector<BoundaryID>::iterator it = boundary_ids.begin();
271  it != boundary_ids.end();
272  ++it)
273  {
274  preBoundary(elem, side, *it, lower_d_elem);
276  onBoundary(elem, side, *it, lower_d_elem);
277  }
278 
279  const Elem * neighbor = elem->neighbor_ptr(side);
280  if (neighbor)
281  {
282  preInternalSide(elem, side);
283 
285  _neighbor_subdomain = neighbor->subdomain_id();
288 
289  if (shouldComputeInternalSide(*elem, *neighbor))
290  onInternalSide(elem, side);
291 
292  if (boundary_ids.size() > 0)
293  for (std::vector<BoundaryID>::iterator it = boundary_ids.begin();
294  it != boundary_ids.end();
295  ++it)
296  onInterface(elem, side, *it);
297 
298  postInternalSide(elem, side);
299  }
300  else
301  onExternalSide(elem, side);
302  } // sides
303 
304  postElement(elem);
305  } // range
306 
307  post();
309  }
310  catch (libMesh::LogicError & e)
311  {
312  mooseException("We caught a libMesh error in ThreadedElementLoopBase:", e.what());
313  }
314  catch (MetaPhysicL::LogicError & e)
315  {
317  }
318  }
319  catch (MooseException & e)
320  {
322  }
323 }
virtual void onExternalSide(const Elem *elem, unsigned int side)
Called when iterating over external sides (no side neighbor)
virtual bool keepGoing()
Whether or not the loop should continue.
void resetExecPrintedSets() const
Resets the set of blocks and boundaries visited.
virtual bool shouldComputeInternalSide(const Elem &elem, const Elem &neighbor) const
Whether to compute the internal side for the provided element-neighbor pair.
virtual void onElement(const Elem *elem)
Assembly of the element (not including surface assembly)
const std::set< SubdomainID > & interiorLowerDBlocks() const
Definition: MooseMesh.h:1403
void translateMetaPhysicLError(const MetaPhysicL::LogicError &)
emit a relatively clear error message when we catch a MetaPhysicL logic error
Definition: MooseError.C:112
virtual void printBoundaryExecutionInformation(const unsigned int) const
Print information about the particular ordering of objects on each boundary.
virtual void pre()
Called before the element range loop.
const Elem * getLowerDElem(const Elem *, unsigned short int) const
Returns a const pointer to a lower dimensional element that corresponds to a side of a higher dimensi...
Definition: MooseMesh.C:1698
virtual void subdomainChanged()
Called every time the current subdomain changes (i.e.
virtual void neighborSubdomainChanged()
Called every time the neighbor subdomain changes (i.e.
virtual void preInternalSide(const Elem *elem, unsigned int side)
Called before evaluations on an element internal side.
virtual void postInternalSide(const Elem *elem, unsigned int side)
Called after evaluations on an element internal side.
const SubdomainID INVALID_BLOCK_ID
Definition: MooseTypes.C:20
virtual void onBoundary(const Elem *elem, unsigned int side, BoundaryID bnd_id, const Elem *lower_d_elem=nullptr)
Called when doing boundary assembling.
virtual void postElement(const Elem *elem)
Called after the element assembly is done (including surface assembling)
virtual void printGeneralExecutionInformation() const
Print information about the loop ordering.
virtual void onInterface(const Elem *elem, unsigned int side, BoundaryID bnd_id)
Called when doing interface assembling.
SubdomainID _old_neighbor_subdomain
The subdomain for the last neighbor.
const std::set< SubdomainID > & boundaryLowerDBlocks() const
Definition: MooseMesh.h:1407
virtual void onInternalSide(const Elem *elem, unsigned int side)
Called when doing internal edge assembling.
const_iterator end() const
Provides a way for users to bail out of the current solve.
virtual void caughtMooseException(MooseException &)
Called if a MooseException is caught anywhere during the computation.
const_iterator begin() const
SubdomainID _subdomain
The subdomain for the current element.
std::vector< BoundaryID > getBoundaryIDs(const Elem *const elem, const unsigned short int side) const
Returns a vector of boundary IDs for the requested element on the requested side. ...
SubdomainID _old_subdomain
The subdomain for the last element.
virtual void post()
Called after the element range loop.
virtual void printBlockExecutionInformation() const
Print information about the particular ordering of objects on each block.
virtual void preElement(const Elem *elem)
Called before the element assembly.
virtual void preBoundary(const Elem *elem, unsigned int side, BoundaryID bnd_id, const Elem *lower_d_elem=nullptr)
Called before the boundary assembly.
SubdomainID _neighbor_subdomain
The subdomain for the current neighbor.

◆ post()

template<typename AuxKernelType >
void ComputeElemAuxVarsThread< AuxKernelType >::post ( )
overridevirtual

Called after the element range loop.

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 126 of file ComputeElemAuxVarsThread.C.

127 {
130 
133 }
void clearActiveMaterialProperties(const THREAD_ID tid)
Clear the active material properties.
virtual void clearActiveFEVariableCoupleableMatrixTags(const THREAD_ID tid) override
virtual void clearActiveFEVariableCoupleableVectorTags(const THREAD_ID tid) override
virtual void clearActiveElementalMooseVariables(const THREAD_ID tid) override
Clear the active elemental MooseVariableFEBase.

◆ postElement()

void ThreadedElementLoopBase< ConstElemRange >::postElement ( const Elem *  elem)
virtualinherited

Called after the element assembly is done (including surface assembling)

Parameters
elem- active element

Reimplemented in ComputeJacobianBlocksThread, NonlinearThread, ComputeIndicatorThread, ComputeJacobianThread, and ComputeMarkerThread.

Definition at line 351 of file ThreadedElementLoopBase.h.

352 {
353 }

◆ postInternalSide()

void ThreadedElementLoopBase< ConstElemRange >::postInternalSide ( const Elem *  elem,
unsigned int  side 
)
virtualinherited

Called after evaluations on an element internal side.

Parameters
elem- Element we are on
side- local side number of the element 'elem'

Reimplemented in ComputeJacobianBlocksThread.

Definition at line 381 of file ThreadedElementLoopBase.h.

382 {
383 }

◆ pre()

void ThreadedElementLoopBase< ConstElemRange >::pre ( )
virtualinherited

Called before the element range loop.

Definition at line 327 of file ThreadedElementLoopBase.h.

Referenced by ComputeJacobianForScalingThread::operator()().

328 {
329 }

◆ preBoundary()

void ThreadedElementLoop< ConstElemRange >::preBoundary ( const Elem *  elem,
unsigned int  side,
BoundaryID  bnd_id,
const Elem *  lower_d_elem = nullptr 
)
overridevirtualinherited

Called before the boundary assembly.

Parameters
elem- The element we are checking is on the boundary.
side- The side of the element in question.
bnd_id- ID of the boundary we are at
lower_d_elem- Lower dimensional element (e.g. Mortar)

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 129 of file ThreadedElementLoop.h.

133 {
135 }
virtual void setCurrentBoundaryID(BoundaryID bid, const THREAD_ID tid) override
sets the current boundary ID in assembly
Base class for assembly-like calculations.

◆ preElement()

void ThreadedElementLoop< ConstElemRange >::preElement ( const Elem *  elem)
overridevirtualinherited

Called before the element assembly.

Parameters
elem- active element

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 115 of file ThreadedElementLoop.h.

Referenced by ComputeJacobianForScalingThread::operator()().

116 {
118 }
virtual void setCurrentSubdomainID(const Elem *elem, const THREAD_ID tid) override
Base class for assembly-like calculations.

◆ preInternalSide()

void ThreadedElementLoop< ConstElemRange >::preInternalSide ( const Elem *  elem,
unsigned int  side 
)
overridevirtualinherited

Called before evaluations on an element internal side.

Parameters
elem- Element we are on
side- local side number of the element 'elem'

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 122 of file ThreadedElementLoop.h.

123 {
125 }
Base class for assembly-like calculations.
virtual void setNeighborSubdomainID(const Elem *elem, unsigned int side, const THREAD_ID tid) override

◆ prepareElement()

void ThreadedElementLoop< ConstElemRange >::prepareElement ( const Elem *  elem)
protectedinherited

Definition at line 186 of file ThreadedElementLoop.h.

Referenced by NonlinearThread::onElement().

187 {
188  _fe_problem.prepare(elem, this->_tid);
189  _fe_problem.reinitElem(elem, this->_tid);
191 }
virtual void prepare(const Elem *elem, const THREAD_ID tid) override
virtual void reinitElem(const Elem *elem, const THREAD_ID tid) override
void reinitMaterials(SubdomainID blk_id, const THREAD_ID tid, bool swap_stateful=true)
SubdomainID _subdomain
The subdomain for the current element.

◆ printBlockExecutionInformation()

template<typename AuxKernelType >
void ComputeElemAuxVarsThread< AuxKernelType >::printBlockExecutionInformation ( ) const
overrideprotectedvirtual

Print list of specific objects executed and in which order.

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 155 of file ComputeElemAuxVarsThread.C.

156 {
159  return;
160 
161  const auto & console = _fe_problem.console();
162  const auto & kernels = _aux_kernels.getActiveBlockObjects(_subdomain, _tid);
163  console << "[DBG] Ordering of AuxKernels on block " << _subdomain << std::endl;
164  printExecutionOrdering<AuxKernelType>(kernels, false);
166 }
bool hasActiveBlockObjects(THREAD_ID tid=0) const
const MooseObjectWarehouse< AuxKernelType > & _aux_kernels
Storage object containing active AuxKernel objects.
const std::map< SubdomainID, std::vector< std::shared_ptr< T > > > & getActiveBlockObjects(THREAD_ID tid=0) const
bool shouldPrintExecution(const THREAD_ID tid) const
Check whether the problem should output execution orders at this time.
std::set< SubdomainID > _blocks_exec_printed
Keep track of which blocks were visited.
const ConsoleStream & console() const
Return console handle.
Definition: Problem.h:48
SubdomainID _subdomain
The subdomain for the current element.

◆ printBoundaryExecutionInformation()

virtual void ThreadedElementLoopBase< ConstElemRange >::printBoundaryExecutionInformation ( const unsigned int  ) const
inlineprotectedvirtualinherited

Print information about the particular ordering of objects on each boundary.

Reimplemented in NonlinearThread.

Definition at line 184 of file ThreadedElementLoopBase.h.

184 {}

◆ printExecutionOrdering() [1/2]

void ThreadedElementLoop< ConstElemRange >::printExecutionOrdering ( const std::vector< T *> &  objs,
const bool  print_header = true,
const std::string &  line_prefix = "[DBG]" 
) const
protectedinherited

Routine to output the ordering of objects within a vector of pointers to these objects.

These objects must implement the name() routine, and it must return a string or compatible type.

Template Parameters
Tthe object type
Parameters
objsthe vector with all the objects (should be pointers)
objects_typethe name of the type of objects. Defaults to the CPP object name
print_headerwhether to print a header about the timing of execution and the type of objects

Definition at line 148 of file ThreadedElementLoop.h.

151 {
152  if (!objs.size())
153  return;
154 
155  auto & console = _fe_problem.console();
156  const auto objects_type = MooseUtils::prettyCppType(objs[0]);
157  std::vector<MooseObject *> moose_objs;
158  for (auto obj_ptr : objs)
159  moose_objs.push_back(dynamic_cast<MooseObject *>(obj_ptr));
160  const auto names = ConsoleUtils::mooseObjectVectorToString(moose_objs);
161 
162  // Print string with a DBG prefix and with sufficient line breaks
163  std::string message = print_header ? "Executing " + objects_type + " on " +
165  : "";
166  message += (print_header ? "Order of execution:\n" : "") + names;
167  console << ConsoleUtils::formatString(message, line_prefix) << std::endl;
168 }
std::string mooseObjectVectorToString(const std::vector< MooseObject *> &objs, const std::string &sep=" ")
Routine to output the name of MooseObjects in a string.
Definition: ConsoleUtils.C:598
const std::string & name() const
Definition: MooseEnumItem.h:35
const ExecFlagType & getCurrentExecuteOnFlag() const
Return/set the current execution flag.
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
const ConsoleStream & console() const
Return console handle.
Definition: Problem.h:48
std::string prettyCppType(const std::string &cpp_type)
Definition: MooseUtils.C:1246

◆ printExecutionOrdering() [2/2]

void ThreadedElementLoop< ConstElemRange >::printExecutionOrdering ( const std::vector< std::shared_ptr< T >> &  objs_ptrs,
const bool  print_header = true,
const std::string &  line_prefix = "[DBG]" 
) const
protectedinherited

Definition at line 173 of file ThreadedElementLoop.h.

177 {
178  std::vector<T *> regular_ptrs;
179  for (auto shared_ptr : objs_ptrs)
180  regular_ptrs.push_back(shared_ptr.get());
181  printExecutionOrdering<T>(regular_ptrs, print_header, line_prefix);
182 }

◆ printGeneralExecutionInformation()

template<typename AuxKernelType >
void ComputeElemAuxVarsThread< AuxKernelType >::printGeneralExecutionInformation ( ) const
overrideprotectedvirtual

Print list of object types executed and in which order.

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 143 of file ComputeElemAuxVarsThread.C.

144 {
146  return;
147 
148  const auto & console = _fe_problem.console();
149  const auto & execute_on = _fe_problem.getCurrentExecuteOnFlag();
150  console << "[DBG] Executing auxiliary kernels on elements on " << execute_on << std::endl;
151 }
const MooseObjectWarehouse< AuxKernelType > & _aux_kernels
Storage object containing active AuxKernel objects.
const ExecFlagType & getCurrentExecuteOnFlag() const
Return/set the current execution flag.
bool shouldPrintExecution(const THREAD_ID tid) const
Check whether the problem should output execution orders at this time.
const ConsoleStream & console() const
Return console handle.
Definition: Problem.h:48
bool hasActiveObjects(THREAD_ID tid=0) const

◆ resetExecPrintedSets()

void ThreadedElementLoopBase< ConstElemRange >::resetExecPrintedSets ( ) const
protectedinherited

Resets the set of blocks and boundaries visited.

Definition at line 444 of file ThreadedElementLoopBase.h.

445 {
446  _blocks_exec_printed.clear();
447  _boundaries_exec_printed.clear();
448 }
std::set< SubdomainID > _blocks_exec_printed
Keep track of which blocks were visited.
std::set< BoundaryID > _boundaries_exec_printed
Keep track of which boundaries were visited.

◆ shouldComputeInternalSide()

bool ThreadedElementLoopBase< ConstElemRange >::shouldComputeInternalSide ( const Elem &  elem,
const Elem &  neighbor 
) const
protectedvirtualinherited

Whether to compute the internal side for the provided element-neighbor pair.

Typically this will return true if the element id is less than the neighbor id when the elements are equal level, or when the element is more refined than the neighbor, and then false otherwise. One type of loop where the logic will be different is when projecting stateful material properties

Reimplemented in NonlinearThread, and FlagElementsThread.

Definition at line 419 of file ThreadedElementLoopBase.h.

421 {
422  auto level = [this](const auto & elem_arg)
423  {
424  if (_mesh.doingPRefinement())
425  return elem_arg.p_level();
426  else
427  return elem_arg.level();
428  };
429  const auto elem_id = elem.id(), neighbor_id = neighbor.id();
430  const auto elem_level = level(elem), neighbor_level = level(neighbor);
431 
432  // When looping over elements and then sides, we need to make sure that we do not duplicate
433  // effort, e.g. if a face is shared by element 1 and element 2, then we do not want to do compute
434  // work both when we are visiting element 1 *and* then later when visiting element 2. Our rule is
435  // to only compute when we are visiting the element that has the lower element id when element and
436  // neighbor are of the same adaptivity level, and then if they are not of the same level, then
437  // we only compute when we are visiting the finer element
438  return (neighbor.active() && (neighbor_level == elem_level) && (elem_id < neighbor_id)) ||
439  (neighbor_level < elem_level);
440 }
void doingPRefinement(bool doing_p_refinement)
Indicate whether the kind of adaptivity we&#39;re doing is p-refinement.
Definition: MooseMesh.h:1347

◆ subdomainChanged()

template<typename AuxKernelType >
void ComputeElemAuxVarsThread< AuxKernelType >::subdomainChanged ( )
overridevirtual

Called every time the current subdomain changes (i.e.

the subdomain of this element is not the same as the subdomain of the last element). Beware of over-using this! You might think that you can do some expensive stuff in here and get away with it... but there are applications that have TONS of subdomains....

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 52 of file ComputeElemAuxVarsThread.C.

53 {
55 
56  std::set<MooseVariableFEBase *> needed_moose_vars;
57  std::unordered_set<unsigned int> needed_mat_props;
58  std::set<TagID> needed_fe_var_matrix_tags;
59  std::set<TagID> needed_fe_var_vector_tags;
60 
62  _subdomain, needed_fe_var_vector_tags, _tid);
64  {
65  const std::vector<std::shared_ptr<AuxKernelType>> & kernels =
67  for (const auto & aux : kernels)
68  {
69  aux->subdomainSetup();
70  const auto & mv_deps = aux->getMooseVariableDependencies();
71  const auto & mp_deps = aux->getMatPropDependencies();
72  needed_moose_vars.insert(mv_deps.begin(), mv_deps.end());
73  needed_mat_props.insert(mp_deps.begin(), mp_deps.end());
74 
75  auto & fe_var_coup_vtags = aux->getFEVariableCoupleableVectorTags();
76  needed_fe_var_vector_tags.insert(fe_var_coup_vtags.begin(), fe_var_coup_vtags.end());
77 
78  auto & fe_var_coup_mtags = aux->getFEVariableCoupleableMatrixTags();
79  needed_fe_var_matrix_tags.insert(fe_var_coup_mtags.begin(), fe_var_coup_mtags.end());
80  }
81  }
82 
84  _fe_problem.prepareMaterials(needed_mat_props, _subdomain, _tid);
85  _fe_problem.setActiveFEVariableCoupleableMatrixTags(needed_fe_var_matrix_tags, _tid);
86  _fe_problem.setActiveFEVariableCoupleableVectorTags(needed_fe_var_vector_tags, _tid);
87 }
bool hasActiveBlockObjects(THREAD_ID tid=0) const
const MooseObjectWarehouse< AuxKernelType > & _aux_kernels
Storage object containing active AuxKernel objects.
virtual void setActiveFEVariableCoupleableMatrixTags(std::set< TagID > &mtags, const THREAD_ID tid) override
const std::map< SubdomainID, std::vector< std::shared_ptr< T > > > & getActiveBlockObjects(THREAD_ID tid=0) const
const MaterialWarehouse & getMaterialWarehouse() const
virtual void setActiveElementalMooseVariables(const std::set< MooseVariableFEBase *> &moose_vars, const THREAD_ID tid) override
Set the MOOSE variables to be reinited on each element.
void updateBlockFEVariableCoupledVectorTagDependency(SubdomainID id, std::set< TagID > &needed_fe_var_vector_tags, THREAD_ID tid=0) const
Update FE variable coupleable vector tag vector.
virtual void subdomainSetup(SubdomainID subdomain, const THREAD_ID tid)
SubdomainID _subdomain
The subdomain for the current element.
virtual void setActiveFEVariableCoupleableVectorTags(std::set< TagID > &vtags, const THREAD_ID tid) override
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...

Member Data Documentation

◆ _aux_kernels

template<typename AuxKernelType>
const MooseObjectWarehouse<AuxKernelType>& ComputeElemAuxVarsThread< AuxKernelType >::_aux_kernels
protected

Storage object containing active AuxKernel objects.

Definition at line 50 of file ComputeElemAuxVarsThread.h.

◆ _aux_sys

template<typename AuxKernelType>
AuxiliarySystem& ComputeElemAuxVarsThread< AuxKernelType >::_aux_sys
protected

Definition at line 47 of file ComputeElemAuxVarsThread.h.

◆ _blocks_exec_printed

std::set<SubdomainID> ThreadedElementLoopBase< ConstElemRange >::_blocks_exec_printed
mutableprotectedinherited

◆ _boundaries_exec_printed

std::set<BoundaryID> ThreadedElementLoopBase< ConstElemRange >::_boundaries_exec_printed
mutableprotectedinherited

Keep track of which boundaries were visited.

Definition at line 190 of file ThreadedElementLoopBase.h.

Referenced by NonlinearThread::printBoundaryExecutionInformation().

◆ _fe_problem

FEProblemBase& ThreadedElementLoop< ConstElemRange >::_fe_problem
protectedinherited

Definition at line 62 of file ThreadedElementLoop.h.

Referenced by ComputeResidualThread::accumulate(), ComputeResidualAndJacobianThread::accumulate(), ComputeResidualThread::accumulateLower(), ComputeResidualAndJacobianThread::accumulateLower(), ComputeJacobianThread::accumulateLower(), ComputeResidualThread::accumulateNeighbor(), ComputeResidualAndJacobianThread::accumulateNeighbor(), ComputeJacobianThread::accumulateNeighbor(), ComputeResidualThread::accumulateNeighborLower(), ComputeResidualAndJacobianThread::accumulateNeighborLower(), ComputeJacobianThread::accumulateNeighborLower(), ComputeJacobianThread::compute(), ComputeFullJacobianThread::computeOnBoundary(), ComputeFullJacobianThread::computeOnElement(), NonlinearThread::computeOnElement(), ComputeFullJacobianThread::computeOnInterface(), ComputeFullJacobianThread::computeOnInternalFace(), ComputeResidualThread::determineObjectWarehouses(), ComputeJacobianThread::determineObjectWarehouses(), ComputeResidualAndJacobianThread::determineObjectWarehouses(), NonlinearThread::onBoundary(), ComputeUserObjectsThread::onBoundary(), ComputeElemDampingThread::onElement(), NonlinearThread::onElement(), ComputeUserObjectsThread::onElement(), NonlinearThread::onInterface(), ComputeUserObjectsThread::onInterface(), NonlinearThread::onInternalSide(), ComputeUserObjectsThread::onInternalSide(), ComputeUserObjectsThread::post(), ComputeJacobianThread::postElement(), ComputeJacobianBlocksThread::postElement(), ComputeJacobianBlocksThread::postInternalSide(), NonlinearThread::prepareFace(), ComputeUserObjectsThread::printBlockExecutionInformation(), NonlinearThread::printBlockExecutionInformation(), NonlinearThread::printBoundaryExecutionInformation(), ComputeElemDampingThread::printGeneralExecutionInformation(), ComputeUserObjectsThread::printGeneralExecutionInformation(), NonlinearThread::printGeneralExecutionInformation(), NonlinearThread::subdomainChanged(), and ComputeUserObjectsThread::subdomainChanged().

◆ _mesh

MooseMesh& ThreadedElementLoopBase< ConstElemRange >::_mesh
protectedinherited

◆ _need_materials

template<typename AuxKernelType>
bool ComputeElemAuxVarsThread< AuxKernelType >::_need_materials
protected

Definition at line 52 of file ComputeElemAuxVarsThread.h.

◆ _neighbor_subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_neighbor_subdomain
protectedinherited

◆ _old_neighbor_subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_old_neighbor_subdomain
protectedinherited

The subdomain for the last neighbor.

Definition at line 175 of file ThreadedElementLoopBase.h.

◆ _old_subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_old_subdomain
protectedinherited

The subdomain for the last element.

Definition at line 169 of file ThreadedElementLoopBase.h.

Referenced by ComputeJacobianForScalingThread::operator()().

◆ _subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_subdomain
protectedinherited

The subdomain for the current element.

Definition at line 166 of file ThreadedElementLoopBase.h.

Referenced by ComputeFullJacobianThread::computeOnBoundary(), ComputeFullJacobianThread::computeOnElement(), NonlinearThread::computeOnElement(), ComputeResidualThread::computeOnInternalFace(), ComputeFullJacobianThread::computeOnInternalFace(), ComputeResidualAndJacobianThread::computeOnInternalFace(), ComputeJacobianThread::computeOnInternalFace(), NonlinearThread::computeOnInternalFace(), ComputeResidualThread::determineObjectWarehouses(), ComputeResidualAndJacobianThread::determineObjectWarehouses(), ComputeMaterialsObjectThread::onBoundary(), ComputeUserObjectsThread::onBoundary(), ComputeMarkerThread::onElement(), ComputeIndicatorThread::onElement(), ComputeMaterialsObjectThread::onElement(), ComputeUserObjectsThread::onElement(), ComputeMaterialsObjectThread::onInterface(), ComputeMaterialsObjectThread::onInternalSide(), NonlinearThread::onInternalSide(), ComputeJacobianForScalingThread::operator()(), ComputeJacobianBlocksThread::postInternalSide(), ComputeMarkerThread::printBlockExecutionInformation(), ComputeIndicatorThread::printBlockExecutionInformation(), ComputeUserObjectsThread::printBlockExecutionInformation(), NonlinearThread::printBlockExecutionInformation(), ComputeUserObjectsThread::querySubdomain(), NonlinearThread::shouldComputeInternalSide(), ComputeMarkerThread::subdomainChanged(), ComputeMaterialsObjectThread::subdomainChanged(), ComputeIndicatorThread::subdomainChanged(), NonlinearThread::subdomainChanged(), and ComputeUserObjectsThread::subdomainChanged().

◆ _tid

THREAD_ID ThreadedElementLoopBase< ConstElemRange >::_tid
protectedinherited

Definition at line 163 of file ThreadedElementLoopBase.h.

Referenced by ComputeResidualThread::accumulate(), ComputeResidualAndJacobianThread::accumulate(), ComputeResidualThread::accumulateLower(), ComputeResidualAndJacobianThread::accumulateLower(), ComputeJacobianThread::accumulateLower(), ComputeResidualThread::accumulateNeighbor(), ComputeResidualAndJacobianThread::accumulateNeighbor(), ComputeJacobianThread::accumulateNeighbor(), ComputeResidualThread::accumulateNeighborLower(), ComputeResidualAndJacobianThread::accumulateNeighborLower(), ComputeJacobianThread::accumulateNeighborLower(), ComputeFullJacobianThread::computeOnBoundary(), NonlinearThread::computeOnBoundary(), ComputeFullJacobianThread::computeOnElement(), NonlinearThread::computeOnElement(), ComputeFullJacobianThread::computeOnInterface(), NonlinearThread::computeOnInterface(), ComputeResidualThread::computeOnInternalFace(), ComputeFullJacobianThread::computeOnInternalFace(), ComputeResidualAndJacobianThread::computeOnInternalFace(), ComputeJacobianThread::computeOnInternalFace(), NonlinearThread::computeOnInternalFace(), ComputeResidualThread::determineObjectWarehouses(), ComputeJacobianThread::determineObjectWarehouses(), ComputeResidualAndJacobianThread::determineObjectWarehouses(), ComputeMaterialsObjectThread::onBoundary(), NonlinearThread::onBoundary(), ComputeUserObjectsThread::onBoundary(), ComputeMarkerThread::onElement(), ComputeElemDampingThread::onElement(), ComputeIndicatorThread::onElement(), ComputeMaterialsObjectThread::onElement(), NonlinearThread::onElement(), ComputeUserObjectsThread::onElement(), ComputeMaterialsObjectThread::onInterface(), NonlinearThread::onInterface(), ComputeUserObjectsThread::onInterface(), ComputeIndicatorThread::onInternalSide(), ComputeMaterialsObjectThread::onInternalSide(), NonlinearThread::onInternalSide(), ComputeUserObjectsThread::onInternalSide(), ComputeJacobianForScalingThread::operator()(), ComputeMaterialsObjectThread::post(), ComputeMarkerThread::post(), ComputeIndicatorThread::post(), ComputeUserObjectsThread::post(), ComputeJacobianThread::postElement(), ComputeJacobianBlocksThread::postElement(), ComputeJacobianBlocksThread::postInternalSide(), NonlinearThread::prepareFace(), ComputeMarkerThread::printBlockExecutionInformation(), ComputeIndicatorThread::printBlockExecutionInformation(), ComputeUserObjectsThread::printBlockExecutionInformation(), NonlinearThread::printBlockExecutionInformation(), NonlinearThread::printBoundaryExecutionInformation(), ComputeElemDampingThread::printGeneralExecutionInformation(), ComputeMarkerThread::printGeneralExecutionInformation(), ComputeIndicatorThread::printGeneralExecutionInformation(), ComputeUserObjectsThread::printGeneralExecutionInformation(), NonlinearThread::printGeneralExecutionInformation(), ComputeUserObjectsThread::queryBoundary(), ComputeUserObjectsThread::querySubdomain(), NonlinearThread::shouldComputeInternalSide(), ComputeMarkerThread::subdomainChanged(), ComputeIndicatorThread::subdomainChanged(), ComputeMaterialsObjectThread::subdomainChanged(), NonlinearThread::subdomainChanged(), and ComputeUserObjectsThread::subdomainChanged().


The documentation for this class was generated from the following files: