www.mooseframework.org
ComputeJacobianThread.h
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 #pragma once
11 
12 #include "NonlinearThread.h"
14 
15 #include "libmesh/elem_range.h"
16 
17 // Forward declarations
18 class FEProblemBase;
20 class IntegratedBCBase;
21 class DGKernelBase;
23 class Kernel;
24 class FVElementalKernel;
25 
27 {
28 public:
29  ComputeJacobianThread(FEProblemBase & fe_problem, const std::set<TagID> & tags);
30 
31  // Splitting Constructor
33 
34  virtual ~ComputeJacobianThread();
35 
36  virtual void postElement(const Elem * /*elem*/) override;
37 
38  void join(const ComputeJacobianThread & /*y*/);
39 
40 protected:
41  const std::set<TagID> & _tags;
42 
43  void determineObjectWarehouses() override;
44  void accumulateNeighbor() override;
45  virtual void accumulateNeighborLower() override;
46  virtual void accumulateLower() override;
47  virtual void accumulate() override{};
48 
49  virtual void compute(ResidualObject &) override { mooseError("Not implemented"); };
50  void compute(KernelBase & kernel) override;
51  void compute(FVElementalKernel & kernel) override;
52  void compute(IntegratedBCBase & bc) override;
53  void compute(DGKernelBase & dg, const Elem * neighbor) override;
54  void compute(InterfaceKernelBase & ik) override;
55 
56  std::string objectType() const override { return "Jacobian"; }
57 };
virtual void postElement(const Elem *) override
Called after the element assembly is done (including surface assembling)
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
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
std::string objectType() const override
Return what the loops is meant to compute.
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.
Nonlinear system to be solved.
FVElemental is used for calculating residual contributions from volume integral terms of a PDE where ...
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
virtual void accumulateLower() override
Add lower-d residual/Jacobian into assembly global data.
tbb::split split
ComputeJacobianThread(FEProblemBase &fe_problem, const std::set< TagID > &tags)
This is the common base class for objects that give residual contributions.
void join(const ComputeJacobianThread &)
Definition: Kernel.h:15
Base class for deriving any boundary condition of a integrated type.
InterfaceKernelBase is the base class for all InterfaceKernel type classes.
virtual void accumulate() override
Add element residual/Jacobian into assembly global data.
virtual void compute(ResidualObject &) override
Will dispatch to computeResidual/computeJacobian/computeResidualAndJacobian based on the derived clas...