Interface to the sparse matrices from PETSc. More...

#include <PETSc.hh>

Inheritance diagram for concepts::PETScMat:
concepts::Matrix< Real >

Public Types

typedef Cmplxtype< Real >::type c_type
 Complex type of data type. More...
 
typedef _Matrix_iterator< Real, const Real &, const Real * > const_iterator
 
typedef std::conditional< std::is_same< typename Realtype< Real >::type, Real >::value, typename Realtype< Real >::type, typename Cmplxtype< Real >::type >::type d_type
 Data type, depending if F is real or complex. More...
 
typedef _Matrix_iterator< Real, Real &, Real * > iterator
 
typedef Realtype< Real >::type r_type
 Real type of data type. More...
 
typedef Real value_type
 

Public Member Functions

virtual void add (const uint i, const uint j, const Real value, const bool use_threshold=false, const Real threshold_value=1e-8)
 Addition operator Add the value to the entry (i,j) if not zero (for use_threshold = false), or if the absolute avlue is bigger than the threshold_value (for use_threshold = true) More...
 
iterator begin (uint r=0)
 Iterator over the elements, standing at position (r,0). More...
 
const_iterator begin (uint r=0) const
 Constant iterator over the elements, standing at position (r,0) More...
 
iterator end ()
 Iterator, standing behind last element. More...
 
const_iterator end () const
 Constant iterator, standing behind last element. More...
 
const uint nofCols () const
 Number of columns. More...
 
const uint nofRows () const
 Number of rows. More...
 
 operator Mat ()
 Conversion operator. Returns the PETSc matrix. More...
 
virtual void operator() (const Function< c_type > &fncY, Function< c_type > &fncX)=0
 
virtual void operator() (const Function< double > &fncY, Function< double > &fncX)
 
virtual void operator() (const Function< r_type > &fncY, Function< Real > &fncX)=0
 Computes fncX = A(fncY) where A is this matrix. More...
 
virtual void operator() (const Function< std::complex< double > > &fncY, Function< std::complex< double > > &fncX)
 
virtual Realoperator() (const uint i, const uint j)
 Should return and allow access to entry with indices i and j. More...
 
virtual Real operator() (const uint i, const uint j) const
 Returns entry with indices i and j. More...
 
void operator() (const Vector< double > &fncY, Vector< double > &fncX)
 
virtual bool operator== (const Matrix< Real > &otherMat) const
 
 PETScMat (const Space< Real > &spc, BilinearForm< Real, Real > &bf)
 Constructor. More...
 
 PETScMat (const Space< Real > &spcX, const Space< Real > &spcY)
 Constructor. More...
 
 PETScMat (const SparseMatrix< Real > &matrix)
 Constructor. Copies the matrix from matrix. More...
 
virtual void set (const uint i, const uint j, const Real value, const bool use_threshold=false, const Real threshold_value=1e-8)
 Affectation operator Affet the value to the entry (i,j) if not zero (for use_threshold = false), or if the absolute value is bigger than the threshold_value (for use_threshold = true) More...
 
virtual const Space< double > & spaceX () const
 
virtual const Space< double > & spaceY () const
 
void storeMatlab (const char *name) const
 
virtual void transpMult (const Vector< c_type > &fncY, Vector< c_type > &fncX)=0
 
virtual void transpMult (const Vector< double > &fncY, Vector< double > &fncX)
 
virtual void transpMult (const Vector< r_type > &fncY, Vector< Real > &fncX)=0
 Computes fncX = AT fncY where A is this matrix. More...
 
virtual void transpMult (const Vector< std::complex< double > > &fncY, Vector< std::complex< double > > &fncX)
 
 ~PETScMat ()
 

Static Public Member Functions

static void assembly (Matrix< Real > &dest, const BilinearForm< Real, G > &bf, const ElementPairList< G > &pairs)
 Assembly operator for dest using the bilinear form bf. More...
 
static void assembly (Matrix< Real > &dest, const Sequence< ElementWithCell< G > * > seq, const BilinearForm< Real, G > &bf, const Real threshold=0.0)
 Assembly operator for dest using the bilinear form bf. More...
 
