DUNE PDELab (git)

Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode > Class Template Reference

The local assembler for DUNE grids. More...

#include <dune/pdelab/gridoperator/default/localassembler.hh>

Public Types

typedef Dune::PDELab::LocalAssemblerTraits< GO > Traits
 The traits class.
 
typedef Traits::Residual::ElementType RangeField
 The local operators type for real numbers e.g. time.
 
typedef Dune::PDELab::LocalAssemblerBase< typename Traits::MatrixBackend, CU, CV > Base
 The base class of this local assembler.
 
typedef LOP LocalOperator
 The local operator.
 
typedef Dune::PDELab::LocalFunctionSpace< GFSU, Dune::PDELab::TrialSpaceTagLFSU
 

Public Member Functions

 DefaultLocalAssembler (LOP &lop, std::shared_ptr< typename GO::BorderDOFExchanger > border_dof_exchanger)
 Constructor with empty constraints.
 
 DefaultLocalAssembler (LOP &lop, const CU &cu, const CV &cv, std::shared_ptr< typename GO::BorderDOFExchanger > border_dof_exchanger)
 Constructor for non trivial constraints.
 
LOP & localOperator ()
 get a reference to the local operator
 
const LOP & localOperator () const
 get a reference to the local operator
 
template<class EG >
bool skipEntity (const EG &eg) const
 
template<class IG >
bool skipIntersection (const IG &ig) const
 
void setTime (Real time_)
 
RangeField weight () const
 Obtain the weight that was set last.
 
void setWeight (RangeField weight)
 Notifies the assembler about the current weight of assembling.
 
bool doPreProcessing () const
 Query whether to do preprocessing in the engines. More...
 
void preProcessing (bool v)
 
bool doPostProcessing () const
 Query whether to do postprocessing in the engines. More...
 
void postProcessing (bool v)
 
const GO::Traits::TrialGridFunctionSpaceConstraints & trialConstraints () const
 get the constraints on the trial grid function space
 
const GO::Traits::TestGridFunctionSpaceConstraints & testConstraints () const
 get the constraints on the test grid function space
 
