www.mooseframework.org
ComputeJacobianThread.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "ComputeJacobianThread.h"
11 
12 #include "DGKernelBase.h"
13 #include "FEProblem.h"
14 #include "IntegratedBCBase.h"
15 #include "InterfaceKernelBase.h"
16 #include "MooseVariableFE.h"
17 #include "NonlinearSystem.h"
18 #include "SwapBackSentinel.h"
19 #include "TimeDerivative.h"
20 #include "FVElementalKernel.h"
21 #include "MaterialBase.h"
22 #include "ConsoleUtils.h"
23 
24 #include "libmesh/threads.h"
25 
27  const std::set<TagID> & tags)
28  : NonlinearThread(fe_problem), _tags(tags)
29 {
30 }
31 
32 // Splitting Constructor
34  : NonlinearThread(x, split), _tags(x._tags)
35 {
36 }
37 
39 
40 void
42 {
43  if (kernel.isImplicit())
44  {
45  kernel.prepareShapes(kernel.variable().number());
46  kernel.computeJacobian();
48  mooseError("Nonlocal kernels only supported for non-diagonal coupling. Please specify an SMP "
49  "preconditioner, with appropriate row-column coupling or specify full = true.");
50  }
51 }
52 
53 void
55 {
56  if (fvkernel.isImplicit())
57  fvkernel.computeJacobian();
58 }
59 
60 void
62 {
63  if (bc.shouldApply() && bc.isImplicit())
64  {
65  bc.prepareShapes(bc.variable().number());
66  bc.computeJacobian();
68  mooseError("Nonlocal boundary conditions only supported for non-diagonal coupling. Please "
69  "specify an SMP preconditioner, with appropriate row-column coupling or specify "
70  "full = true.");
71  }
72 }
73 
74 void
75 ComputeJacobianThread::compute(DGKernelBase & dg, const Elem * neighbor)
76 {
77  if (dg.isImplicit())
78  {
79  dg.prepareShapes(dg.variable().number());
81  if (dg.hasBlocks(neighbor->subdomain_id()))
82  dg.computeJacobian();
83  }
84 }
85 
86 void
88 {
89  if (intk.isImplicit())
90  {
91  intk.prepareShapes(intk.variable().number());
93  intk.computeJacobian();
94  }
95 }
96 
97 void
99 {
100  // If users pass a empty vector or a full size of vector,
101  // we take all kernels
102  if (!_tags.size() || _tags.size() == _fe_problem.numMatrixTags())
103  {
108  }
109  // If we have one tag only,
110  // We call tag based storage
111  else if (_tags.size() == 1)
112  {
117  }
118  // This one may be expensive, and hopefully we do not use it so often
119  else
120  {
125  }
126 }
127 
128 void
130 {
132 }
133 
134 void
136 {
138 }
139 
140 void
142 {
144 }
145 
146 void
148 {
150  _num_cached++;
151 
152  if (_num_cached % 20 == 0)
153  {
154  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
156  }
157 }
158 
159 void
161 {
162 }
virtual void postElement(const Elem *) override
Called after the element assembly is done (including surface assembling)
virtual void addJacobianLowerD(const THREAD_ID tid) override
MooseObjectWarehouse< InterfaceKernelBase > * _ik_warehouse
void determineObjectWarehouses() override
Determine the objects we will actually compute based on vector/matrix tag information.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299
unsigned int number() const
Get variable number coming from libMesh.
MooseObjectTagWarehouse< IntegratedBCBase > & _integrated_bcs
Reference to BC storage structures.
unsigned int _num_cached
virtual void prepareShapes(unsigned int var_num)
Prepare shape functions.
virtual const MooseVariableFEBase & neighborVariable() const =0
The neighbor variable number that this interface kernel operates on.
virtual void computeJacobian() override
Computes the jacobian for the current side.
Definition: DGKernelBase.C:138
void computingScalingJacobian(bool computing_scaling_jacobian)
Setter for whether we&#39;re computing the scaling jacobian.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
Serves as a base class for DGKernel and ADDGKernel.
Definition: DGKernelBase.h:32
const std::set< TagID > & _tags
virtual void accumulateNeighborLower() override
Add neighbor and lower residual/Jacobian into assembly global data.
void prepareShapes(unsigned int var_num) override final
Prepare shape functions.
FVElemental is used for calculating residual contributions from volume integral terms of a PDE where ...
void prepareNeighborShapes(unsigned int var_num)
Prepare neighbor shape functions.
void accumulateNeighbor() override
Add neighbor residual/Jacobian into assembly global data.
This is the common base class for the three main kernel types implemented in MOOSE, Kernel, VectorKernel and ArrayKernel.
Definition: KernelBase.h:22
void prepareShapes(unsigned int var_num) override final
Prepare shape functions.
Definition: DGKernelBase.C:208
MooseObjectTagWarehouse< KernelBase > & _kernels
virtual void addJacobianNeighbor(const THREAD_ID tid) override
virtual bool checkNonlocalCouplingRequirement()
Definition: SubProblem.h:88
MooseObjectWarehouse< T > & getMatrixTagObjectWarehouse(TagID tag_id, THREAD_ID tid)
Retrieve a moose object warehouse in which every moose object has the given matrix tag...
virtual void computeJacobian()=0
Compute this object&#39;s contribution to the diagonal Jacobian entries.
void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
void prepareShapes(unsigned int var_num) override final
Prepare shape functions.
virtual void accumulateLower() override
Add lower-d residual/Jacobian into assembly global data.
virtual unsigned int numMatrixTags() const
The total number of tags.
Definition: SubProblem.h:220
tbb::split split
MooseObjectWarehouse< IntegratedBCBase > * _ibc_warehouse
ComputeJacobianThread(FEProblemBase &fe_problem, const std::set< TagID > &tags)
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...
void join(const ComputeJacobianThread &)
MooseObjectTagWarehouse< InterfaceKernelBase > & _interface_kernels
Reference to interface kernel storage structure.
Base class for deriving any boundary condition of a integrated type.
virtual bool shouldApply()
Hook for turning the boundary condition on and off.
InterfaceKernelBase is the base class for all InterfaceKernel type classes.
virtual const MooseVariableBase & variable() const =0
Returns the variable that this object operates on.
MooseObjectWarehouse< KernelBase > * _tag_kernels
MooseObjectWarehouse< DGKernelBase > * _dg_warehouse
bool hasBlocks(const SubdomainName &name) const
Test if the supplied block name is valid for this object.
virtual void cacheJacobian(const THREAD_ID tid) override
virtual void addJacobianNeighborLowerD(const THREAD_ID tid) override
MooseObjectTagWarehouse< DGKernelBase > & _dg_kernels
Reference to DGKernel storage structure.
virtual void compute(ResidualObject &) override
Will dispatch to computeResidual/computeJacobian/computeResidualAndJacobian based on the derived clas...
virtual void addCachedJacobian(const THREAD_ID tid) override