static void assembly (Matrix< Real > &dest, const Sequence< ElementWithCell< G > * > seqX, const Space< G > &spcY, const BilinearForm< Real, G > &bf, const Real threshold=0.0)
 Assembly operator for dest using the bilinear form bf. More...
 
static void assembly (Matrix< Real > &dest, const Space< G > &spc, const BilinearForm< Real, G > &bf, const Real threshold=0.0)
 Assembly operator for dest using the bilinear form bf. More...
 
static void assembly (Matrix< Real > &dest, const Space< G > &spc, const Sequence< bool > &seq, const BilinearForm< Real, G > &bf, const Real threshold=0.0)
 Assembly operator for dest using the bilinear form bf on space spc. More...
 
static void assembly (Matrix< Real > &dest, const Space< G > &spcX, const Space< G > &spcY, const BilinearForm< Real, G > &bf, const Real threshold=0.0, const bool single=false)
 Assembly operator for dest using the bilinear form bf. More...
 
static void assembly (Matrix< Real > &dest, const Space< G > &spcX, const Space< G > &spcY, const Sequence< bool > &seq, const BilinearForm< Real, G > &bf, const Real threshold=0.0, const bool single=false)
 Assembly operator for dest using the bilinear form bf. More...
 
static void assembly (Matrix< Real > &dest, Scan< Element< G > > *sc, const BilinearForm< Real, G > &bf, const Real threshold=0.0)
 Assembly operator for dest using the bilinear form bf. More...
 

Protected Member Functions

std::ostream & info (std::ostream &os) const
 

Private Attributes

Mat A_
 PETSc matrix. More...
 
uint nX_
 Dimension of image space (spcX_) More...
 
uint nY_
 Dimension of source space (spcY_) More...
 
const Space< Real > & spcX_
 Image space. More...
 
const Space< Real > & spcY_
 Source space. More...
 
Vec x_
 PETSc vector handle. More...
 
Vec y_
 

Timing Interface

These functions are used to get timings from class internal computations.

The values are stored in a user defined concepts::InOutParameters structure in different arrays (see setTimings). These arrays can be grouped into a table for easier postprocessing with

table.addMap(concepts::ResultsTable::DOUBLE, "jacobian", output);
table.addMap(concepts::ResultsTable::DOUBLE, "whole_sumfact", output);
std::ofstream ofs("table.gnuplot");
ofs << std::setprecision(20);
static InOutParameterstimings_
 Place to store timing values. More...
 
static uint timeCntr_
 Counter for timing table. More...
 
static void setTimings (InOutParameters *timings)
 Sets the class to store the timing values in. More...
 
static bool timings ()
 Returns true if the class is able to do timings. More...
 

Detailed Description

Interface to the sparse matrices from PETSc.

Very usefull if you intend to use PETSc to solve your systems since this makes it possible to use the preconditioners of PETSc.

It is not possible tough, to combine different PETSc matrices. If you want to compute a linear combination of bilinear forms, then insert the linear combination of the bilinear forms into the constructor.

PETSc is a suite of data structures and routines for the scalable (parallel) solution of scientific applications modeled by partial differential equations. It employs the MPI standard for all message-passing communication.

Our interface to PETSc uses only its serial capablities.

See also
PETSc solvers in Concepts
documentation of PETSc on matrices
Satish Balay, Kris Buschelman, William D. Gropp, Dinesh Kaushik, Lois Curfman McInnes, and Barry F. Smith, PETSc home page, 2001.
Satish Balay, William D. Gropp, Lois Curfman McInnes, and Barry F. Smith. Efficient Management of Parallelism in Object Oriented Numerical Software Libraries. In E. Arge, A. M. Bruaset, and H. P. Langtangen, editors, Modern Software Tools in Scientific Computing, pages 163-202. Birkhauser Press, 1997.
Satish Balay, William D. Gropp, Lois Curfman McInnes, and Barry F. Smith. PETSc Users Manual. Technical Report ANL-95/11 - Revision 2.1.0, Argonne National Laboratory, 2001.
Author
Philipp Frauenfelder, 2001
Examples
hpFEM2d.cc.

Definition at line 148 of file PETSc.hh.

Member Typedef Documentation

◆ c_type

typedef Cmplxtype<Real >::type concepts::Matrix< Real >::c_type
inherited

Complex type of data type.

Definition at line 45 of file matrix.hh.

