https://mooseframework.inl.gov
ComputeJacobianThread.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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 "HDGKernel.h"
17 #include "MooseVariableFE.h"
18 #include "NonlinearSystem.h"
19 #include "SwapBackSentinel.h"
20 #include "TimeDerivative.h"
21 #include "FVElementalKernel.h"
22 #include "MaterialBase.h"
23 #include "ConsoleUtils.h"
24 #include "Assembly.h"
25 
26 #include "libmesh/threads.h"
27 
29  const std::set<TagID> & tags)
30  : NonlinearThread(fe_problem), _tags(tags)
31 {
32 }
33 
34 // Splitting Constructor
36  : NonlinearThread(x, split), _tags(x._tags)
37 {
38 }
39 
41 
42 void
44 {
45  if (kernel.isImplicit())
46  {
47  kernel.prepareShapes(kernel.variable().number());
48  kernel.computeJacobian();
50  mooseError("Nonlocal kernels only supported for non-diagonal coupling. Please specify an SMP "
51  "preconditioner, with appropriate row-column coupling or specify full = true.");
52  }
53 }
54 
55 void
57 {
58  if (fvkernel.isImplicit())
59  fvkernel.computeJacobian();
60 }
61 
62 void
64 {
65  if (bc.shouldApply() && bc.isImplicit())
66  {
67  bc.prepareShapes(bc.variable().number());
68  bc.computeJacobian();
70  mooseError("Nonlocal boundary conditions only supported for non-diagonal coupling. Please "
71  "specify an SMP preconditioner, with appropriate row-column coupling or specify "
72  "full = true.");
73  }
74 }
75 
76 void
77 ComputeJacobianThread::compute(DGKernelBase & dg, const Elem * neighbor)
78 {
79  if (dg.isImplicit())
80  {
81  dg.prepareShapes(dg.variable().number());
83  if (dg.hasBlocks(neighbor->subdomain_id()))
84  dg.computeJacobian();
85  }
86 }
87 
88 void
90 {
91  if (intk.isImplicit())
92  {
93  intk.prepareShapes(intk.variable().number());
95  intk.computeJacobian();
96  }
97 }
98 
99 void
101 {
102  // If users pass a empty vector or a full size of vector,
103  // we take all kernels
104  if (!_tags.size() || _tags.size() == _fe_problem.numMatrixTags())
105  {
111  }
112  // If we have one tag only,
113  // We call tag based storage
114  else if (_tags.size() == 1)
115  {
121  }
122  // This one may be expensive, and hopefully we do not use it so often
123  else
124  {
130  }
131 }
132 
133 void
135 {
137 }
138 
139 void
141 {
143 }
144 
145 void
147 {
149 }
150 
151 void
153 {
155  _num_cached++;
156 
157  if (_num_cached % 20 == 0 || _fe_problem.assembly(_tid, _nl.number()).hasStaticCondensation())
159 }
160 
161 void
163 {
164 }
165 
166 void
168 {
170  "We should not be called if we have no active HDG kernels");
171  mooseAssert(false, "HDGKernels must compute the full Jacobian");
172 }
virtual void postElement(const Elem *) override
Called after the element assembly is done (including surface assembling)
virtual void addJacobianLowerD(const THREAD_ID tid) override
bool hasActiveBlockObjects(THREAD_ID tid=0) const
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:302
unsigned int number() const
Get variable number coming from libMesh.
MooseObjectTagWarehouse< IntegratedBCBase > & _integrated_bcs
Reference to BC storage structures.
unsigned int _num_cached
virtual bool checkNonlocalCouplingRequirement() const override
MooseObjectTagWarehouse< HDGKernel > & _hdg_kernels
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:139
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 Assembly & assembly(const THREAD_ID tid, const unsigned int sys_num) override
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.
MooseObjectWarehouse< HDGKernel > * _hdg_warehouse
This is the common base class for the three main kernel types implemented in MOOSE, Kernel, VectorKernel and ArrayKernel.
Definition: KernelBase.h:23
virtual void computeOnInternalFace() override
void prepareShapes(unsigned int var_num) override final
Prepare shape functions.
Definition: DGKernelBase.C:209
MooseObjectTagWarehouse< KernelBase > & _kernels
virtual void addJacobianNeighbor(const THREAD_ID tid) override
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.
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1159
NonlinearSystemBase & _nl
virtual void accumulateLower() override
Add lower-d residual/Jacobian into assembly global data.
virtual bool shouldApply() const override
Hook for turning the boundary condition on and off.
virtual unsigned int numMatrixTags() const
The total number of tags.
Definition: SubProblem.h:248
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.
InterfaceKernelBase is the base class for all InterfaceKernel type classes.
SubdomainID _subdomain
The subdomain for the current element.
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