DUNE-FEM (unstable)

Dune::Fem::NumpyLinearOperator< DomainFunction, RangeFunction > Struct Template Reference

NumpyLinearOperator. More...

#include <dune/fem/operator/linear/numpyoperator.hh>

Public Types

typedef DomainFunction DomainFunctionType
 type of discrete function in the operator's domain
 
typedef RangeFunction RangeFunctionType
 type of discrete function in the operator's range
 
typedef DomainFunction::RangeFieldType DomainFieldType
 field type of the operator's domain
 
typedef RangeFunction::RangeFieldType RangeFieldType
 field type of the operator's range
 

Public Member Functions

virtual void operator() (const DomainFunction &arg, RangeFunction &dest) const
 application operator More...
 
template<class DomainFunction , class RangeFunction >
void apply (const DomainFunction &arg, RangeFunction &dest) const
 apply matrix to discrete function
 
const DomainSpaceType & domainSpace () const
 get domain space (i.e. space that builds the rows)
 
const RangeSpaceType & rangeSpace () const
 get range space (i.e. space that builds the columns)
 
MatrixTypeexportMatrix () const
 get reference to storage object
 
ObjectType * newObject () const
 interface method from LocalMatrixFactory
 
LocalMatrixType localMatrix (const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity) const
 
LocalMatrixType localMatrix () const
 
LocalColumnObjectType localColumn (const DomainEntityType &domainEntity) const
 get local column
 
void compress ()
 compress matrix to a real CRS format
 
void reserve (const Stencil &stencil, bool verbose=false)
 reserve memory
 
void extractDiagonal (DiscreteFunctionType &diag) const
 
void resort ()
 resort row numbering in matrix to have ascending numbering
 
virtual void flushAssembly ()
 commit intermediate states of linear operator assembly
 
void beginAssemble ()
 Initiate the assemble of values using the LocalContribution concept. More...
 
void endAssemble ()
 Finalize the assemble of values using the LocalContribution concept. More...
 
virtual bool symmetric () const
 
virtual bool positiveDefinite () const
 
virtual void finalize ()
 finalization of operator More...
 

Protected Member Functions

MatrixTypematrix () const
 get reference to storage object, for internal use
 

Detailed Description

template<class DomainFunction, class RangeFunction>
struct Dune::Fem::NumpyLinearOperator< DomainFunction, RangeFunction >

NumpyLinearOperator.

Member Function Documentation

◆ beginAssemble()

void Dune::Fem::AssembledOperator< DomainFunction, RangeFunction >::beginAssemble ( )
inlineinherited

Initiate the assemble of values using the LocalContribution concept.

Template Parameters
AssembleOperationthe specific operation (Add, Set, ...)

◆ endAssemble()

void Dune::Fem::AssembledOperator< DomainFunction, RangeFunction >::endAssemble ( )
inlineinherited

Finalize the assemble of values using the LocalContribution concept.

Template Parameters
AssembleOperationthe specific operation (Add, Set, ...)

◆ extractDiagonal()

void Dune::Fem::SparseRowMatrixObject< DomainFunction::DiscreteFunctionSpaceType , RangeFunction::DiscreteFunctionSpaceType , SparseRowMatrix< double, size_t, pybind11::array_t< double >, pybind11::array_t< size_t > > >::extractDiagonal ( DiscreteFunctionType diag) const
inlineinherited

extract diagonal entries from matrix into discrete function this only works for square matrices

◆ finalize()

virtual void Dune::Fem::Operator< DomainFunction, RangeFunction >::finalize ( )
inlinevirtualinherited

finalization of operator

Note
The default implementation is empty.

Reimplemented in Dune::Fem::SparseRowLinearOperator< DomainFunction, RangeFunction, Matrix >.

◆ localMatrix() [1/2]

LocalMatrixType Dune::Fem::SparseRowMatrixObject< DomainFunction::DiscreteFunctionSpaceType , RangeFunction::DiscreteFunctionSpaceType , SparseRowMatrix< double, size_t, pybind11::array_t< double >, pybind11::array_t< size_t > > >::localMatrix ( ) const
inlineinherited
Deprecated:
Use TemporaryLocalMatrix in combination with {add,set,get}LocalMatrix on matrix object return local matrix object

◆ localMatrix() [2/2]

LocalMatrixType Dune::Fem::SparseRowMatrixObject< DomainFunction::DiscreteFunctionSpaceType , RangeFunction::DiscreteFunctionSpaceType , SparseRowMatrix< double, size_t, pybind11::array_t< double >, pybind11::array_t< size_t > > >::localMatrix ( const DomainEntityType &  domainEntity,
const RangeEntityType &  rangeEntity 
) const
inlineinherited
Deprecated:
Use TemporaryLocalMatrix in combination with {add,set,get}LocalMatrix on matrix object return local matrix object

◆ nonlinear()

virtual bool Dune::Fem::LinearOperator< DomainFunction, RangeFunction >::nonlinear ( ) const
inlinevirtualinherited

Return true if the Operator is nonlinear and false otherwise (default is true).

Reimplemented from Dune::Fem::Operator< DomainFunction, RangeFunction >.

◆ operator()()

template<class DomainFunction , class RangeFunction >
virtual void Dune::Fem::NumpyLinearOperator< DomainFunction, RangeFunction >::operator() ( const DomainFunction &  u,
RangeFunction &  w 
) const
inlinevirtual

application operator

Parameters
[in]uargument discrete function
[out]wdestination discrete function
Note
This method has to be implemented by all derived classes.

Implements Dune::Fem::Operator< DomainFunction, RangeFunction >.

References Dune::Fem::NumpyLinearOperator< DomainFunction, RangeFunction >::apply().

◆ positiveDefinite()

virtual bool Dune::Fem::LinearOperator< DomainFunction, RangeFunction >::positiveDefinite ( ) const
inlinevirtualinherited

Return true if the Operator is positive definite.

◆ symmetric()

virtual bool Dune::Fem::LinearOperator< DomainFunction, RangeFunction >::symmetric ( ) const
inlinevirtualinherited

Return true if the Operator is symmetric.


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