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

#include <ComputeJacobianThread.h>

Inheritance diagram for ComputeJacobianThread:
[legend]

Public Member Functions

 ComputeJacobianThread (FEProblemBase &fe_problem, SparseMatrix< Number > &jacobian, Moose::KernelType kernel_type=Moose::KT_ALL)
 
 ComputeJacobianThread (ComputeJacobianThread &x, Threads::split split)
 
virtual ~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...
 
void join (const ComputeJacobianThread &)
 
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 ()
 
virtual void computeFaceJacobian (BoundaryID bnd_id)
 
virtual void computeInternalFaceJacobian (const Elem *neighbor)
 
virtual void computeInternalInterFaceJacobian (BoundaryID bnd_id)
 

Protected Attributes

SparseMatrix< Number > & _jacobian
 
NonlinearSystemBase_nl
 
unsigned int _num_cached
 
const MooseObjectWarehouse< IntegratedBC > & _integrated_bcs
 
const MooseObjectWarehouse< DGKernel > & _dg_kernels
 
const MooseObjectWarehouse< InterfaceKernel > & _interface_kernels
 
const KernelWarehouse_kernels
 
Moose::KernelType _kernel_type
 
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 30 of file ComputeJacobianThread.h.

Constructor & Destructor Documentation

ComputeJacobianThread::ComputeJacobianThread ( FEProblemBase fe_problem,
SparseMatrix< Number > &  jacobian,
Moose::KernelType  kernel_type = Moose::KT_ALL 
)

Definition at line 31 of file ComputeJacobianThread.C.

35  _jacobian(jacobian),
36  _nl(fe_problem.getNonlinearSystemBase()),
37  _num_cached(0),
42  _kernel_type(kernel_type)
43 {
44 }
Moose::KernelType _kernel_type
SparseMatrix< Number > & _jacobian
NonlinearSystemBase & getNonlinearSystemBase()
const MooseObjectWarehouse< IntegratedBC > & getIntegratedBCWarehouse()
const MooseObjectWarehouse< DGKernel > & getDGKernelWarehouse()
NonlinearSystemBase & _nl
const MooseObjectWarehouse< DGKernel > & _dg_kernels
const KernelWarehouse & _kernels
const MooseObjectWarehouse< InterfaceKernel > & getInterfaceKernelWarehouse()
const MooseObjectWarehouse< InterfaceKernel > & _interface_kernels
const MooseObjectWarehouse< IntegratedBC > & _integrated_bcs
const KernelWarehouse & getKernelWarehouse()
Access functions to Warehouses from outside NonlinearSystemBase.
ComputeJacobianThread::ComputeJacobianThread ( ComputeJacobianThread x,
Threads::split  split 
)

Definition at line 47 of file ComputeJacobianThread.C.

49  _jacobian(x._jacobian),
50  _nl(x._nl),
51  _num_cached(x._num_cached),
52  _integrated_bcs(x._integrated_bcs),
53  _dg_kernels(x._dg_kernels),
54  _interface_kernels(x._interface_kernels),
55  _kernels(x._kernels),
56  _kernel_type(x._kernel_type)
57 {
58 }
Moose::KernelType _kernel_type
SparseMatrix< Number > & _jacobian
static PetscErrorCode Vec x
NonlinearSystemBase & _nl
const MooseObjectWarehouse< DGKernel > & _dg_kernels
const KernelWarehouse & _kernels
const MooseObjectWarehouse< InterfaceKernel > & _interface_kernels
const MooseObjectWarehouse< IntegratedBC > & _integrated_bcs
ComputeJacobianThread::~ComputeJacobianThread ( )
virtual

Definition at line 60 of file ComputeJacobianThread.C.

60 {}

Member Function Documentation

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

void ComputeJacobianThread::computeFaceJacobian ( BoundaryID  bnd_id)
protectedvirtual

done only when nonlocal integrated_bcs exist in the system

Reimplemented in ComputeFullJacobianThread.

Definition at line 114 of file ComputeJacobianThread.C.

Referenced by onBoundary().

