DUNE PDELab (git)

Dune::PDELab::NonoverlappingOperator< GFS, M, X, Y > Class Template Reference

Operator for the non-overlapping parallel case. More...

#include <dune/pdelab/backend/istl/novlpistlsolverbackend.hh>

Public Types

using matrix_type = Backend::Native< M >
 export type of matrix
 
using domain_type = Backend::Native< X >
 export type of vectors the matrix is applied to
 
using range_type = Backend::Native< Y >
 export type of result vectors
 
typedef X::field_type field_type
 export type of the entries for x
 

Public Member Functions

 NonoverlappingOperator (const GFS &gfs_, const M &A)
 Construct a non-overlapping operator. More...
 
virtual void apply (const X &x, Y &y) const override
 apply operator More...
 
virtual void applyscaleadd (field_type alpha, const X &x, Y &y) const override
 apply operator to x, scale and add: \( y = y + \alpha A(x) \) More...
 
SolverCategory::Category category () const override
 Category of the linear operator (see SolverCategory::Category)
 
virtual const M & getmat () const override
 extract the matrix
 

Detailed Description

template<typename GFS, typename M, typename X, typename Y>
class Dune::PDELab::NonoverlappingOperator< GFS, M, X, Y >

Operator for the non-overlapping parallel case.

Calculate \(y:=Ax\).

Template Parameters
GFSThe GridFunctionSpace the vectors apply to.
MType of the matrix. Should be one of the ISTL matrix types.
XType of the vectors the matrix is applied to.
YType of the result vectors.

Constructor & Destructor Documentation

◆ NonoverlappingOperator()

template<typename GFS , typename M , typename X , typename Y >
Dune::PDELab::NonoverlappingOperator< GFS, M, X, Y >::NonoverlappingOperator ( const GFS &  gfs_,
const M &  A 
)
inline

Construct a non-overlapping operator.

Parameters
gfs_GridFunctionsSpace for the vectors.
AMatrix for this operator. This should be the locally assembled matrix.
Note
The constructed object stores references to all the objects given as parameters here. They should be valid for as long as the constructed object is used. They are not needed to destruct the constructed object.

Member Function Documentation

◆ apply()

template<typename GFS , typename M , typename X , typename Y >
virtual void Dune::PDELab::NonoverlappingOperator< GFS, M, X, Y >::apply ( const X &  x,
Y &  y 
) const
inlineoverridevirtual

apply operator

Compute \(y:=A(x)\) on this process, then make y consistent (sum up corresponding entries of y on the different processes and store the result back in y on each process).

Implements Dune::LinearOperator< X, Y >.

References Dune::ForwardCommunication, and Dune::InteriorBorder_InteriorBorder_Interface.

◆ applyscaleadd()

template<typename GFS , typename M , typename X , typename Y >
virtual void Dune::PDELab::NonoverlappingOperator< GFS, M, X, Y >::applyscaleadd ( field_type  alpha,
const X &  x,
Y &  y 
) const
inlineoverridevirtual

apply operator to x, scale and add: \( y = y + \alpha A(x) \)

Compute \(y:=\alpha A(x)\) on this process, then make y consistent (sum up corresponding entries of y on the different processes and store the result back in y on each process).

Implements Dune::LinearOperator< X, Y >.

References Dune::ForwardCommunication, and Dune::InteriorBorder_InteriorBorder_Interface.


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 (May 16, 22:29, 2024)