https://mooseframework.inl.gov
ComputeResidualAndJacobianThread.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 
11 #include "NonlinearSystem.h"
12 #include "Problem.h"
13 #include "FEProblem.h"
14 #include "KernelBase.h"
15 #include "IntegratedBCBase.h"
16 #include "DGKernelBase.h"
17 #include "InterfaceKernelBase.h"
18 #include "Material.h"
19 #include "TimeKernel.h"
20 #include "SwapBackSentinel.h"
21 #include "FVTimeKernel.h"
22 #include "HDGKernel.h"
23 
24 #include "libmesh/threads.h"
25 
27  FEProblemBase & fe_problem,
28  const std::set<TagID> & vector_tags,
29  const std::set<TagID> & matrix_tags)
30  : NonlinearThread(fe_problem), _vector_tags(vector_tags), _matrix_tags(matrix_tags)
31 {
32 }
33 
34 // Splitting Constructor
37  : NonlinearThread(x, split), _vector_tags(x._vector_tags), _matrix_tags(x._matrix_tags)
38 {
39 }
40 
42 
43 void
45 {
47 }
48 
49 void
51 {
54 }
55 
56 void
58 {
62 }
63 
64 void
66 {
69 }
70 
71 void
73 {
76  _num_cached++;
77 
78  if (_num_cached % 20 == 0)
79  {
80  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
83  }
84 }
85 
86 void
88 {
89 }
90 
91 void
93 {
94  // We need to filter out vector tags that don't belong to the current nonlinear system
95  const auto & residual_vector_tags = _fe_problem.getVectorTags(Moose::VECTOR_TAG_RESIDUAL);
96 
97  // We would only like to consider the tags that belong to the current system
98  std::set<TagID> filtered_residual_tags;
99  _fe_problem.selectVectorTagsFromSystem(_nl, residual_vector_tags, filtered_residual_tags);
100 
101  if (_vector_tags.size() && _vector_tags.size() != filtered_residual_tags.size())
102  mooseError("Can only currently compute the residual and Jacobian together if we are computing "
103  "the full suite of residual tags");
104 
105  if (_matrix_tags.size() && _matrix_tags.size() != _fe_problem.numMatrixTags())
106  mooseError("Can only currently compute the residual and Jacobian together if we are computing "
107  "the full suite of Jacobian tags");
108 
114 
115  if (_fe_problem.haveFV())
116  {
117  _fv_kernels.clear();
119  .query()
120  .template condition<AttribSysNum>(_nl.number())
121  .template condition<AttribSystem>("FVElementalKernel")
122  .template condition<AttribSubdomains>(_subdomain)
123  .template condition<AttribThread>(_tid)
125  }
126 }
127 
128 void
130 {
132  "We should not be called if we have no active HDG kernels");
133  for (const auto & hdg_kernel : _hdg_warehouse->getActiveBlockObjects(_subdomain, _tid))
134  if (hdg_kernel->hasBlocks(_subdomain))
135  hdg_kernel->computeResidualAndJacobianOnSide();
136 }
void accumulateNeighborLower() override
Add neighbor and lower residual/Jacobian into assembly global data.
virtual void computeResidualAndJacobian()
Compute this object&#39;s contribution to the residual and Jacobian simultaneously.
void determineObjectWarehouses() override
Determine the objects we will actually compute based on vector/matrix tag information.
virtual void addResidualLower(const THREAD_ID tid) override
virtual void addJacobianLowerD(const THREAD_ID tid) override
bool hasActiveBlockObjects(THREAD_ID tid=0) const
const std::map< SubdomainID, std::vector< std::shared_ptr< T > > > & getActiveBlockObjects(THREAD_ID tid=0) const
MooseObjectWarehouse< InterfaceKernelBase > * _ik_warehouse
virtual bool haveFV() const override
returns true if this problem includes/needs finite volume functionality.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
MooseObjectTagWarehouse< IntegratedBCBase > & _integrated_bcs
Reference to BC storage structures.
unsigned int _num_cached
std::vector< T * > & queryInto(std::vector< T *> &results, Args &&... args)
queryInto executes the query and stores the results in the given vector.
Definition: TheWarehouse.h:311
MooseObjectTagWarehouse< HDGKernel > & _hdg_kernels
static void selectVectorTagsFromSystem(const SystemBase &system, const std::vector< VectorTag > &input_vector_tags, std::set< TagID > &selected_tags)
Select the vector tags which belong to a specific system.
Definition: SubProblem.C:289
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
void accumulate() override
Add element residual/Jacobian into assembly global data.
MooseObjectWarehouse< HDGKernel > * _hdg_warehouse
void accumulateNeighbor() override
Add neighbor residual/Jacobian into assembly global data.
TheWarehouse & theWarehouse() const
MooseObjectTagWarehouse< KernelBase > & _kernels
virtual void addJacobianNeighbor(const THREAD_ID tid) override
std::vector< VectorTag > getVectorTags(const std::set< TagID > &tag_ids) const
Definition: SubProblem.C:172
virtual void cacheResidual(const THREAD_ID tid) override
void join(const ComputeResidualAndJacobianThread &)
virtual void addResidualNeighbor(const THREAD_ID tid) override
void compute(ResidualObject &ro) override
Will dispatch to computeResidual/computeJacobian/computeResidualAndJacobian based on the derived clas...
unsigned int number() const
Gets the number of this system.
Definition: SystemBase.C:1159
NonlinearSystemBase & _nl
void accumulateLower() override
Add lower-d residual/Jacobian into assembly global data.
const std::set< TagID > & _vector_tags
the tags denoting the vectors we want our residual objects to fill
virtual unsigned int numMatrixTags() const
The total number of tags.
Definition: SubProblem.h:248
tbb::split split
MooseObjectWarehouse< IntegratedBCBase > * _ibc_warehouse
This is the common base class for objects that give residual contributions.
Query query()
query creates and returns an initialized a query object for querying objects from the warehouse...
Definition: TheWarehouse.h:466
MooseObjectTagWarehouse< InterfaceKernelBase > & _interface_kernels
Reference to interface kernel storage structure.
ComputeResidualAndJacobianThread(FEProblemBase &fe_problem, const std::set< TagID > &vector_tags, const std::set< TagID > &matrix_tags)
const std::set< TagID > & _matrix_tags
the tags denoting the matrices we want our residual objects to fill
SubdomainID _subdomain
The subdomain for the current element.
MooseObjectWarehouse< KernelBase > * _tag_kernels
MooseObjectWarehouse< DGKernelBase > * _dg_warehouse
std::vector< FVElementalKernel * > _fv_kernels
Current subdomain FVElementalKernels.
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 addCachedResidual(const THREAD_ID tid) override
virtual void addCachedJacobian(const THREAD_ID tid) override