115 {
116  const std::vector<std::shared_ptr<IntegratedBC>> & bcs =
118  for (const auto & bc : bcs)
119  if (bc->shouldApply() && bc->isImplicit())
120  {
121  bc->subProblem().prepareFaceShapes(bc->variable().number(), _tid);
122  bc->computeJacobian();
125  {
126  std::shared_ptr<NonlocalIntegratedBC> nonlocal_integrated_bc =
127  std::dynamic_pointer_cast<NonlocalIntegratedBC>(bc);
128  if (nonlocal_integrated_bc)
130  }
131  }
132 }
virtual void computeNonlocalJacobian()
computeNonlocalJacobian and computeNonlocalOffDiagJacobian methods are introduced for providing the j...
const std::map< BoundaryID, std::vector< std::shared_ptr< T > > > & getActiveBoundaryObjects(THREAD_ID tid=0) const
NonlocalIntegratedBC is used for solving integral terms in integro-differential equations.
virtual bool checkNonlocalCouplingRequirement()
Definition: SubProblem.h:62
const MooseObjectWarehouse< IntegratedBC > & _integrated_bcs
void ComputeJacobianThread::computeInternalFaceJacobian ( const Elem *  neighbor)
protectedvirtual

Reimplemented in ComputeFullJacobianThread.

Definition at line 135 of file ComputeJacobianThread.C.

Referenced by onInternalSide().

136 {
137  // No need to call hasActiveObjects, this is done in the calling method (see onInternalSide)
138  const std::vector<std::shared_ptr<DGKernel>> & dgks =
140  for (const auto & dg : dgks)
141  if (dg->isImplicit())
142  {
143  dg->subProblem().prepareFaceShapes(dg->variable().number(), _tid);
144  dg->subProblem().prepareNeighborShapes(dg->variable().number(), _tid);
145  if (dg->hasBlocks(neighbor->subdomain_id()))
146  dg->computeJacobian();
147  }
148 }
const std::map< SubdomainID, std::vector< std::shared_ptr< T > > > & getActiveBlockObjects(THREAD_ID tid=0) const
const MooseObjectWarehouse< DGKernel > & _dg_kernels
SubdomainID _subdomain
The subdomain for the current element.
void ComputeJacobianThread::computeInternalInterFaceJacobian ( BoundaryID  bnd_id)
protectedvirtual

Reimplemented in ComputeFullJacobianThread.

Definition at line 151 of file ComputeJacobianThread.C.

Referenced by onInterface().

152 {
153  // No need to call hasActiveObjects, this is done in the calling method (see onInterface)
154  const std::vector<std::shared_ptr<InterfaceKernel>> & intks =
156  for (const auto & intk : intks)
157  if (intk->isImplicit())
158  {
159  intk->subProblem().prepareFaceShapes(intk->variable().number(), _tid);
160  intk->subProblem().prepareNeighborShapes(intk->neighborVariable().number(), _tid);
161  intk->computeJacobian();
162  }
163 }
const std::map< BoundaryID, std::vector< std::shared_ptr< T > > > & getActiveBoundaryObjects(THREAD_ID tid=0) const
const MooseObjectWarehouse< InterfaceKernel > & _interface_kernels
void ComputeJacobianThread::computeJacobian ( )
protectedvirtual

done only when nonlocal kernels exist in the system

Reimplemented in ComputeFullJacobianThread.

Definition at line 63 of file ComputeJacobianThread.C.

Referenced by onElement().