◆ const_iterator

typedef _Matrix_iterator<Real , const Real &, const Real *> concepts::Matrix< Real >::const_iterator
inherited

Definition at line 51 of file matrix.hh.

◆ d_type

typedef std::conditional<std::is_same<typename Realtype<Real >::type, Real >::value , typename Realtype<Real >::type, typename Cmplxtype<Real >::type >::type concepts::Matrix< Real >::d_type
inherited

Data type, depending if F is real or complex.

Definition at line 48 of file matrix.hh.

◆ iterator

typedef _Matrix_iterator<Real , Real &, Real *> concepts::Matrix< Real >::iterator
inherited

Definition at line 50 of file matrix.hh.

◆ r_type

typedef Realtype<Real >::type concepts::Matrix< Real >::r_type
inherited

Real type of data type.

Definition at line 43 of file matrix.hh.

◆ value_type

typedef Real concepts::Matrix< Real >::value_type
inherited

Definition at line 41 of file matrix.hh.

Constructor & Destructor Documentation

◆ PETScMat() [1/3]

concepts::PETScMat::PETScMat ( const Space< Real > &  spcX,
const Space< Real > &  spcY 
)

Constructor.

Creates an empty matrix.

◆ PETScMat() [2/3]

concepts::PETScMat::PETScMat ( const Space< Real > &  spc,
BilinearForm< Real, Real > &  bf 
)

Constructor.

Computes the global matrix by assembling the element matrices.

◆ PETScMat() [3/3]

concepts::PETScMat::PETScMat ( const SparseMatrix< Real > &  matrix)

Constructor. Copies the matrix from matrix.

◆ ~PETScMat()

concepts::PETScMat::~PETScMat ( )

Member Function Documentation

◆ add()

virtual void concepts::Matrix< Real >::add ( const uint  i,
const uint  j,
const Real  value,
const bool  use_threshold = false,
const Real  threshold_value = 1e-8 
)
virtualinherited

Addition operator Add the value to the entry (i,j) if not zero (for use_threshold = false), or if the absolute avlue is bigger than the threshold_value (for use_threshold = true)

◆ assembly() [1/8]

static void concepts::Matrix< Real >::assembly ( Matrix< Real > &  dest,
const BilinearForm< Real , G > &  bf,
const ElementPairList< G > &  pairs 
)
staticinherited

Assembly operator for dest using the bilinear form bf.

This assembly operator uses the element pairs taken from pairs. For every two elements found in a ElementPair in pairs, the bilinear form is evaluated and the result assembled into dest.

◆ assembly() [2/8]

static void concepts::Matrix< Real >::assembly ( Matrix< Real > &  dest,
const Sequence< ElementWithCell< G > * >  seq,
const BilinearForm< Real , G > &  bf,
const Real  threshold = 0.0 
)
staticinherited

Assembly operator for dest using the bilinear form bf.

This assembly operator does not compute element matrices for two different elements. The elements are taken from the element sequence seq.

◆ assembly() [3/8]

static void concepts::Matrix< Real >::assembly ( Matrix< Real > &  dest,
const Sequence< ElementWithCell< G > * >  seqX,
const Space< G > &  spcY,
const BilinearForm< Real , G > &  bf,
const Real  threshold = 0.0 
)
staticinherited

Assembly operator for dest using the bilinear form bf.

This assembly operator computes also the element matrices for two different elements. The elements are taken from the element sequence seqX for the test space, and the trial space spcY is fully taken into account.

◆ assembly() [4/8]

static void concepts::Matrix< Real >::assembly ( Matrix< Real > &  dest,
const Space< G > &  spc,
const BilinearForm< Real , G > &  bf,
const Real  threshold = 0.0 
)
staticinherited

Assembly operator for dest using the bilinear form bf.

This assembly operator does not compute element matrices for two different elements. The elements are taken from the space spc.

◆ assembly() [5/8]

static void concepts::Matrix< Real >::assembly ( Matrix< Real > &  dest,
const Space< G > &  spc,
const Sequence< bool > &  seq,
const BilinearForm< Real , G > &  bf,
const Real  threshold = 0.0 
)
staticinherited

Assembly operator for dest using the bilinear form bf on space spc.

