DUNE PDELab (2.8)

Dune::PDELab::FastDGGridOperator< GFSU, GFSV, LOP, MB, DF, RF, JF, CU, CV > Class Template Reference

#include <dune/pdelab/gridoperator/fastdg.hh>

Classes

struct  SetupGridOperator
 

Public Types

typedef FastDGAssembler< GFSU, GFSV, CU, CV > Assembler
 The global assembler type.
 
using Domain = Dune::PDELab::Backend::Vector< GFSU, DF >
 The type of the domain (solution).
 
using Range = Dune::PDELab::Backend::Vector< GFSV, RF >
 The type of the range (residual).
 
using Jacobian = Dune::PDELab::Backend::Matrix< MB, Domain, Range, JF >
 The type of the jacobian.
 
typedef MB::template Pattern< Jacobian, GFSV, GFSU > Pattern
 The sparsity pattern container for the jacobian matrix.
 
typedef FastDGLocalAssembler< FastDGGridOperator, LOP, GFSU::Traits::EntitySet::Partitions::partitionIterator()==InteriorBorder_PartitionLocalAssembler
 The local assembler type.
 
typedef Dune::PDELab::GridOperatorTraits< GFSU, GFSV, MB, DF, RF, JF, CU, CV, Assembler, LocalAssemblerTraits
 The grid operator traits.
 

Public Member Functions

 FastDGGridOperator (const GFSU &gfsu_, const CU &cu_, const GFSV &gfsv_, const CV &cv_, LOP &lop_, const MB &mb_=MB())
 Constructor for non trivial constraints.
 
 FastDGGridOperator (const GFSU &gfsu_, const GFSV &gfsv_, LOP &lop_, const MB &mb_=MB())
 Constructor for empty constraints.
 
const GFSU & trialGridFunctionSpace () const
 Get the trial grid function space.
 
const GFSV & testGridFunctionSpace () const
 Get the test grid function space.
 
GFSU::Traits::SizeType globalSizeU () const
 Get dimension of space u.
 
GFSV::Traits::SizeType globalSizeV () const
 Get dimension of space v.
 
template<typename F , typename X >
void interpolate (const X &xold, F &f, X &x) const
 Interpolate the constrained dofs from given function.
 
void fill_pattern (Pattern &p) const
 Fill pattern of jacobian matrix.
 
void residual (const Domain &x, Range &r) const
 Assemble residual.
 
void jacobian (const Domain &x, Jacobian &a) const
 Assembler jacobian.
 
void jacobian_apply (const Domain &update, Range &result) const
 Apply jacobian matrix to the vector update without explicitly assembling it.
 
void jacobian_apply (const Domain &solution, const Domain &update, Range &result) const
 Apply jacobian matrix to the vector update without explicitly assembling it.
 
void nonlinear_jacobian_apply (const Domain &solution, const Domain &update, Range &result) const
 Apply jacobian matrix to the vector update without explicitly assembling it.
 
const Traits::MatrixBackendmatrixBackend () const
 Get the matrix backend for this grid operator.
 

Static Public Member Functions

template<typename GridOperatorTuple >
static void setupGridOperators (GridOperatorTuple tuple)
 

Detailed Description

template<typename GFSU, typename GFSV, typename LOP, typename MB, typename DF, typename RF, typename JF, typename CU = Dune::PDELab::EmptyTransformation, typename CV = Dune::PDELab::EmptyTransformation>
class Dune::PDELab::FastDGGridOperator< GFSU, GFSV, LOP, MB, DF, RF, JF, CU, CV >

Fast D(iscontinuous) G(alerkin) grid operator implementation.

Template Parameters
GFSUGridFunctionSpace for ansatz functions
GFSVGridFunctionSpace for test functions
MBThe matrix backend to be used for representation of the jacobian
DFThe domain field type of the operator
RFThe range field type of the operator
JFThe jacobian field type
CUConstraints maps for the individual dofs (trial space)
CVConstraints maps for the individual dofs (test space)

Member Function Documentation

◆ setupGridOperators()

template<typename GFSU , typename GFSV , typename LOP , typename MB , typename DF , typename RF , typename JF , typename CU = Dune::PDELab::EmptyTransformation, typename CV = Dune::PDELab::EmptyTransformation>
template<typename GridOperatorTuple >
static void Dune::PDELab::FastDGGridOperator< GFSU, GFSV, LOP, MB, DF, RF, JF, CU, CV >::setupGridOperators ( GridOperatorTuple  tuple)
inlinestatic

Method to set up a number of grid operators which are used in a joint assembling. It is assumed that all operators are specializations of the same template type

References Dune::Hybrid::forEach().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)