64 {
65  const MooseObjectWarehouse<KernelBase> * warehouse;
66  switch (_kernel_type)
67  {
68  case Moose::KT_ALL:
69  warehouse = &_nl.getKernelWarehouse();
70  break;
71 
72  case Moose::KT_TIME:
73  warehouse = &_nl.getTimeKernelWarehouse();
74  break;
75 
76  case Moose::KT_NONTIME:
77  warehouse = &_nl.getNonTimeKernelWarehouse();
78  break;
79 
80  case Moose::KT_EIGEN:
81  warehouse = &_nl.getEigenKernelWarehouse();
82  break;
83 
84  case Moose::KT_NONEIGEN:
85  warehouse = &_nl.getNonEigenKernelWarehouse();
86  break;
87 
88  default:
89  mooseError("Unknown kernel type \n");
90  }
91 
92  if (warehouse->hasActiveBlockObjects(_subdomain, _tid))
93  {
94  const std::vector<std::shared_ptr<KernelBase>> & kernels =
96  for (const auto & kernel : kernels)
97  if (kernel->isImplicit())
98  {
99  kernel->subProblem().prepareShapes(kernel->variable().number(), _tid);
100  kernel->computeJacobian();
103  {
104  std::shared_ptr<NonlocalKernel> nonlocal_kernel =
105  std::dynamic_pointer_cast<NonlocalKernel>(kernel);
106  if (nonlocal_kernel)
107  kernel->computeNonlocalJacobian();
108  }
109  }
110  }
111 }
Moose::KernelType _kernel_type
const KernelWarehouse & getNonEigenKernelWarehouse()
void mooseError(Args &&...args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:182
const KernelWarehouse & getEigenKernelWarehouse()
virtual bool checkNonlocalCouplingRequirement()
Definition: SubProblem.h:62
NonlinearSystemBase & _nl
const std::map< SubdomainID, std::vector< std::shared_ptr< T > > > & getActiveBlockObjects(THREAD_ID tid=0) const
NonlocalKernel is used for solving integral terms in integro-differential equations.
bool hasActiveBlockObjects(THREAD_ID tid=0) const
SubdomainID _subdomain
The subdomain for the current element.
const KernelWarehouse & getTimeKernelWarehouse()
virtual void computeNonlocalJacobian()
computeNonlocalJacobian and computeNonlocalOffDiagJacobian methods are introduced for providing the j...
const KernelWarehouse & getNonTimeKernelWarehouse()
const KernelWarehouse & getKernelWarehouse()
Access functions to Warehouses from outside NonlinearSystemBase.
void ComputeJacobianThread::join ( const ComputeJacobianThread )

Definition at line 310 of file ComputeJacobianThread.C.

311 {
312 }
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 50 of file ThreadedElementLoop.h.

50 { return !_fe_problem.hasException(); }
virtual bool hasException()
Whether or not an exception has occurred.
virtual 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 >.

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

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

209 {
211  {
212  _fe_problem.reinitElemFace(elem, side, bnd_id, _tid);
213 
214  // Set up Sentinel class so that, even if reinitMaterials() throws, we
215  // still remember to swap back during stack unwinding.
217 
218  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
220 
221  computeFaceJacobian(bnd_id);
222  }
223 }
virtual void reinitMaterialsBoundary(BoundaryID boundary_id, THREAD_ID tid, bool swap_stateful=true)
bool hasActiveBoundaryObjects(THREAD_ID tid=0) const
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)
const MooseObjectWarehouse< IntegratedBC > & _integrated_bcs
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}.
void ComputeJacobianThread::onElement ( const Elem *  elem)
overridevirtual

Assembly of the element (not including surface assembly)

Parameters
elem- active element

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 190 of file ComputeJacobianThread.C.

191 {
192  _fe_problem.prepare(elem, _tid);
193 
194  _fe_problem.reinitElem(elem, _tid);
195 
196  // Set up Sentinel class so that, even if reinitMaterials() throws, we
197  // still remember to swap back during stack unwinding.
200 
201  if (_nl.getScalarVariables(_tid).size() > 0)
203 
204  computeJacobian();
205 }
virtual void prepare(const Elem *elem, THREAD_ID tid) override
const std::vector< MooseVariableScalar * > & getScalarVariables(THREAD_ID tid)
Definition: SystemBase.h:422
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}.
void ComputeJacobianThread::onInterface ( const Elem *  elem,
unsigned int  side,
BoundaryID  bnd_id 
)
overridevirtual

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