This assembly operator does not compute element matrices for two different elements. The elements are computing on the cells that are flagged by seq.

◆ assembly() [6/8]

static void concepts::Matrix< Real >::assembly ( Matrix< Real > &  dest,
const Space< G > &  spcX,
const Space< G > &  spcY,
const BilinearForm< Real , G > &  bf,
const Real  threshold = 0.0,
const bool  single = false 
)
staticinherited

Assembly operator for dest using the bilinear form bf.

This assembly operator computes also the element matrices for two different elements (coming from test space spcX and trial space spcY).

◆ assembly() [7/8]

static void concepts::Matrix< Real >::assembly ( Matrix< Real > &  dest,
const Space< G > &  spcX,
const Space< G > &  spcY,
const Sequence< bool > &  seq,
const BilinearForm< Real , G > &  bf,
const Real  threshold = 0.0,
const bool  single = false 
)
staticinherited

Assembly operator for dest using the bilinear form bf.

This assembly operator computes also the element matrices for two different elements (coming from test space spcX and trial space spcY). The elements are computing on the cells that are flagged by seq.

◆ assembly() [8/8]

static void concepts::Matrix< Real >::assembly ( Matrix< Real > &  dest,
Scan< Element< G > > *  sc,
const BilinearForm< Real , G > &  bf,
const Real  threshold = 0.0 
)
staticinherited

Assembly operator for dest using the bilinear form bf.

This assembly operator does not compute element matrices for two different elements. The elements are taken from the space scanner sc.

◆ begin() [1/2]

iterator concepts::Matrix< Real >::begin ( uint  r = 0)
inlineinherited

Iterator over the elements, standing at position (r,0).

Might be implemented differently for derived classes.

Definition at line 64 of file matrix.hh.

◆ begin() [2/2]

const_iterator concepts::Matrix< Real >::begin ( uint  r = 0) const
inlineinherited

Constant iterator over the elements, standing at position (r,0)

Might be implemented differently for derived classes.

Definition at line 71 of file matrix.hh.

◆ end() [1/2]

iterator concepts::Matrix< Real >::end
inlineinherited

Iterator, standing behind last element.

Definition at line 66 of file matrix.hh.

◆ end() [2/2]

const_iterator concepts::Matrix< Real >::end
inlineinherited

Constant iterator, standing behind last element.

Definition at line 73 of file matrix.hh.

◆ info()

std::ostream& concepts::PETScMat::info ( std::ostream &  os) const
protected

◆ nofCols()

const uint concepts::Matrix< Real >::nofCols
inlineinherited

Number of columns.

Definition at line 58 of file matrix.hh.

◆ nofRows()

const uint concepts::Matrix< Real >::nofRows
inlineinherited

Number of rows.

Definition at line 56 of file matrix.hh.

◆ operator Mat()

concepts::PETScMat::operator Mat ( )
inline

Conversion operator. Returns the PETSc matrix.

Definition at line 187 of file PETSc.hh.

◆ operator()() [1/7]

virtual void concepts::Matrix< Real >::operator() ( const Function< c_type > &  fncY,
Function< c_type > &  fncX 
)
pure virtualinherited

◆ operator()() [2/7]

virtual void concepts::PETScMat::operator() ( const Function< double > &  fncY,
Function< double > &  fncX 
)
virtual

◆ operator()() [3/7]

virtual void concepts::Matrix< Real >::operator() ( const Function< r_type > &  fncY,
Function< Real > &  fncX 
)
pure virtualinherited

Computes fncX = A(fncY) where A is this matrix.

Implemented in concepts::Permutation< Real >.

◆ operator()() [4/7]

virtual void concepts::PETScMat::operator() ( const Function< std::complex< double > > &  fncY,
Function< std::complex< double > > &  fncX 
)
virtual

◆ operator()() [5/7]

virtual Real& concepts::PETScMat::operator() ( const uint  i,
const uint  j 
)
virtual

Should return and allow access to entry with indices i and j.

This is not possible with the interface of PETSc. Therefore, this method throws an exception.

Exceptions
MissingFeature

Implements concepts::Matrix< Real >.

◆ operator()() [6/7]

virtual Real concepts::PETScMat::operator() ( const uint  i,
const uint  j 
) const
virtual

Returns entry with indices i and j.

