2 * SparseMatrixPetsc.hxx
4 * Created on: 03/11/2017. 2017
8 #ifndef SOURCE_DIRECTORY__BASE_INC_SparseMatrixPetsc_HXX_
9 #define SOURCE_DIRECTORY__BASE_INC_SparseMatrixPetsc_HXX_
14 #include "MEDCouplingUMesh.hxx"
20 class SparseMatrixPetsc: public GenericMatrix {
24 virtual ~SparseMatrixPetsc();
27 * constructor for a block sparse matrix
28 * @param numberOfRows : The number of rows
29 * @param numberOfColumns : The number of columns
31 SparseMatrixPetsc( int numberOfRows, int numberOfColumns ) ;
34 * constructor for a block sparse matrix with number of non zero coefficients given
35 * @param numberOfRows : The number of rows
36 * @param numberOfColumns : The number of columns
37 * @param nnz : The maximum number of nonzeros coefficients per line (or an upper bound)
39 SparseMatrixPetsc( int numberOfRows, int numberOfColumns, int nnz ) ;
42 * constructor for a sparse matrix with block structure
43 * @param blockSize : The block size
44 * @param numberOfRows : The number of rows
45 * @param numberOfColumns : The number of columns
46 * @param nnz : The maximum number of nonzeros coefficients per line (or an upper bound)
47 * @comment blockSize should always divide numberOfRows and numberOfColumns
49 SparseMatrixPetsc( int blockSize, int numberOfRows, int numberOfColumns, int nnz );
53 * @param SparseMatrixPetsc : The SparseMatrixPetsc object to be copied
55 SparseMatrixPetsc ( const SparseMatrixPetsc& sparseMatrixPetsc ) ;
58 * constructor with data
59 * @param Petsc matris : mat
61 SparseMatrixPetsc(Mat mat);
63 SparseMatrixPetsc transpose() const ;
65 double operator()( int i, int j ) const ;
67 void setValue( int i, int j, double value ) ;
68 void addValue( int i, int j, double value ) ;
70 void setValue( int i, int j, Matrix M ) ;
71 void addValue( int i, int j, Matrix M ) ;
73 void setValuesBlocked( int i, int j, Matrix M ) ;
74 void addValuesBlocked( int i, int j, Matrix M ) ;
76 SparseMatrixPetsc& operator+= (const SparseMatrixPetsc& SparseMatrixPetsc) ;
78 SparseMatrixPetsc& operator-= (const SparseMatrixPetsc& SparseMatrixPetsc) ;
80 SparseMatrixPetsc& operator*= (double value) ;
82 SparseMatrixPetsc& operator/= (double value) ;
84 SparseMatrixPetsc& operator*= (const SparseMatrixPetsc& matrix) ;
86 Vector operator* (const Vector& vector) const ;
88 const SparseMatrixPetsc& operator= ( const SparseMatrixPetsc& SparseMatrixPetsc ) ;
90 friend SparseMatrixPetsc operator+ (const SparseMatrixPetsc& SparseMatrixPetsc1, const SparseMatrixPetsc& SparseMatrixPetsc2);
92 friend SparseMatrixPetsc operator- (const SparseMatrixPetsc& SparseMatrixPetsc1, const SparseMatrixPetsc& SparseMatrixPetsc2);
94 friend SparseMatrixPetsc operator* (double value , const SparseMatrixPetsc& SparseMatrixPetsc ) ;
96 friend SparseMatrixPetsc operator* (const SparseMatrixPetsc& SparseMatrixPetsc, double value ) ;
98 friend SparseMatrixPetsc operator/ (const SparseMatrixPetsc& SparseMatrixPetsc, double value) ;
100 friend SparseMatrixPetsc operator*(const SparseMatrixPetsc& M, const SparseMatrixPetsc& N) ;
102 void viewMatrix() const ;
103 double getMatrixCoeff(int i, int j) const;
105 bool containsPetscMatrix() const;
106 Mat getPetscMatrix() const;
108 void diagonalShift(double lambda);
109 void zeroEntries();//sets the matrix coefficients to zero
111 std::vector< double > getEigenvalues(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
112 std::vector< Vector > getEigenvectors(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
113 MEDCoupling::DataArrayDouble * getEigenvectorsDataArrayDouble(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
114 std::vector< double > getSingularValues(int nsv, SVDWhich which=SVD_SMALLEST, double tol=1e-6) const;
115 double getConditionNumber(bool isSingular=false, double tol=1e-6) const;
117 bool isSymmetric(double tol=1.e-6) const ;
119 void leftDiagonalScale(Vector v);
120 void rightDiagonalScale(Vector v);
125 int _numberOfNonZeros ;//The maximum number of nonzeros coefficients per line (or an upper bound)
127 int computeSpectrum(int nev, double ** valP, double ***vecP, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
128 int computeSVD (int nsv, double ** valS, SVDWhich which=SVD_SMALLEST , double tol=1e-6) const;
130 Vector vecToVector(const Vec& vec) const ;
131 Vec vectorToVec( const Vector& myVector ) const ;
135 #endif /* SOURCE_DIRECTORY__BASE_INC_SparseMatrixPetsc_HXX_ */