std::enable_if< AlwaysTrue< X >::value &&!std::is_same< GO::Traits::TestGridFunctionSpaceConstraints, EmptyTransformation >::value >::type forwardtransform (X &x, const bool postrestrict=false) const
 Transforms a vector \( \boldsymbol{x} \) from \( V\) to \( V'\). If postrestrict == true then \(\boldsymbol{R}^T_{\boldsymbol{\tilde U}', \boldsymbol{U}'} \boldsymbol{S}_{\boldsymbol{\tilde V}}\) is applied instead of the full transformation.

 
std::enable_if< AlwaysTrue< X >::value &&!std::is_same< GO::Traits::TestGridFunctionSpaceConstraints, EmptyTransformation >::value >::type backtransform (X &x, const bool prerestrict=false) const
 Transforms a vector \( \boldsymbol{x} \) from \( V'\) to \( V\). If prerestrict == true then \(\boldsymbol{S}^T_{\boldsymbol{\tilde U}}\) is applied instead of the full transformation.

 
void preStage (Real time_, int r_)
 
LocalPatternAssemblerEnginelocalPatternAssemblerEngine (typename Traits::MatrixPattern &p)
 
LocalResidualAssemblerEnginelocalResidualAssemblerEngine (typename Traits::Residual &r, const typename Traits::Solution &x)
 
LocalJacobianAssemblerEnginelocalJacobianAssemblerEngine (typename Traits::Jacobian &a, const typename Traits::Solution &x)
 
LocalJacobianApplyAssemblerEnginelocalJacobianApplyAssemblerEngine (const typename Traits::Domain &update, typename Traits::Range &result)
 
LocalJacobianApplyAssemblerEnginelocalJacobianApplyAssemblerEngine (const typename Traits::Domain &solution, const typename Traits::Domain &update, typename Traits::Range &result)
 

Static Public Member Functions

static constexpr bool doAlphaVolume ()
 Query methods for the assembler engines. Theses methods do not belong to the assembler interface, but simplify the implementations of query methods in the engines;.
 

Protected Member Functions

void eread (const GCView &globalcontainer_view, LocalMatrix< T > &localcontainer) const
 read local stiffness matrix for entity
 
void ewrite (const LocalMatrix< T > &localcontainer, GCView &globalcontainer_view) const
 write local stiffness matrix for entity
 
void eadd (const LocalMatrix< T > &localcontainer, GCView &globalcontainer_view) const
 write local stiffness matrix for entity
 
std::enable_if< AlwaysTrue< M >::value &&!std::is_same< GO::Traits::TestGridFunctionSpaceConstraints, EmptyTransformation >::value >::type scatter_jacobian (M &local_container, GCView &global_container_view, bool symmetric_mode) const
 Scatter local jacobian to global container.
 
void etadd_symmetric (M &localcontainer, GCView &globalcontainer_view) const
 Add local matrix to global matrix, and apply Dirichlet constraints in a symmetric fashion. Apart from that, identical to etadd().
 
void add_entry (P &globalpattern, const LFSVIndices &lfsv_indices, Index i, const LFSUIndices &lfsu_indices, Index j) const
 Adding matrix entry to pattern with respect to the constraints contributions. This assembles the entries addressed by etadd(..). See the documentation there for more information about the matrix pattern.
 
void set_trivial_rows (const GFSV &gfsv, GC &globalcontainer, const C &c) const
 insert dirichlet constraints for row and assemble T^T_U in constrained rows
 
typedef DefaultLocalPatternAssemblerEngine< DefaultLocalAssemblerLocalPatternAssemblerEngine
 

Detailed Description

template<typename GO, typename LOP, bool nonoverlapping_mode = false>
class Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >

The local assembler for DUNE grids.

Template Parameters
GFSUGridFunctionSpace for ansatz functions
GFSVGridFunctionSpace for test functions
XThe solution vector representation type
RThe residual vector representation type
AThe jacobian matrix representation type
BThe matrix backend
PThe matrix pattern representation type
CUConstraints maps for the individual dofs (trial space)
CVConstraints maps for the individual dofs (test space)

Member Typedef Documentation

◆ LFSU

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
typedef Dune::PDELab::LocalFunctionSpace<GFSU, Dune::PDELab::TrialSpaceTag> Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::LFSU

The local function spaces

◆ LocalPatternAssemblerEngine

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
typedef DefaultLocalPatternAssemblerEngine<DefaultLocalAssembler> Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::LocalPatternAssemblerEngine

The local assembler engines

Member Function Documentation

◆ doPostProcessing()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
bool Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::doPostProcessing ( ) const
inline

Query whether to do postprocessing in the engines.

This method is used by the engines.

◆ doPreProcessing()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
bool Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::doPreProcessing ( ) const
inline

Query whether to do preprocessing in the engines.

This method is used by the engines.

◆ localJacobianApplyAssemblerEngine() [1/2]

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
LocalJacobianApplyAssemblerEngine& Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::localJacobianApplyAssemblerEngine ( const typename Traits::Domain solution,
const typename Traits::Domain update,
typename Traits::Range result 
)
inline

Returns a reference to the requested engine. This engine is completely configured and ready to use.

◆ localJacobianApplyAssemblerEngine() [2/2]

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
LocalJacobianApplyAssemblerEngine& Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::localJacobianApplyAssemblerEngine ( const typename Traits::Domain update,
typename Traits::Range result 
)
inline

Returns a reference to the requested engine. This engine is completely configured and ready to use.

Referenced by Dune::PDELab::GridOperator< GFSU, GFSV, LOP, MB, DF, RF, JF, CU, CV >::jacobian_apply(), and Dune::PDELab::GridOperator< GFSU, GFSV, LOP, MB, DF, RF, JF, CU, CV >::nonlinear_jacobian_apply().

◆ localJacobianAssemblerEngine()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
LocalJacobianAssemblerEngine& Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::localJacobianAssemblerEngine ( typename Traits::Jacobian a,
const typename Traits::Solution x 
)
inline

◆ localPatternAssemblerEngine()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
LocalPatternAssemblerEngine& Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::localPatternAssemblerEngine ( typename Traits::MatrixPattern p)
inline

Access methods which provid "ready to use" engines Returns a reference to the requested engine. This engine is completely configured and ready to use.

References Dune::PDELab::DefaultLocalPatternAssemblerEngine< LA >::setPattern().

Referenced by Dune::PDELab::GridOperator< GFSU, GFSV, LOP, MB, DF, RF, JF, CU, CV >::fill_pattern().

◆ localResidualAssemblerEngine()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
LocalResidualAssemblerEngine& Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::localResidualAssemblerEngine ( typename Traits::Residual r,
const typename Traits::Solution x 
)
inline

◆ postProcessing()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
void Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::postProcessing ( bool  v)
inline

This method allows to set the behavior with regard to any postprocessing within the engines. It is called by the setupGridOperators() method of the GridOperator and should not be called directly.

◆ preProcessing()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
void Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::preProcessing ( bool  v)
inline

This method allows to set the behavior with regard to any preprocessing within the engines. It is called by the setupGridOperators() method of the GridOperator and should not be called directly.

◆ preStage()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
void Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::preStage ( Real  time_,
int  r_ 
)
inline

Time stepping interface

◆ setTime()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
void Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::setTime ( Real  time_)
inline

Notifies the local assembler about the current time of assembling. Should be called before assembling if the local operator has time dependencies.

◆ skipEntity()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
template<class EG >
bool Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::skipEntity ( const EG &  eg) const
inline

Assemble on a given cell without function spaces.

Returns
If true, the assembling for this cell is assumed to be complete and the assembler continues with the next grid cell.

◆ skipIntersection()

template<typename GO , typename LOP , bool nonoverlapping_mode = false>
template<class IG >
bool Dune::PDELab::DefaultLocalAssembler< GO, LOP, nonoverlapping_mode >::skipIntersection ( const IG &  ig) const
inline

Assemble on a given intersection without function spaces.

Returns
If true, the assembling for this intersection is assumed to be complete and the assembler continues with the next grid intersection.

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.80.0 (Apr 26, 22:29, 2024)