Implements concepts::Matrix< Real >.

◆ operator()() [7/7]

void concepts::PETScMat::operator() ( const Vector< double > &  fncY,
Vector< double > &  fncX 
)

◆ operator==()

virtual bool concepts::Matrix< Real >::operator== ( const Matrix< Real > &  otherMat) const
inlinevirtualinherited

Definition at line 179 of file matrix.hh.

◆ set()

virtual void concepts::Matrix< Real >::set ( const uint  i,
const uint  j,
const Real  value,
const bool  use_threshold = false,
const Real  threshold_value = 1e-8 
)
virtualinherited

Affectation operator Affet the value to the entry (i,j) if not zero (for use_threshold = false), or if the absolute value is bigger than the threshold_value (for use_threshold = true)

◆ setTimings()

static void concepts::Matrix< Real >::setTimings ( InOutParameters timings)
staticinherited

Sets the class to store the timing values in.

Additionally, the timeCntr_ is reset to 0. This counter is used to fill in the values into the arrays listed below in subsequent calls. The following timings are taken and stored in timings:

  • evaluation of bilinear form in bilinear_form
  • application of T matrix in tmatrix_apply
  • assembling into global matrix in global_assembly

◆ spaceX()

virtual const Space<double>& concepts::PETScMat::spaceX ( ) const
inlinevirtual

Definition at line 183 of file PETSc.hh.

◆ spaceY()

virtual const Space<double>& concepts::PETScMat::spaceY ( ) const
inlinevirtual

Definition at line 184 of file PETSc.hh.

◆ storeMatlab()

void concepts::PETScMat::storeMatlab ( const char *  name) const

◆ timings()

static bool concepts::Matrix< Real >::timings
staticinherited

Returns true if the class is able to do timings.

The ability to do timings depends on a compiler switch in matrix.cc file.

◆ transpMult() [1/4]

virtual void concepts::Matrix< Real >::transpMult ( const Vector< c_type > &  fncY,
Vector< c_type > &  fncX 
)
pure virtualinherited

◆ transpMult() [2/4]

virtual void concepts::PETScMat::transpMult ( const Vector< double > &  fncY,
Vector< double > &  fncX 
)
virtual

◆ transpMult() [3/4]

virtual void concepts::Matrix< Real >::transpMult ( const Vector< r_type > &  fncY,
Vector< Real > &  fncX 
)
pure virtualinherited

Computes fncX = AT fncY where A is this matrix.

Implemented in concepts::Permutation< Real >.

◆ transpMult() [4/4]

virtual void concepts::PETScMat::transpMult ( const Vector< std::complex< double > > &  fncY,
Vector< std::complex< double > > &  fncX 
)
virtual

Member Data Documentation

◆ A_

Mat concepts::PETScMat::A_
private

PETSc matrix.

Definition at line 209 of file PETSc.hh.

◆ nX_

uint concepts::PETScMat::nX_
private

Dimension of image space (spcX_)

Definition at line 197 of file PETSc.hh.

◆ nY_

uint concepts::PETScMat::nY_
private

Dimension of source space (spcY_)

Definition at line 203 of file PETSc.hh.

◆ spcX_

const Space<Real>& concepts::PETScMat::spcX_
private

Image space.

Definition at line 194 of file PETSc.hh.

◆ spcY_

const Space<Real>& concepts::PETScMat::spcY_
private

Source space.

Definition at line 200 of file PETSc.hh.

◆ timeCntr_

uint concepts::Matrix< Real >::timeCntr_
staticprivateinherited

Counter for timing table.

Definition at line 236 of file matrix.hh.

◆ timings_

InOutParameters* concepts::Matrix< Real >::timings_
staticprivateinherited

Place to store timing values.

Definition at line 234 of file matrix.hh.

◆ x_

Vec concepts::PETScMat::x_
private

PETSc vector handle.

Definition at line 206 of file PETSc.hh.

◆ y_

Vec concepts::PETScMat::y_
private

Definition at line 206 of file PETSc.hh.


The documentation for this class was generated from the following file:
void addMap(enum mapTypes type, const char *name, const InOutParameters &holder)
Organizes the results in the hashes from InOutParameters in a nice table.
Definition: resultsTable.hh:23
void print(std::ostream &os) const
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich