www.mooseframework.org
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ComputeJacobianBlocksThread Class Reference

Specialization for filling multiple "small" preconditioning matrices simulatenously. More...

#include <ComputeJacobianBlocksThread.h>

Inheritance diagram for ComputeJacobianBlocksThread:
[legend]

Public Member Functions

 ComputeJacobianBlocksThread (FEProblemBase &fe_problem, std::vector< JacobianBlock *> &blocks, const std::set< TagID > &tags)
 
 ComputeJacobianBlocksThread (ComputeJacobianBlocksThread &x, Threads::split split)
 
virtual ~ComputeJacobianBlocksThread ()
 
void join (const ComputeJacobianThread &)
 
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 onBoundary (const Elem *elem, unsigned int side, BoundaryID bnd_id) override
 Called when doing boundary assembling. More...
 
virtual void onInternalSide (const Elem *elem, unsigned int side) override
 Called when doing internal edge assembling. More...
 
virtual void onInterface (const Elem *elem, unsigned int side, BoundaryID bnd_id) override
 Called when doing interface assembling. More...
 
virtual void post () override
 Called after the element range loop. More...
 
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 neighborSubdomainChanged () override
 Called every time the neighbor subdomain changes (i.e. More...
 
void operator() (const ConstElemRange &range, bool bypass_threading=false)
 
virtual void pre ()
 Called before the element range loop. More...
 
virtual void postInternalSide (const Elem *elem, unsigned int side)
 Called after evaluations on an element internal side. More...
 

Protected Member Functions

virtual void postElement (const Elem *elem) override
 Called after the element assembly is done (including surface assembling) More...
 
virtual void computeJacobian () override
 
virtual void computeFaceJacobian (BoundaryID bnd_id) override
 
virtual void computeInternalFaceJacobian (const Elem *neighbor) override
 
virtual void computeInternalInterFaceJacobian (BoundaryID bnd_id) override
 

Protected Attributes

std::vector< JacobianBlock * > _blocks
 
NonlinearSystemBase_nl
 
unsigned int _num_cached
 
MooseObjectTagWarehouse< IntegratedBCBase > & _integrated_bcs
 
MooseObjectWarehouse< IntegratedBCBase > * _ibc_warehouse
 
MooseObjectTagWarehouse< DGKernel > & _dg_kernels
 
MooseObjectWarehouse< DGKernel > * _dg_warehouse
 
MooseObjectTagWarehouse< InterfaceKernel > & _interface_kernels
 
MooseObjectWarehouse< InterfaceKernel > * _ik_warehouse
 
MooseObjectTagWarehouse< KernelBase > & _kernels
 
MooseObjectWarehouse< KernelBase > * _warehouse
 
const std::set< TagID > & _tags
 
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...
 

Detailed Description

Specialization for filling multiple "small" preconditioning matrices simulatenously.

Definition at line 40 of file ComputeJacobianBlocksThread.h.

Constructor & Destructor Documentation

◆ ComputeJacobianBlocksThread() [1/2]

ComputeJacobianBlocksThread::ComputeJacobianBlocksThread ( FEProblemBase fe_problem,
std::vector< JacobianBlock *> &  blocks,
const std::set< TagID > &  tags 
)

Definition at line 21 of file ComputeJacobianBlocksThread.C.

24  : ComputeFullJacobianThread(fe_problem, tags), _blocks(blocks)
25 {
26 }
ComputeFullJacobianThread(FEProblemBase &fe_problem, const std::set< TagID > &tags)
std::vector< JacobianBlock * > _blocks

◆ ComputeJacobianBlocksThread() [2/2]

ComputeJacobianBlocksThread::ComputeJacobianBlocksThread ( ComputeJacobianBlocksThread x,
Threads::split  split 
)

Definition at line 29 of file ComputeJacobianBlocksThread.C.

31  : ComputeFullJacobianThread(x, split), _blocks(x._blocks)
32 {
33 }
ComputeFullJacobianThread(FEProblemBase &fe_problem, const std::set< TagID > &tags)
static PetscErrorCode Vec x
std::vector< JacobianBlock * > _blocks

◆ ~ComputeJacobianBlocksThread()

ComputeJacobianBlocksThread::~ComputeJacobianBlocksThread ( )
virtual

Definition at line 35 of file ComputeJacobianBlocksThread.C.

35 {}

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 77 of file ThreadedElementLoop.h.

78 {
79  Threads::spin_mutex::scoped_lock lock(threaded_element_mutex);
80 
81  std::string what(e.what());
83 }
virtual const char * what() const
Get out the error message.
virtual void setException(const std::string &message)
Set an exception.
static Threads::spin_mutex threaded_element_mutex
This mutex is used by all derived classes of the ThreadedElementLoop.

◆ computeFaceJacobian()

void ComputeFullJacobianThread::computeFaceJacobian ( BoundaryID  bnd_id)
overrideprotectedvirtualinherited

done only when nonlocal integrated_bcs exist in the system

Reimplemented from ComputeJacobianThread.

Definition at line 127 of file ComputeFullJacobianThread.C.

128 {
129  std::vector<std::pair<MooseVariableFEBase *, MooseVariableFEBase *>> & ce =
131  for (const auto & it : ce)
132  {
133  MooseVariableFEBase & ivar = *(it.first);
134  MooseVariableFEBase & jvar = *(it.second);
137  {
138  // only if there are dofs for j-variable (if it is subdomain restricted var, there may not be
139  // any)
140  const auto & bcs = _ibc_warehouse->getBoundaryObjects(bnd_id, _tid);
141  for (const auto & bc : bcs)
142  if (bc->shouldApply() && bc->variable().number() == ivar.number() && bc->isImplicit())
143  {
144  bc->subProblem().prepareFaceShapes(jvar.number(), _tid);
145  bc->computeJacobianBlock(jvar);
146  }
147  }
148  }
149 
152  {
153  std::vector<std::pair<MooseVariableFEBase *, MooseVariableFEBase *>> & cne =
155  for (const auto & it : cne)
156  {
157  MooseVariableFEBase & ivariable = *(it.first);
158  MooseVariableFEBase & jvariable = *(it.second);
159 
160  unsigned int ivar = ivariable.number();
161  unsigned int jvar = jvariable.number();
162 
163  if (ivariable.activeOnSubdomain(_subdomain) && jvariable.activeOnSubdomain(_subdomain) &&
165  {
166  const std::vector<std::shared_ptr<IntegratedBCBase>> & integrated_bcs =
168  for (const auto & integrated_bc : integrated_bcs)
169  {
170  std::shared_ptr<NonlocalIntegratedBC> nonlocal_integrated_bc =
172  if (nonlocal_integrated_bc)
173  if ((integrated_bc->variable().number() == ivar) && integrated_bc->isImplicit())
174  {
175  integrated_bc->subProblem().prepareFaceShapes(jvar, _tid);
176  integrated_bc->computeNonlocalOffDiagJacobian(jvar);
177  }
178  }
179  }
180  }
181  }
182 
183  const std::vector<MooseVariableScalar *> & scalar_vars = _nl.getScalarVariables(_tid);
184  if (scalar_vars.size() > 0)
185  {
186  // go over nl-variables (non-scalar)
187  const std::vector<MooseVariableFEBase *> & vars = _nl.getVariables(_tid);
188  for (const auto & ivar : vars)
189  if (ivar->activeOnSubdomain(_subdomain) > 0 &&
191  {
192  // for each variable get the list of active kernels
193  const auto & bcs = _ibc_warehouse->getActiveBoundaryObjects(bnd_id, _tid);
194  for (const auto & bc : bcs)
195  if (bc->variable().number() == ivar->number() && bc->isImplicit())
196  {
197  // now, get the list of coupled scalar vars and compute their off-diag jacobians
198  const std::vector<MooseVariableScalar *> coupled_scalar_vars =
199  bc->getCoupledMooseScalarVars();
200 
201  // Do: dvar / dscalar_var, only want to process only nl-variables (not aux ones)
202  for (const auto & jvar : coupled_scalar_vars)
203  if (_nl.hasScalarVariable(jvar->name()))
204  bc->computeJacobianBlockScalar(jvar->number());
205  }
206  }
207  }
208 }
SubProblem & subProblem()
Get a reference to the subproblem.
const std::vector< MooseVariableScalar * > & getScalarVariables(THREAD_ID tid)
Definition: SystemBase.h:555
virtual void prepareFaceShapes(unsigned int var, THREAD_ID tid)=0
unsigned int number() const
Get variable number coming from libMesh.
NonlocalIntegratedBC is used for solving integral terms in integro-differential equations.
virtual bool activeOnSubdomain(SubdomainID subdomain) const =0
Is the variable active on the subdomain?
const std::vector< MooseVariableFEBase * > & getVariables(THREAD_ID tid)
Definition: SystemBase.h:550
std::unique_ptr< T_DEST, T_DELETER > dynamic_pointer_cast(std::unique_ptr< T_SRC, T_DELETER > &src)
std::vector< std::pair< MooseVariableFEBase *, MooseVariableFEBase * > > & couplingEntries(THREAD_ID tid)
bool hasActiveBoundaryObjects(THREAD_ID tid=0) const
virtual bool checkNonlocalCouplingRequirement()
Definition: SubProblem.h:62
const std::map< BoundaryID, std::vector< std::shared_ptr< T > > > & getBoundaryObjects(THREAD_ID tid=0) const
MooseObjectWarehouse< IntegratedBCBase > * _ibc_warehouse
NonlinearSystemBase & _nl
const std::map< BoundaryID, std::vector< std::shared_ptr< T > > > & getActiveBoundaryObjects(THREAD_ID tid=0) const
std::vector< std::pair< MooseVariableFEBase *, MooseVariableFEBase * > > & nonlocalCouplingEntries(THREAD_ID tid)
SubdomainID _subdomain
The subdomain for the current element.
virtual bool hasScalarVariable(const std::string &var_name) const
Definition: SystemBase.C:646

◆ computeInternalFaceJacobian()

void ComputeFullJacobianThread::computeInternalFaceJacobian ( const Elem *  neighbor)
overrideprotectedvirtualinherited

Reimplemented from ComputeJacobianThread.

Definition at line 211 of file ComputeFullJacobianThread.C.

212 {
214  {
215  const auto & ce = _fe_problem.couplingEntries(_tid);
216  for (const auto & it : ce)
217  {
218  const auto & dgks = _dg_warehouse->getActiveBlockObjects(_subdomain, _tid);
219  for (const auto & dg : dgks)
220  {
221  MooseVariableFEBase & ivariable = *(it.first);
222  MooseVariableFEBase & jvariable = *(it.second);
223 
224  unsigned int ivar = ivariable.number();
225  unsigned int jvar = jvariable.number();
226 
227  if (dg->variable().number() == ivar && dg->isImplicit() &&
228  dg->hasBlocks(neighbor->subdomain_id()) && jvariable.activeOnSubdomain(_subdomain))
229  {
230  dg->subProblem().prepareFaceShapes(jvar, _tid);
231  dg->subProblem().prepareNeighborShapes(jvar, _tid);
232  dg->computeOffDiagJacobian(jvar);
233  }
234  }
235  }
236  }
237 }
bool hasActiveBlockObjects(THREAD_ID tid=0) const
const std::map< SubdomainID, std::vector< std::shared_ptr< T > > > & getActiveBlockObjects(THREAD_ID tid=0) const
unsigned int number() const
Get variable number coming from libMesh.
virtual bool activeOnSubdomain(SubdomainID subdomain) const =0
Is the variable active on the subdomain?
std::vector< std::pair< MooseVariableFEBase *, MooseVariableFEBase * > > & couplingEntries(THREAD_ID tid)
MooseObjectWarehouse< DGKernel > * _dg_warehouse
SubdomainID _subdomain
The subdomain for the current element.

◆ computeInternalInterFaceJacobian()

void ComputeFullJacobianThread::computeInternalInterFaceJacobian ( BoundaryID  bnd_id)
overrideprotectedvirtualinherited

Reimplemented from ComputeJacobianThread.

Definition at line 240 of file ComputeFullJacobianThread.C.

241 {
243  {
244  const auto & ce = _fe_problem.couplingEntries(_tid);
245  for (const auto & it : ce)
246  {
247  const std::vector<std::shared_ptr<InterfaceKernel>> & int_ks =
249  for (const auto & interface_kernel : int_ks)
250  {
251  if (!interface_kernel->isImplicit())
252  continue;
253 
254  unsigned int ivar = it.first->number();
255  unsigned int jvar = it.second->number();
256 
257  interface_kernel->subProblem().prepareFaceShapes(jvar, _tid);
258  interface_kernel->subProblem().prepareNeighborShapes(jvar, _tid);
259 
260  if (interface_kernel->variable().number() == ivar)
261  interface_kernel->computeElementOffDiagJacobian(jvar);
262 
263  if (interface_kernel->neighborVariable().number() == ivar)
264  interface_kernel->computeNeighborOffDiagJacobian(jvar);
265  }
266  }
267  }
268 }
std::vector< std::pair< MooseVariableFEBase *, MooseVariableFEBase * > > & couplingEntries(THREAD_ID tid)
MooseObjectWarehouse< InterfaceKernel > * _ik_warehouse
bool hasActiveBoundaryObjects(THREAD_ID tid=0) const
const std::map< BoundaryID, std::vector< std::shared_ptr< T > > > & getActiveBoundaryObjects(THREAD_ID tid=0) const

◆ computeJacobian()

void ComputeFullJacobianThread::computeJacobian ( )
overrideprotectedvirtualinherited

done only when nonlocal kernels exist in the system

Reimplemented from ComputeJacobianThread.

Definition at line 39 of file ComputeFullJacobianThread.C.

40 {
41  std::vector<std::pair<MooseVariableFEBase *, MooseVariableFEBase *>> & ce =
43  for (const auto & it : ce)
44  {
45  MooseVariableFEBase & ivariable = *(it.first);
46  MooseVariableFEBase & jvariable = *(it.second);
47 
48  unsigned int ivar = ivariable.number();
49  unsigned int jvar = jvariable.number();
50 
51  if (ivariable.activeOnSubdomain(_subdomain) && jvariable.activeOnSubdomain(_subdomain) &&
53  {
54  // only if there are dofs for j-variable (if it is subdomain restricted var, there may not be
55  // any)
56  const auto & kernels = _warehouse->getActiveVariableBlockObjects(ivar, _subdomain, _tid);
57  for (const auto & kernel : kernels)
58  if ((kernel->variable().number() == ivar) && kernel->isImplicit())
59  {
60  kernel->subProblem().prepareShapes(jvar, _tid);
61  kernel->computeOffDiagJacobian(jvariable);
62  }
63  }
64  }
65 
68  {
69  std::vector<std::pair<MooseVariableFEBase *, MooseVariableFEBase *>> & cne =
71  for (const auto & it : cne)
72  {
73  MooseVariableFEBase & ivariable = *(it.first);
74  MooseVariableFEBase & jvariable = *(it.second);
75 
76  unsigned int ivar = ivariable.number();
77  unsigned int jvar = jvariable.number();
78 
79  if (ivariable.activeOnSubdomain(_subdomain) && jvariable.activeOnSubdomain(_subdomain) &&
81  {
82  const auto & kernels = _warehouse->getActiveVariableBlockObjects(ivar, _subdomain, _tid);
83  for (const auto & kernel : kernels)
84  {
85  std::shared_ptr<NonlocalKernel> nonlocal_kernel =
87  if (nonlocal_kernel)
88  if ((kernel->variable().number() == ivar) && kernel->isImplicit())
89  {
90  kernel->subProblem().prepareShapes(jvar, _tid);
91  kernel->computeNonlocalOffDiagJacobian(jvar);
92  }
93  }
94  }
95  }
96  }
97 
98  const std::vector<MooseVariableScalar *> & scalar_vars = _nl.getScalarVariables(_tid);
99  if (scalar_vars.size() > 0)
100  {
101  // go over nl-variables (non-scalar)
102  const std::vector<MooseVariableFEBase *> & vars = _nl.getVariables(_tid);
103  for (const auto & ivariable : vars)
104  if (ivariable->activeOnSubdomain(_subdomain) > 0 &&
106  {
107  // for each variable get the list of active kernels
108  const auto & kernels =
110  for (const auto & kernel : kernels)
111  if (kernel->isImplicit())
112  {
113  // now, get the list of coupled scalar vars and compute their off-diag jacobians
114  const std::vector<MooseVariableScalar *> coupled_scalar_vars =
115  kernel->getCoupledMooseScalarVars();
116 
117  // Do: dvar / dscalar_var, only want to process only nl-variables (not aux ones)
118  for (const auto & jvariable : coupled_scalar_vars)
119  if (_nl.hasScalarVariable(jvariable->name()))
120  kernel->computeOffDiagJacobianScalar(jvariable->number());
121  }
122  }
123  }
124 }
const std::vector< std::shared_ptr< T > > & getActiveVariableBlockObjects(unsigned int variable_id, SubdomainID block_id, THREAD_ID tid=0) const
const std::vector< MooseVariableScalar * > & getScalarVariables(THREAD_ID tid)
Definition: SystemBase.h:555
unsigned int number() const
Get variable number coming from libMesh.
virtual void prepareShapes(unsigned int var, THREAD_ID tid)=0
virtual bool activeOnSubdomain(SubdomainID subdomain) const =0
Is the variable active on the subdomain?
const std::vector< MooseVariableFEBase * > & getVariables(THREAD_ID tid)
Definition: SystemBase.h:550
std::unique_ptr< T_DEST, T_DELETER > dynamic_pointer_cast(std::unique_ptr< T_SRC, T_DELETER > &src)
std::vector< std::pair< MooseVariableFEBase *, MooseVariableFEBase * > > & couplingEntries(THREAD_ID tid)
bool hasActiveVariableBlockObjects(unsigned int variable_id, SubdomainID block_id, THREAD_ID tid=0) const
Methods for checking/getting variable kernels for a variable and SubdomainID.
virtual bool checkNonlocalCouplingRequirement()
Definition: SubProblem.h:62
MooseObjectWarehouse< KernelBase > * _warehouse
SubProblem & subProblem()
Returns a reference to the SubProblem for which this Kernel is active.
Definition: KernelBase.h:101
NonlinearSystemBase & _nl
std::vector< std::pair< MooseVariableFEBase *, MooseVariableFEBase * > > & nonlocalCouplingEntries(THREAD_ID tid)
const std::string & name() const
Get the variable name.
NonlocalKernel is used for solving integral terms in integro-differential equations.
SubdomainID _subdomain
The subdomain for the current element.
virtual bool hasScalarVariable(const std::string &var_name) const
Definition: SystemBase.C:646

◆ join()

void ComputeJacobianBlocksThread::join ( const ComputeJacobianThread )
inline

Definition at line 52 of file ComputeJacobianBlocksThread.h.

52 {}

◆ 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.

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 101 of file ThreadedElementLoop.h.

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

◆ onBoundary()

void ComputeJacobianThread::onBoundary ( const Elem *  elem,
unsigned int  side,
BoundaryID  bnd_id 
)
overridevirtualinherited

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

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 199 of file ComputeJacobianThread.C.

200 {
202  {
203  _fe_problem.reinitElemFace(elem, side, bnd_id, _tid);
204 
205  // Set up Sentinel class so that, even if reinitMaterials() throws, we
206  // still remember to swap back during stack unwinding.
208 
209  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
211 
212  computeFaceJacobian(bnd_id);
213  }
214 }
virtual void reinitMaterialsBoundary(BoundaryID boundary_id, THREAD_ID tid, bool swap_stateful=true)
bool hasActiveBoundaryObjects(THREAD_ID tid=0) const
MooseObjectWarehouse< IntegratedBCBase > * _ibc_warehouse
virtual void reinitElemFace(const Elem *elem, unsigned int side, BoundaryID bnd_id, THREAD_ID tid) override
virtual void reinitMaterialsFace(SubdomainID blk_id, THREAD_ID tid, bool swap_stateful=true)
virtual void computeFaceJacobian(BoundaryID bnd_id)
virtual void swapBackMaterialsFace(THREAD_ID tid)
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}.

◆ onElement()

void ComputeJacobianThread::onElement ( const Elem *  elem)
overridevirtualinherited

Assembly of the element (not including surface assembly)

Parameters
elem- active element

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 181 of file ComputeJacobianThread.C.

182 {
183  _fe_problem.prepare(elem, _tid);
184 
185  _fe_problem.reinitElem(elem, _tid);
186 
187  // Set up Sentinel class so that, even if reinitMaterials() throws, we
188  // still remember to swap back during stack unwinding.
191 
192  if (_nl.getScalarVariables(_tid).size() > 0)
194 
195  computeJacobian();
196 }
virtual void prepare(const Elem *elem, THREAD_ID tid) override
const std::vector< MooseVariableScalar * > & getScalarVariables(THREAD_ID tid)
Definition: SystemBase.h:555
virtual void reinitElem(const Elem *elem, THREAD_ID tid) override
virtual void reinitMaterials(SubdomainID blk_id, THREAD_ID tid, bool swap_stateful=true)
virtual void swapBackMaterials(THREAD_ID tid)
NonlinearSystemBase & _nl
SubdomainID _subdomain
The subdomain for the current element.
virtual void reinitOffDiagScalars(THREAD_ID tid) override
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}.

◆ onInterface()

void ComputeJacobianThread::onInterface ( const Elem *  elem,
unsigned int  side,
BoundaryID  bnd_id 
)
overridevirtualinherited

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 from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 251 of file ComputeJacobianThread.C.

252 {
254  {
255  // Pointer to the neighbor we are currently working on.
256  const Elem * neighbor = elem->neighbor_ptr(side);
257 
258  if (neighbor->active())
259  {
260  _fe_problem.reinitNeighbor(elem, side, _tid);
261 
262  // Set up Sentinels so that, even if one of the reinitMaterialsXXX() calls throws, we
263  // still remember to swap back during stack unwinding.
265  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
267 
269  _fe_problem.reinitMaterialsNeighbor(neighbor->subdomain_id(), _tid);
270 
272 
273  {
274  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
276  }
277  }
278  }
279 }
virtual void reinitMaterialsNeighbor(SubdomainID blk_id, THREAD_ID tid, bool swap_stateful=true)
virtual void addJacobianNeighbor(THREAD_ID tid) override
virtual void reinitMaterialsBoundary(BoundaryID boundary_id, THREAD_ID tid, bool swap_stateful=true)
MooseObjectWarehouse< InterfaceKernel > * _ik_warehouse
bool hasActiveBoundaryObjects(THREAD_ID tid=0) const
virtual void computeInternalInterFaceJacobian(BoundaryID bnd_id)
virtual void swapBackMaterialsNeighbor(THREAD_ID tid)
virtual void reinitMaterialsFace(SubdomainID blk_id, THREAD_ID tid, bool swap_stateful=true)
virtual void reinitNeighbor(const Elem *elem, unsigned int side, THREAD_ID tid) override
virtual void swapBackMaterialsFace(THREAD_ID tid)
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}.

◆ onInternalSide()

void ComputeJacobianThread::onInternalSide ( const Elem *  elem,
unsigned int  side 
)
overridevirtualinherited

Called when doing internal edge assembling.

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

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 217 of file ComputeJacobianThread.C.

218 {
220  {
221  // Pointer to the neighbor we are currently working on.
222  const Elem * neighbor = elem->neighbor_ptr(side);
223 
224  // Get the global id of the element and the neighbor
225  const dof_id_type elem_id = elem->id(), neighbor_id = neighbor->id();
226 
227  if ((neighbor->active() && (neighbor->level() == elem->level()) && (elem_id < neighbor_id)) ||
228  (neighbor->level() < elem->level()))
229  {
230  _fe_problem.reinitNeighbor(elem, side, _tid);
231 
232  // Set up Sentinels so that, even if one of the reinitMaterialsXXX() calls throws, we
233  // still remember to swap back during stack unwinding.
235  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
236 
238  _fe_problem.reinitMaterialsNeighbor(neighbor->subdomain_id(), _tid);
239 
240  computeInternalFaceJacobian(neighbor);
241 
242  {
243  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
245  }
246  }
247  }
248 }
bool hasActiveBlockObjects(THREAD_ID tid=0) const
virtual void reinitMaterialsNeighbor(SubdomainID blk_id, THREAD_ID tid, bool swap_stateful=true)
virtual void addJacobianNeighbor(THREAD_ID tid) override
virtual void swapBackMaterialsNeighbor(THREAD_ID tid)
virtual void reinitMaterialsFace(SubdomainID blk_id, THREAD_ID tid, bool swap_stateful=true)
MooseObjectWarehouse< DGKernel > * _dg_warehouse
virtual void reinitNeighbor(const Elem *elem, unsigned int side, THREAD_ID tid) override
virtual void computeInternalFaceJacobian(const Elem *neighbor)
SubdomainID _subdomain
The subdomain for the current element.
virtual void swapBackMaterialsFace(THREAD_ID tid)
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}.

◆ operator()()

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

Definition at line 171 of file ThreadedElementLoopBase.h.

172 {
173  try
174  {
175  ParallelUniqueId puid;
176  _tid = bypass_threading ? 0 : puid.id;
177 
178  pre();
179 
182  typename RangeType::const_iterator el = range.begin();
183  for (el = range.begin(); el != range.end(); ++el)
184  {
185  if (!keepGoing())
186  break;
187 
188  const Elem * elem = *el;
189 
190  preElement(elem);
191 
193  _subdomain = elem->subdomain_id();
194  if (_subdomain != _old_subdomain)
196 
197  onElement(elem);
198 
199  for (unsigned int side = 0; side < elem->n_sides(); side++)
200  {
201  std::vector<BoundaryID> boundary_ids = _mesh.getBoundaryIDs(elem, side);
202 
203  if (boundary_ids.size() > 0)
204  for (std::vector<BoundaryID>::iterator it = boundary_ids.begin();
205  it != boundary_ids.end();
206  ++it)
207  onBoundary(elem, side, *it);
208 
209  const Elem * neighbor = elem->neighbor_ptr(side);
210  if (neighbor != nullptr)
211  {
212  preInternalSide(elem, side);
213 
215  _neighbor_subdomain = neighbor->subdomain_id();
218 
219  onInternalSide(elem, side);
220 
221  if (boundary_ids.size() > 0)
222  for (std::vector<BoundaryID>::iterator it = boundary_ids.begin();
223  it != boundary_ids.end();
224  ++it)
225  onInterface(elem, side, *it);
226 
227  postInternalSide(elem, side);
228  }
229  } // sides
230  postElement(elem);
231 
232  } // range
233 
234  post();
235  }
236  catch (MooseException & e)
237  {
239  }
240 }
virtual bool keepGoing()
Whether or not the loop should continue.
virtual void onElement(const Elem *elem)
Assembly of the element (not including surface assembly)
virtual void pre()
Called before the element range loop.
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.h:320
virtual void postElement(const Elem *elem)
Called after the element assembly is done (including surface assembling)
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.
virtual void onInternalSide(const Elem *elem, unsigned int side)
Called when doing internal edge assembling.
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.
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. ...
Definition: MooseMesh.C:2105
SubdomainID _subdomain
The subdomain for the current element.
SubdomainID _old_subdomain
The subdomain for the last element.
virtual void post()
Called after the element range loop.
virtual void preElement(const Elem *elem)
Called before the element assembly.
SubdomainID _neighbor_subdomain
The subdomain for the current neighbor.
virtual void onBoundary(const Elem *elem, unsigned int side, BoundaryID bnd_id)
Called when doing boundary assembling.

◆ post()

void ComputeJacobianThread::post ( )
overridevirtualinherited

Called after the element range loop.

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 295 of file ComputeJacobianThread.C.

296 {
299 }
virtual void clearActiveMaterialProperties(THREAD_ID tid) override
Clear the active material properties.
virtual void clearActiveElementalMooseVariables(THREAD_ID tid) override
Clear the active elemental MooseVariableFEBase.

◆ postElement()

void ComputeJacobianBlocksThread::postElement ( const Elem *  elem)
overrideprotectedvirtual

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

Parameters
elem- active element

Reimplemented from ComputeJacobianThread.

Definition at line 38 of file ComputeJacobianBlocksThread.C.

39 {
40  std::vector<dof_id_type>
41  dof_indices; // Do this out here to avoid creating and destroying it thousands of times
42 
43  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
44 
45  for (const auto & block : _blocks)
46  {
47  const DofMap & dof_map = block->_precond_system.get_dof_map();
48  dof_map.dof_indices(elem, dof_indices);
49 
51  block->_jacobian, block->_ivar, block->_jvar, dof_map, dof_indices, _tid);
52  }
53 }
DofMap & dof_map
std::vector< JacobianBlock * > _blocks
virtual void addJacobianBlock(SparseMatrix< Number > &jacobian, unsigned int ivar, unsigned int jvar, const DofMap &dof_map, std::vector< dof_id_type > &dof_indices, THREAD_ID tid) override

◆ 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'

Definition at line 288 of file ThreadedElementLoopBase.h.

289 {
290 }

◆ pre()

void ThreadedElementLoopBase< ConstElemRange >::pre ( )
virtualinherited

Called before the element range loop.

Definition at line 244 of file ThreadedElementLoopBase.h.

245 {
246 }

◆ 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 87 of file ThreadedElementLoop.h.

88 {
90 }
virtual void setCurrentSubdomainID(const Elem *elem, 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 94 of file ThreadedElementLoop.h.

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

◆ subdomainChanged()

void ComputeJacobianThread::subdomainChanged ( )
overridevirtualinherited

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 127 of file ComputeJacobianThread.C.

128 {
130 
131  // Update variable Dependencies
132  std::set<MooseVariableFEBase *> needed_moose_vars;
137 
138  // Update material dependencies
139  std::set<unsigned int> needed_mat_props;
144 
146  _fe_problem.setActiveMaterialProperties(needed_mat_props, _tid);
148 
149  // If users pass a empty vector or a full size of vector,
150  // we take all kernels
151  if (!_tags.size() || _tags.size() == _fe_problem.numMatrixTags())
152  {
153  _warehouse = &_kernels;
157  }
158  // If we have one tag only,
159  // We call tag based storage
160  else if (_tags.size() == 1)
161  {
166  }
167  // This one may be expensive, and hopefully we do not use it so often
168  else
169  {
174  }
175 
176  for (auto & var : needed_moose_vars)
177  var->computingJacobian(true);
178 }
void updateBoundaryVariableDependency(std::set< MooseVariableFEBase *> &needed_moose_vars, THREAD_ID tid=0) const
void updateBlockMatPropDependency(SubdomainID id, std::set< unsigned int > &needed_mat_props, THREAD_ID tid=0) const
virtual void prepareMaterials(SubdomainID blk_id, THREAD_ID tid)
Add the MooseVariables that the current materials depend on to the dependency list.
virtual void setActiveElementalMooseVariables(const std::set< MooseVariableFEBase *> &moose_vars, THREAD_ID tid) override
Set the MOOSE variables to be reinited on each element.
MooseObjectTagWarehouse< KernelBase > & _kernels
void updateBlockVariableDependency(SubdomainID id, std::set< MooseVariableFEBase *> &needed_moose_vars, THREAD_ID tid=0) const
MooseObjectTagWarehouse< IntegratedBCBase > & _integrated_bcs
virtual void setActiveMaterialProperties(const std::set< unsigned int > &mat_prop_ids, THREAD_ID tid) override
Record and set the material properties required by the current computing thread.
virtual void subdomainSetup(SubdomainID subdomain, THREAD_ID tid)
MooseObjectWarehouse< InterfaceKernel > * _ik_warehouse
const std::set< TagID > & _tags
MooseObjectWarehouse< T > & getMatrixTagObjectWarehouse(TagID tag_id, THREAD_ID tid)
Retrieve a moose object warehouse in which every moose object has the given matrix tag...
MooseObjectWarehouse< KernelBase > * _warehouse
MooseObjectWarehouse< IntegratedBCBase > * _ibc_warehouse
MooseObjectWarehouse< T > & getMatrixTagsObjectWarehouse(const std::set< TagID > &tags, THREAD_ID tid)
Retrieve a moose object warehouse in which every moose object has one of the given matrix tags...
MooseObjectWarehouse< DGKernel > * _dg_warehouse
SubdomainID _subdomain
The subdomain for the current element.
virtual unsigned int numMatrixTags()
The total number of tags.
Definition: SubProblem.h:145
MooseObjectTagWarehouse< InterfaceKernel > & _interface_kernels
MooseObjectTagWarehouse< DGKernel > & _dg_kernels
void updateBoundaryMatPropDependency(std::set< unsigned int > &needed_mat_props, THREAD_ID tid=0) const

Member Data Documentation

◆ _blocks

std::vector<JacobianBlock *> ComputeJacobianBlocksThread::_blocks
protected

Definition at line 57 of file ComputeJacobianBlocksThread.h.

Referenced by postElement().

◆ _dg_kernels

MooseObjectTagWarehouse<DGKernel>& ComputeJacobianThread::_dg_kernels
protectedinherited

Definition at line 57 of file ComputeJacobianThread.h.

Referenced by ComputeJacobianThread::subdomainChanged().

◆ _dg_warehouse

MooseObjectWarehouse<DGKernel>* ComputeJacobianThread::_dg_warehouse
protectedinherited

◆ _fe_problem

FEProblemBase& ThreadedElementLoop< ConstElemRange >::_fe_problem
protectedinherited

◆ _ibc_warehouse

MooseObjectWarehouse<IntegratedBCBase>* ComputeJacobianThread::_ibc_warehouse
protectedinherited

◆ _ik_warehouse

MooseObjectWarehouse<InterfaceKernel>* ComputeJacobianThread::_ik_warehouse
protectedinherited

◆ _integrated_bcs

MooseObjectTagWarehouse<IntegratedBCBase>& ComputeJacobianThread::_integrated_bcs
protectedinherited

Definition at line 52 of file ComputeJacobianThread.h.

Referenced by ComputeJacobianThread::subdomainChanged().

◆ _interface_kernels

MooseObjectTagWarehouse<InterfaceKernel>& ComputeJacobianThread::_interface_kernels
protectedinherited

Definition at line 62 of file ComputeJacobianThread.h.

Referenced by ComputeJacobianThread::subdomainChanged().

◆ _kernels

MooseObjectTagWarehouse<KernelBase>& ComputeJacobianThread::_kernels
protectedinherited

Definition at line 67 of file ComputeJacobianThread.h.

Referenced by ComputeJacobianThread::subdomainChanged().

◆ _mesh

MooseMesh& ThreadedElementLoopBase< ConstElemRange >::_mesh
protectedinherited

Definition at line 136 of file ThreadedElementLoopBase.h.

◆ _neighbor_subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_neighbor_subdomain
protectedinherited

The subdomain for the current neighbor.

Definition at line 146 of file ThreadedElementLoopBase.h.

◆ _nl

NonlinearSystemBase& ComputeJacobianThread::_nl
protectedinherited

◆ _num_cached

unsigned int ComputeJacobianThread::_num_cached
protectedinherited

Definition at line 49 of file ComputeJacobianThread.h.

Referenced by ComputeJacobianThread::postElement().

◆ _old_neighbor_subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_old_neighbor_subdomain
protectedinherited

The subdomain for the last neighbor.

Definition at line 149 of file ThreadedElementLoopBase.h.

◆ _old_subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_old_subdomain
protectedinherited

The subdomain for the last element.

Definition at line 143 of file ThreadedElementLoopBase.h.

◆ _subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_subdomain
protectedinherited

◆ _tags

const std::set<TagID>& ComputeJacobianThread::_tags
protectedinherited

Definition at line 72 of file ComputeJacobianThread.h.

Referenced by ComputeJacobianThread::subdomainChanged().

◆ _tid

THREAD_ID ThreadedElementLoopBase< ConstElemRange >::_tid
protectedinherited

Definition at line 137 of file ThreadedElementLoopBase.h.

Referenced by ComputeFullJacobianThread::computeFaceJacobian(), ComputeJacobianThread::computeFaceJacobian(), ComputeFullJacobianThread::computeInternalFaceJacobian(), ComputeJacobianThread::computeInternalFaceJacobian(), ComputeFullJacobianThread::computeInternalInterFaceJacobian(), ComputeJacobianThread::computeInternalInterFaceJacobian(), ComputeFullJacobianThread::computeJacobian(), ComputeJacobianThread::computeJacobian(), ComputeJacobianThread::onBoundary(), ComputeResidualThread::onBoundary(), ComputeMaterialsObjectThread::onBoundary(), ComputeUserObjectsThread::onBoundary(), ComputeMarkerThread::onElement(), ComputeElemDampingThread::onElement(), ComputeElemAuxVarsThread::onElement(), ComputeJacobianThread::onElement(), ComputeIndicatorThread::onElement(), ComputeResidualThread::onElement(), ComputeMaterialsObjectThread::onElement(), ComputeUserObjectsThread::onElement(), ComputeJacobianThread::onInterface(), ComputeResidualThread::onInterface(), ComputeJacobianThread::onInternalSide(), ComputeIndicatorThread::onInternalSide(), ComputeResidualThread::onInternalSide(), ComputeMaterialsObjectThread::onInternalSide(), ComputeUserObjectsThread::onInternalSide(), ComputeMarkerThread::post(), ComputeElemAuxVarsThread::post(), ComputeMaterialsObjectThread::post(), ComputeIndicatorThread::post(), ComputeJacobianThread::post(), ComputeResidualThread::post(), ComputeUserObjectsThread::post(), ComputeJacobianThread::postElement(), ComputeResidualThread::postElement(), postElement(), ComputeUserObjectsThread::queryBoundary(), ComputeUserObjectsThread::querySubdomain(), ComputeMarkerThread::subdomainChanged(), ComputeElemAuxVarsThread::subdomainChanged(), ComputeJacobianThread::subdomainChanged(), ComputeIndicatorThread::subdomainChanged(), ComputeResidualThread::subdomainChanged(), ComputeMaterialsObjectThread::subdomainChanged(), and ComputeUserObjectsThread::subdomainChanged().

◆ _warehouse

MooseObjectWarehouse<KernelBase>* ComputeJacobianThread::_warehouse
protectedinherited

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