261 {
263  {
264  // Pointer to the neighbor we are currently working on.
265  const Elem * neighbor = elem->neighbor_ptr(side);
266 
267  if (neighbor->active())
268  {
269  _fe_problem.reinitNeighbor(elem, side, _tid);
270 
271  // Set up Sentinels so that, even if one of the reinitMaterialsXXX() calls throws, we
272  // still remember to swap back during stack unwinding.
274  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
275 
277  _fe_problem.reinitMaterialsNeighbor(neighbor->subdomain_id(), _tid);
278 
280 
281  {
282  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
284  }
285  }
286  }
287 }
SparseMatrix< Number > & _jacobian
virtual void reinitMaterialsNeighbor(SubdomainID blk_id, THREAD_ID tid, bool swap_stateful=true)
virtual void computeInternalInterFaceJacobian(BoundaryID bnd_id)
bool hasActiveBoundaryObjects(THREAD_ID tid=0) const
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 addJacobianNeighbor(SparseMatrix< Number > &jacobian, THREAD_ID tid) override
const MooseObjectWarehouse< InterfaceKernel > & _interface_kernels
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}.
void ComputeJacobianThread::onInternalSide ( const Elem *  elem,
unsigned int  side 
)
overridevirtual

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

227 {
229  {
230  // Pointer to the neighbor we are currently working on.
231  const Elem * neighbor = elem->neighbor_ptr(side);
232 
233  // Get the global id of the element and the neighbor
234  const dof_id_type elem_id = elem->id(), neighbor_id = neighbor->id();
235 
236  if ((neighbor->active() && (neighbor->level() == elem->level()) && (elem_id < neighbor_id)) ||
237  (neighbor->level() < elem->level()))
238  {
239  _fe_problem.reinitNeighbor(elem, side, _tid);
240 
241  // Set up Sentinels so that, even if one of the reinitMaterialsXXX() calls throws, we
242  // still remember to swap back during stack unwinding.
244  _fe_problem.reinitMaterialsFace(elem->subdomain_id(), _tid);
245 
247  _fe_problem.reinitMaterialsNeighbor(neighbor->subdomain_id(), _tid);
248 
249  computeInternalFaceJacobian(neighbor);
250 
251  {
252  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
254  }
255  }
256  }
257 }
SparseMatrix< Number > & _jacobian
virtual void reinitMaterialsNeighbor(SubdomainID blk_id, THREAD_ID tid, bool swap_stateful=true)
virtual void swapBackMaterialsNeighbor(THREAD_ID tid)
const MooseObjectWarehouse< DGKernel > & _dg_kernels
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
bool hasActiveBlockObjects(THREAD_ID tid=0) const
virtual void computeInternalFaceJacobian(const Elem *neighbor)
virtual void addJacobianNeighbor(SparseMatrix< Number > &jacobian, THREAD_ID tid) override
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}.
void ThreadedElementLoopBase< ConstElemRange >::operator() ( const ConstElemRange &  range,
bool  bypass_threading = false 
)
inherited
void ComputeJacobianThread::post ( )
overridevirtual

Called after the element range loop.

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Definition at line 303 of file ComputeJacobianThread.C.

304 {
307 }
virtual void clearActiveMaterialProperties(THREAD_ID tid) override
Clear the active material properties.
virtual void clearActiveElementalMooseVariables(THREAD_ID tid) override
Clear the active elemental MooseVariable.
void ComputeJacobianThread::postElement ( const Elem *  elem)
overridevirtual

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

Parameters
elem- active element

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

Reimplemented in ComputeJacobianBlocksThread.

Definition at line 290 of file ComputeJacobianThread.C.

291 {
293  _num_cached++;
294 
295  if (_num_cached % 20 == 0)
296  {
297  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
299  }
300 }
SparseMatrix< Number > & _jacobian
virtual void cacheJacobian(THREAD_ID tid) override
virtual void addCachedJacobian(SparseMatrix< Number > &jacobian, THREAD_ID tid) override
virtual 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'
virtual void ThreadedElementLoopBase< ConstElemRange >::pre ( )
virtualinherited

Called before the element range loop.

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

Called before the element assembly.

Parameters
elem- active element

Reimplemented from ThreadedElementLoopBase< ConstElemRange >.

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

void ComputeJacobianThread::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 166 of file ComputeJacobianThread.C.

167 {
169 
170  // Update variable Dependencies
171  std::set<MooseVariable *> needed_moose_vars;
176 
177  // Update material dependencies
178  std::set<unsigned int> needed_mat_props;
183 
185  _fe_problem.setActiveMaterialProperties(needed_mat_props, _tid);
187 }
void updateBlockVariableDependency(SubdomainID id, std::set< MooseVariable * > &needed_moose_vars, THREAD_ID tid=0) const
virtual void setActiveElementalMooseVariables(const std::set< MooseVariable * > &moose_vars, THREAD_ID tid) override
Set the MOOSE variables to be reinited on each element.
virtual void prepareMaterials(SubdomainID blk_id, THREAD_ID tid)
Add the MooseVariables that the current materials depend on to the dependency list.
void updateBlockMatPropDependency(SubdomainID id, std::set< unsigned int > &needed_mat_props, THREAD_ID tid=0) const
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)
void updateBoundaryVariableDependency(std::set< MooseVariable * > &needed_moose_vars, THREAD_ID tid=0) const
void updateBoundaryMatPropDependency(std::set< unsigned int > &needed_mat_props, THREAD_ID tid=0) const
const MooseObjectWarehouse< DGKernel > & _dg_kernels
SubdomainID _subdomain
The subdomain for the current element.
const KernelWarehouse & _kernels
const MooseObjectWarehouse< InterfaceKernel > & _interface_kernels
const MooseObjectWarehouse< IntegratedBC > & _integrated_bcs

Member Data Documentation

const MooseObjectWarehouse<DGKernel>& ComputeJacobianThread::_dg_kernels
protected
FEProblemBase& ThreadedElementLoop< ConstElemRange >::_fe_problem
protectedinherited
const MooseObjectWarehouse<IntegratedBC>& ComputeJacobianThread::_integrated_bcs
protected

Definition at line 59 of file ComputeJacobianThread.h.

Referenced by computeFaceJacobian(), onBoundary(), and subdomainChanged().

const MooseObjectWarehouse<InterfaceKernel>& ComputeJacobianThread::_interface_kernels
protected
SparseMatrix<Number>& ComputeJacobianThread::_jacobian
protected

Definition at line 53 of file ComputeJacobianThread.h.

Referenced by onInterface(), onInternalSide(), and postElement().

Moose::KernelType ComputeJacobianThread::_kernel_type
protected

Definition at line 70 of file ComputeJacobianThread.h.

Referenced by computeJacobian().

const KernelWarehouse& ComputeJacobianThread::_kernels
protected

Definition at line 68 of file ComputeJacobianThread.h.

Referenced by subdomainChanged().

MooseMesh& ThreadedElementLoopBase< ConstElemRange >::_mesh
protectedinherited

Definition at line 141 of file ThreadedElementLoopBase.h.

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_neighbor_subdomain
protectedinherited

The subdomain for the current neighbor.

Definition at line 151 of file ThreadedElementLoopBase.h.

NonlinearSystemBase& ComputeJacobianThread::_nl
protected

Definition at line 54 of file ComputeJacobianThread.h.

Referenced by computeJacobian(), and onElement().

unsigned int ComputeJacobianThread::_num_cached
protected

Definition at line 56 of file ComputeJacobianThread.h.

Referenced by postElement().

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_old_neighbor_subdomain
protectedinherited

The subdomain for the last neighbor.

Definition at line 154 of file ThreadedElementLoopBase.h.

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_old_subdomain
protectedinherited

The subdomain for the last element.

Definition at line 148 of file ThreadedElementLoopBase.h.

SubdomainID ThreadedElementLoopBase< ConstElemRange >::_subdomain
protectedinherited
THREAD_ID ThreadedElementLoopBase< ConstElemRange >::_tid
protectedinherited

Definition at line 142 of file ThreadedElementLoopBase.h.

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


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