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

#include <ComputeFullJacobianThread.h>

Inheritance diagram for ComputeFullJacobianThread:
[legend]

Public Member Functions

 ComputeFullJacobianThread (FEProblemBase &fe_problem, const std::set< TagID > &tags)
 
 ComputeFullJacobianThread (ComputeFullJacobianThread &x, Threads::split split)
 
virtual ~ComputeFullJacobianThread ()
 
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 postElement (const Elem *) override
 Called after the element assembly is done (including surface 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 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

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
 
MooseObjectTagWarehouse< KernelBase > & _ad_jacobian_kernels
 Reference to ADKernel<JACOBIAN> storage structure. More...
 
MooseObjectWarehouse< KernelBase > * _adjk_warehouse
 Pointer to tags. More...
 
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

Definition at line 20 of file ComputeFullJacobianThread.h.

Constructor & Destructor Documentation

◆ ComputeFullJacobianThread() [1/2]

ComputeFullJacobianThread::ComputeFullJacobianThread ( FEProblemBase fe_problem,
const std::set< TagID > &  tags 
)

Definition at line 23 of file ComputeFullJacobianThread.C.

25  : ComputeJacobianThread(fe_problem, tags)
26 {
27 }
ComputeJacobianThread(FEProblemBase &fe_problem, const std::set< TagID > &tags)

◆ ComputeFullJacobianThread() [2/2]

ComputeFullJacobianThread::ComputeFullJacobianThread ( ComputeFullJacobianThread x,
Threads::split  split 
)

Definition at line 30 of file ComputeFullJacobianThread.C.

33 {
34 }
static PetscErrorCode Vec x
std::vector< std::string > split(const std::string &str, const std::string &delimiter)
Python like split function for strings.
Definition: MooseUtils.C:782
ComputeJacobianThread(FEProblemBase &fe_problem, const std::set< TagID > &tags)

◆ ~ComputeFullJacobianThread()

ComputeFullJacobianThread::~ComputeFullJacobianThread ( )
virtual

Definition at line 36 of file ComputeFullJacobianThread.C.

36 {}

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

77 {
78  Threads::spin_mutex::scoped_lock lock(threaded_element_mutex);
79 
80  std::string what(e.what());
82 }
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)
overrideprotectedvirtual

done only when nonlocal integrated_bcs exist in the system

Reimplemented from ComputeJacobianThread.

Definition at line 135 of file ComputeFullJacobianThread.C.

136 {
137  std::vector<std::pair<MooseVariableFEBase *, MooseVariableFEBase *>> & ce =
139  for (const auto & it : ce)
140  {
141  MooseVariableFEBase & ivar = *(it.first);
142  MooseVariableFEBase & jvar = *(it.second);
145  {
146  // only if there are dofs for j-variable (if it is subdomain restricted var, there may not be
147  // any)
148  const auto & bcs = _ibc_warehouse->getActiveBoundaryObjects(bnd_id, _tid);
149  for (const auto & bc : bcs)
150  if (bc->shouldApply() && bc->variable().number() == ivar.number() && bc->isImplicit())
151  {
152  bc->subProblem().prepareFaceShapes(jvar.number(), _tid);
153  bc->computeJacobianBlock(jvar);
154  }
155  }
156  }
157 
160  {
161  std::vector<std::pair<MooseVariableFEBase *, MooseVariableFEBase *>> & cne =
163  for (const auto & it : cne)
164  {
165  MooseVariableFEBase & ivariable = *(it.first);
166  MooseVariableFEBase & jvariable = *(it.second);
167 
168  unsigned int ivar = ivariable.number();
169  unsigned int jvar = jvariable.number();
170 
171  if (ivariable.activeOnSubdomain(_subdomain) && jvariable.activeOnSubdomain(_subdomain) &&
173  {
174  const std::vector<std::shared_ptr<IntegratedBCBase>> & integrated_bcs =
176  for (const auto & integrated_bc : integrated_bcs)
177  {
178  std::shared_ptr<NonlocalIntegratedBC> nonlocal_integrated_bc =
180  if (nonlocal_integrated_bc)
181  if ((integrated_bc->variable().number() == ivar) && integrated_bc->isImplicit())
182  {
183  integrated_bc->subProblem().prepareFaceShapes(jvar, _tid);
184  integrated_bc->computeNonlocalOffDiagJacobian(jvar);
185  }
186  }
187  }
188  }
189  }
190 
191  const std::vector<MooseVariableScalar *> & scalar_vars = _nl.getScalarVariables(_tid);
192  if (scalar_vars.size() > 0)
193  {
194  // go over nl-variables (non-scalar)
195  const std::vector<MooseVariableFEBase *> & vars = _nl.getVariables(_tid);
196  for (const auto & ivar : vars)
197  if (ivar->activeOnSubdomain(_subdomain) > 0 &&
199  {
200  // for each variable get the list of active kernels
201  const auto & bcs = _ibc_warehouse->getActiveBoundaryObjects(bnd_id, _tid);
202  for (const auto & bc : bcs)
203  if (bc->variable().number() == ivar->number() && bc->isImplicit())
204  {
205  // now, get the list of coupled scalar vars and compute their off-diag jacobians
206  const std::vector<MooseVariableScalar *> coupled_scalar_vars =
207  bc->getCoupledMooseScalarVars();
208 
209  // Do: dvar / dscalar_var, only want to process only nl-variables (not aux ones)
210  for (const auto & jvar : coupled_scalar_vars)
211  if (_nl.hasScalarVariable(jvar->name()))
212  bc->computeJacobianBlockScalar(jvar->number());
213  }
214  }
215  }
216 }
SubProblem & subProblem()
Get a reference to the subproblem.
const std::vector< MooseVariableScalar * > & getScalarVariables(THREAD_ID tid)
Definition: SystemBase.h:574
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:569
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:69
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:661

◆ computeInternalFaceJacobian()

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

Reimplemented from ComputeJacobianThread.

Definition at line 219 of file ComputeFullJacobianThread.C.

220 {
222  {
223  const auto & ce = _fe_problem.couplingEntries(_tid);
224  for (const auto & it : ce)
225  {
226  const auto & dgks = _dg_warehouse->getActiveBlockObjects(_subdomain, _tid);
227  for (const auto & dg : dgks)
228  {
229  MooseVariableFEBase & ivariable = *(it.first);
230  MooseVariableFEBase & jvariable = *(it.second);
231 
232  unsigned int ivar = ivariable.number();
233  unsigned int jvar = jvariable.number();
234 
235  if (dg->variable().number() == ivar && dg->isImplicit() &&
236  dg->hasBlocks(neighbor->subdomain_id()) && jvariable.activeOnSubdomain(_subdomain))
237  {
238  dg->subProblem().prepareFaceShapes(jvar, _tid);
239  dg->subProblem().prepareNeighborShapes(jvar, _tid);
240  dg->computeOffDiagJacobian(jvar);
241  }
242  }
243  }
244  }
245 }
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)
overrideprotectedvirtual

Reimplemented from ComputeJacobianThread.

Definition at line 248 of file ComputeFullJacobianThread.C.

249 {
251  {
252  const auto & ce = _fe_problem.couplingEntries(_tid);
253  for (const auto & it : ce)
254  {
255  const std::vector<std::shared_ptr<InterfaceKernel>> & int_ks =
257  for (const auto & interface_kernel : int_ks)
258  {
259  if (!interface_kernel->isImplicit())
260  continue;
261 
262  unsigned int ivar = it.first->number();
263  unsigned int jvar = it.second->number();
264 
265  interface_kernel->subProblem().prepareFaceShapes(jvar, _tid);
266  interface_kernel->subProblem().prepareNeighborShapes(jvar, _tid);
267 
268  if (interface_kernel->variable().number() == ivar)
269  interface_kernel->computeElementOffDiagJacobian(jvar);
270 
271  if (interface_kernel->neighborVariable().number() == ivar)
272  interface_kernel->computeNeighborOffDiagJacobian(jvar);
273  }
274  }
275  }
276 }
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 ( )
overrideprotectedvirtual

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 
67  {
69  for (const auto & kernel : kernels)
70  if (kernel->isImplicit())
71  kernel->computeADOffDiagJacobian();
72  }
73 
76  {
77  std::vector<std::pair<MooseVariableFEBase *, MooseVariableFEBase *>> & cne =
79  for (const auto & it : cne)
80  {
81  MooseVariableFEBase & ivariable = *(it.first);
82  MooseVariableFEBase & jvariable = *(it.second);
83 
84  unsigned int ivar = ivariable.number();
85  unsigned int jvar = jvariable.number();
86 
87  if (ivariable.activeOnSubdomain(_subdomain) && jvariable.activeOnSubdomain(_subdomain) &&
89  {
90  const auto & kernels = _warehouse->getActiveVariableBlockObjects(ivar, _subdomain, _tid);
91  for (const auto & kernel : kernels)
92  {
93  std::shared_ptr<NonlocalKernel> nonlocal_kernel =
95  if (nonlocal_kernel)
96  if ((kernel->variable().number() == ivar) && kernel->isImplicit())
97  {
98  kernel->subProblem().prepareShapes(jvar, _tid);
99  kernel->computeNonlocalOffDiagJacobian(jvar);
100  }
101  }
102  }
103  }
104  }
105 
106  const std::vector<MooseVariableScalar *> & scalar_vars = _nl.getScalarVariables(_tid);
107  if (scalar_vars.size() > 0)
108  {
109  // go over nl-variables (non-scalar)
110  const std::vector<MooseVariableFEBase *> & vars = _nl.getVariables(_tid);
111  for (const auto & ivariable : vars)
112  if (ivariable->activeOnSubdomain(_subdomain) > 0 &&
114  {
115  // for each variable get the list of active kernels
116  const auto & kernels =
118  for (const auto & kernel : kernels)
119  if (kernel->isImplicit())
120  {
121  // now, get the list of coupled scalar vars and compute their off-diag jacobians
122  const std::vector<MooseVariableScalar *> coupled_scalar_vars =
123  kernel->getCoupledMooseScalarVars();
124 
125  // Do: dvar / dscalar_var, only want to process only nl-variables (not aux ones)
126  for (const auto & jvariable : coupled_scalar_vars)
127  if (_nl.hasScalarVariable(jvariable->name()))
128  kernel->computeOffDiagJacobianScalar(jvariable->number());
129  }
130  }
131  }
132 }
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:574
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 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:569
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.
MooseObjectWarehouse< KernelBase > * _adjk_warehouse
Pointer to tags.
virtual bool checkNonlocalCouplingRequirement()
Definition: SubProblem.h:69
MooseObjectWarehouse< KernelBase > * _warehouse
SubProblem & subProblem()
Returns a reference to the SubProblem for which this Kernel is active.
Definition: KernelBase.h:106
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:661

◆ join()

void ComputeFullJacobianThread::join ( const ComputeJacobianThread )
inline

Definition at line 30 of file ComputeFullJacobianThread.h.

30 {}

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

44 { 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 100 of file ThreadedElementLoop.h.

101 {
104 }
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 210 of file ComputeJacobianThread.C.

211 {
213  {
214  _fe_problem.reinitElemFace(elem, side, bnd_id, _tid);
215 
216  // Set up Sentinel class so that, even if reinitMaterials() throws, we
217  // still remember to swap back during stack unwinding.
219 
220  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
222 
223  computeFaceJacobian(bnd_id);
224  }
225 }
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 192 of file ComputeJacobianThread.C.

193 {
194  _fe_problem.prepare(elem, _tid);
195 
196  _fe_problem.reinitElem(elem, _tid);
197 
198  // Set up Sentinel class so that, even if reinitMaterials() throws, we
199  // still remember to swap back during stack unwinding.
202 
203  if (_nl.getScalarVariables(_tid).size() > 0)
205 
206  computeJacobian();
207 }
virtual void prepare(const Elem *elem, THREAD_ID tid) override
const std::vector< MooseVariableScalar * > & getScalarVariables(THREAD_ID tid)
Definition: SystemBase.h:574
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 262 of file ComputeJacobianThread.C.

263 {
265  {
266  // Pointer to the neighbor we are currently working on.
267  const Elem * neighbor = elem->neighbor_ptr(side);
268 
269  if (neighbor->active())
270  {
271  _fe_problem.reinitNeighbor(elem, side, _tid);
272 
273  // Set up Sentinels so that, even if one of the reinitMaterialsXXX() calls throws, we
274  // still remember to swap back during stack unwinding.
276  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
278 
280  _fe_problem.reinitMaterialsNeighbor(neighbor->subdomain_id(), _tid);
281 
283 
284  {
285  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
287  }
288  }
289  }
290 }
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 228 of file ComputeJacobianThread.C.

229 {
231  {
232  // Pointer to the neighbor we are currently working on.
233  const Elem * neighbor = elem->neighbor_ptr(side);
234 
235  // Get the global id of the element and the neighbor
236  const dof_id_type elem_id = elem->id(), neighbor_id = neighbor->id();
237 
238  if ((neighbor->active() && (neighbor->level() == elem->level()) && (elem_id < neighbor_id)) ||
239  (neighbor->level() < elem->level()))
240  {
241  _fe_problem.reinitNeighbor(elem, side, _tid);
242 
243  // Set up Sentinels so that, even if one of the reinitMaterialsXXX() calls throws, we
244  // still remember to swap back during stack unwinding.
246  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
247 
249  _fe_problem.reinitMaterialsNeighbor(neighbor->subdomain_id(), _tid);
250 
251  computeInternalFaceJacobian(neighbor);
252 
253  {
254  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
256  }
257  }
258  }
259 }
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 173 of file ThreadedElementLoopBase.h.

174 {
175  try
176  {
177  try
178  {
179  ParallelUniqueId puid;
180  _tid = bypass_threading ? 0 : puid.id;
181 
182  pre();
183 
186  typename RangeType::const_iterator el = range.begin();
187  for (el = range.begin(); el != range.end(); ++el)
188  {
189  if (!keepGoing())
190  break;
191 
192  const Elem * elem = *el;
193 
194  preElement(elem);
195 
197  _subdomain = elem->subdomain_id();
198  if (_subdomain != _old_subdomain)
200 
201  onElement(elem);
202 
203  for (unsigned int side = 0; side < elem->n_sides(); side++)
204  {
205  std::vector<BoundaryID> boundary_ids = _mesh.getBoundaryIDs(elem, side);
206 
207  if (boundary_ids.size() > 0)
208  for (std::vector<BoundaryID>::iterator it = boundary_ids.begin();
209  it != boundary_ids.end();
210  ++it)
211  onBoundary(elem, side, *it);
212 
213  const Elem * neighbor = elem->neighbor_ptr(side);
214  if (neighbor != nullptr)
215  {
216  preInternalSide(elem, side);
217 
219  _neighbor_subdomain = neighbor->subdomain_id();
222 
223  onInternalSide(elem, side);
224 
225  if (boundary_ids.size() > 0)
226  for (std::vector<BoundaryID>::iterator it = boundary_ids.begin();
227  it != boundary_ids.end();
228  ++it)
229  onInterface(elem, side, *it);
230 
231  postInternalSide(elem, side);
232  }
233  } // sides
234  postElement(elem);
235 
236  } // range
237 
238  post();
239  }
240  catch (libMesh::LogicError & e)
241  {
242  throw MooseException("We caught a libMesh error");
243  }
244  }
245  catch (MooseException & e)
246  {
248  }
249 }
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.C:16
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:2148
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 306 of file ComputeJacobianThread.C.

307 {
310 }
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 ComputeJacobianThread::postElement ( const Elem *  elem)
overridevirtualinherited

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

Parameters
elem- active element

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Reimplemented in ComputeJacobianBlocksThread.

Definition at line 293 of file ComputeJacobianThread.C.

294 {
296  _num_cached++;
297 
298  if (_num_cached % 20 == 0)
299  {
300  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
302  }
303 }
virtual void addCachedJacobian(THREAD_ID tid) override
virtual void cacheJacobian(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 297 of file ThreadedElementLoopBase.h.

298 {
299 }

◆ pre()

void ThreadedElementLoopBase< ConstElemRange >::pre ( )
virtualinherited

Called before the element range loop.

Definition at line 253 of file ThreadedElementLoopBase.h.

254 {
255 }

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

87 {
89 }
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 93 of file ThreadedElementLoop.h.

94 {
96 }
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 136 of file ComputeJacobianThread.C.

137 {
139 
140  // Update variable Dependencies
141  std::set<MooseVariableFEBase *> needed_moose_vars;
147 
148  // Update material dependencies
149  std::set<unsigned int> needed_mat_props;
155 
157  _fe_problem.setActiveMaterialProperties(needed_mat_props, _tid);
159 
160  // If users pass a empty vector or a full size of vector,
161  // we take all kernels
162  if (!_tags.size() || _tags.size() == _fe_problem.numMatrixTags())
163  {
164  _warehouse = &_kernels;
169  }
170  // If we have one tag only,
171  // We call tag based storage
172  else if (_tags.size() == 1)
173  {
179  }
180  // This one may be expensive, and hopefully we do not use it so often
181  else
182  {
188  }
189 }
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
MooseObjectTagWarehouse< KernelBase > & _ad_jacobian_kernels
Reference to ADKernel<JACOBIAN> storage structure.
MooseObjectWarehouse< KernelBase > * _adjk_warehouse
Pointer to 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:152
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

◆ _ad_jacobian_kernels

MooseObjectTagWarehouse<KernelBase>& ComputeJacobianThread::_ad_jacobian_kernels
protectedinherited

Reference to ADKernel<JACOBIAN> storage structure.

Definition at line 73 of file ComputeJacobianThread.h.

Referenced by ComputeJacobianThread::subdomainChanged().

◆ _adjk_warehouse

MooseObjectWarehouse<KernelBase>* ComputeJacobianThread::_adjk_warehouse
protectedinherited

◆ _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 138 of file ThreadedElementLoopBase.h.

◆ _neighbor_subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_neighbor_subdomain
protectedinherited

The subdomain for the current neighbor.

Definition at line 148 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 151 of file ThreadedElementLoopBase.h.

◆ _old_subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_old_subdomain
protectedinherited

The subdomain for the last element.

Definition at line 145 of file ThreadedElementLoopBase.h.

◆ _subdomain

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_subdomain
protectedinherited

◆ _tags

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

Definition at line 78 of file ComputeJacobianThread.h.

Referenced by ComputeJacobianThread::subdomainChanged().

◆ _tid

THREAD_ID ThreadedElementLoopBase< ConstElemRange >::_tid
protectedinherited

Definition at line 139 of file ThreadedElementLoopBase.h.

Referenced by computeFaceJacobian(), ComputeJacobianThread::computeFaceJacobian(), computeInternalFaceJacobian(), ComputeJacobianThread::computeInternalFaceJacobian(), computeInternalInterFaceJacobian(), ComputeJacobianThread::computeInternalInterFaceJacobian(), 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(), ComputeJacobianBlocksThread::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: