53 mtrx(NULL), symmFlag(false), leqs(0), geqs(0), blocksize(1), di(0), kspInit(false), newValues(true), localIS(NULL), globalIS(NULL) { }
61 MatDestroy(& this->
mtrx);
63 KSPDestroy(& this->
ksp);
76 MatDuplicate( this->
mtrx, MAT_COPY_VALUES, & ( answer->
mtrx ) );
81 answer->
di = this->
di;
104 VecDuplicate(globX, & globY);
106 MatMult(this->
mtrx, globX, globY);
119 VecCreate(PETSC_COMM_SELF, & globY);
120 VecSetType(globY, VECSEQ);
121 VecSetSizes(globY, PETSC_DECIDE, this->
nRows);
123 MatMult(this->
mtrx, globX, globY);
125 VecGetArray(globY, & ptr);
127 for (
int i = 0; i < this->
nRows; i++ ) {
128 answer(i) = ptr [ i ];
131 VecRestoreArray(globY, & ptr);
150 VecCreate(PETSC_COMM_SELF, & globY);
151 VecSetType(globY, VECSEQ);
152 VecSetSizes(globY, PETSC_DECIDE, this->
nColumns);
154 MatMultTranspose(this->
mtrx, globX, globY);
156 VecGetArray(globY, & ptr);
158 for (
int i = 0; i < this->
nColumns; i++ ) {
159 answer(i) = ptr [ i ];
162 VecRestoreArray(globY, & ptr);
189 VecCreate(PETSC_COMM_SELF, & globY);
190 VecSetType(globY, VECSEQ);
191 VecSetSizes(globY, PETSC_DECIDE, nr);
195 for (
int k = 0; k < nc; ++k ) {
197 VecCreateSeqWithArray(PETSC_COMM_SELF, 1, nrB, colVals.
givePointer(), & globX);
198 MatMult(this->
mtrx, globX, globY);
199 VecGetArray(globY, & resultsPtr);
200 for (
int i = 0; i < nr; ++i ) {
201 answer(i, k) = resultsPtr [ i ];
203 VecRestoreArray(globY, & resultsPtr);
211 MatMatMult(this->
mtrx, globB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, & globC);
213 for (
int r = 0; r < nr; r++ ) {
214 MatGetRow(globC, r, NULL, NULL, & vals);
215 for (
int i = 0, i2 = r; i < nc; i++, i2 += nr ) {
216 aptr [ i2 ] = vals [ i ];
218 MatRestoreRow(globC, r, NULL, NULL, & vals);
249 MatMatMult(globB, this->
mtrx, MAT_INITIAL_MATRIX, PETSC_DEFAULT, & globC);
251 for (
int r = 0; r < nc; r++ ) {
252 MatGetRow(globC, r, NULL, NULL, & vals);
253 for (
int i = 0; i < nr; i++ ) {
254 * aptr++ = vals [ i ];
256 MatRestoreRow(globC, r, NULL, NULL, & vals);
266 MatScale(this->
mtrx, x);
281 OOFEM_WARNING(
"Calling function that has not been tested (remove this message when it is verified)");
289 MatDiagonalSet(this->
mtrx, globM, ADD_VALUES);
309 OOFEM_ERROR(
"Parallel is not supported in manual matrix creation");
320 std :: vector< IntArray > rows_upper(
nRows), rows_lower(
nRows);
327 rows_upper [ ii - 1 ].insertSortedOnce(jj - 1);
329 rows_lower [ ii - 1 ].insertSortedOnce(jj - 1);
337 for (
int i = 0; i <
leqs; i++ ) {
338 d_nnz(i) = rows_upper [ i ].giveSize() + rows_lower [ i ].giveSize();
343 MatCreate(PETSC_COMM_SELF, &
mtrx);
345 MatSetFromOptions(
mtrx);
348 MatSetType(
mtrx, MATDENSE);
350 MatSetType(
mtrx, MATSEQAIJ);
355 MatSeqAIJSetPreallocation(
mtrx, 0, d_nnz.givePointer() );
357 MatSetOption(
mtrx, MAT_ROW_ORIENTED, PETSC_FALSE);
358 MatSetOption(
mtrx, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
399 std :: vector< IntArray > rows_upper(
nRows), rows_lower(
nRows);
402 elem->giveLocationArray(r_loc, r_s);
403 elem->giveLocationArray(c_loc, c_s);
404 for (
int ii : r_loc ) {
406 for (
int jj : c_loc ) {
409 rows_upper [ ii - 1 ].insertSortedOnce(jj - 1, c_loc.giveSize() / 2);
411 rows_lower [ ii - 1 ].insertSortedOnce(jj - 1, c_loc.giveSize() / 2);
419 std :: vector< IntArray >r_locs, c_locs;
420 for (
auto &gbc : domain->
giveBcs() ) {
424 for ( std :: size_t k = 0; k < r_locs.size(); k++ ) {
427 for (
int ii : krloc ) {
429 for (
int jj : kcloc ) {
434 rows_lower [ ii - 1 ].insertSortedOnce(jj - 1, kcloc.giveSize() / 2);
445 for (
int i = 0; i <
leqs; i++ ) {
446 d_nnz(i) = rows_upper [ i ].giveSize() + rows_lower [ i ].giveSize();
447 total_nnz += d_nnz(i);
452 MatCreate(PETSC_COMM_SELF, &
mtrx);
454 MatSetFromOptions(
mtrx);
457 MatSetType(
mtrx, MATDENSE);
459 MatSetType(
mtrx, MATSEQAIJ);
464 MatSeqAIJSetPreallocation(
mtrx, 0, d_nnz.givePointer() );
466 MatSetOption(
mtrx, MAT_ROW_ORIENTED, PETSC_FALSE);
467 MatSetOption(
mtrx, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
501 #ifdef __PARALLEL_MODE 509 n2l->
init(eModel, di, s);
510 n2g->
init(eModel, di, s);
512 #ifdef __VERBOSE_PARALLEL 521 #ifdef __VERBOSE_PARALLEL 525 IntArray d_nnz(leqs), o_nnz(leqs), d_nnz_sym(leqs), o_nnz_sym(leqs);
531 std :: vector< IntArray > d_rows_upper(leqs), d_rows_lower(leqs);
532 std :: vector< IntArray > o_rows_upper(leqs), o_rows_lower(leqs);
541 elem->giveLocationArray(loc, s);
546 for (
int i = 1; i <= lloc.
giveSize(); i++ ) {
547 if ( ( ii = lloc.
at(i) ) ) {
548 for (
int j = 1; j <= lloc.
giveSize(); j++ ) {
549 if ( ( jj = gloc.
at(j) ) >= 0 ) {
551 if ( jj >= ( ii - 1 ) ) {
552 d_rows_upper [ ii - 1 ].insertSortedOnce(jj, loc.
giveSize() / 2);
554 d_rows_lower [ ii - 1 ].insertSortedOnce(jj, loc.
giveSize() / 2);
557 if ( jj >= ( ii - 1 ) ) {
558 o_rows_upper [ ii - 1 ].insertSortedOnce(jj, loc.
giveSize() / 2);
560 o_rows_lower [ ii - 1 ].insertSortedOnce(jj, loc.
giveSize() / 2);
573 for (
int n = 1; n <= n2lmap->
giveSize(); ++n ) {
574 if ( n2lmap->
at(n) ) {
579 for (
int i = 0; i <
leqs; i++ ) {
580 d_nnz(i) = d_rows_upper [ i ].giveSize() + d_rows_lower [ i ].giveSize();
581 o_nnz(i) = o_rows_upper [ i ].giveSize() + o_rows_lower [ i ].giveSize();
583 d_nnz_sym(i) = d_rows_upper [ i ].giveSize();
584 o_nnz_sym(i) = o_rows_upper [ i ].giveSize();
594 MatSetSizes(
mtrx, leqs, leqs, geqs, geqs);
595 MatSetType(
mtrx, MATMPIAIJ);
596 MatSetFromOptions(
mtrx);
597 MatMPIAIJSetPreallocation(
mtrx, 0, d_nnz.givePointer(), 0, o_nnz.givePointer() );
606 #ifdef __VERBOSE_PARALLEL 613 IntArray d_nnz, d_nnz_sym, indices, rowstart;
615 MatType type = MATSEQBAIJ;
617 MatCreate(PETSC_COMM_SELF, &
mtrx);
619 MatSetType(
mtrx, type);
620 MatSetFromOptions(
mtrx);
621 MatGetType(
mtrx, &type );
623 if ( strcmp(type, MATDENSE) != 0 ) {
625 std :: vector< IntArray > rows_upper(
leqs), blocks;
628 elem->giveLocationArray(loc, s);
629 for (
int ii : loc ) {
631 for (
int jj : loc ) {
633 rows_upper [ ii - 1 ].insertSortedOnce( jj - 1, loc.giveSize() );
641 std :: vector< IntArray > r_locs, c_locs;
642 for (
auto &gbc : domain->
giveBcs() ) {
646 for ( std :: size_t k = 0; k < r_locs.size(); k++ ) {
647 for (
int ii : r_locs [ k ] ) {
649 for (
int jj : c_locs [ k ] ) {
651 rows_upper [ ii - 1 ].insertSortedOnce( jj - 1, loc.
giveSize() );
662 for (
auto &row_upper : rows_upper ) {
663 nnz += row_upper.giveSize();
666 if ( strcmp(type, MATSEQAIJ) != 0 ) {
669 for (
int bsize = maxblock; bsize > 1; --bsize ) {
670 int nblocks = ceil(rows_upper.size() / (double)bsize);
672 blocks.resize(nblocks);
673 for (
int i = 0; i <
leqs; i++ ) {
674 int blockrow = i / bsize;
675 for (
int j : rows_upper [ i ] ) {
676 int blockcol = j / bsize;
677 blocks [ blockrow ].insertSortedOnce( blockcol );
682 for (
auto &block : blocks ) {
683 bnnz += block.giveSize() * bsize * bsize;
685 OOFEM_LOG_DEBUG(
"Block size %d resulted in nnz = %d -> %d using %d^2 blocks\n", bsize, nnz, bnnz, nblocks);
686 if ( bnnz <= nnz * 1.2 ) {
697 int nblocks = ceil(rows_upper.size() / (double)
blocksize);
698 d_nnz_sym.
resize(nblocks);
700 for (
int i = 0; i < nblocks; i++ ) {
701 d_nnz_sym[i] = blocks [ i ].
giveSize();
703 d_nnz[i] += d_nnz_sym[i];
704 for (
int jj : blocks [ i ] ) {
713 if ( strcmp(type, MATSEQAIJ) == 0 ) {
715 }
else if ( strcmp(type, MATSEQBAIJ) == 0 ) {
717 }
else if ( strcmp(type, MATSEQSBAIJ) == 0 ) {
722 MatSetOption(
mtrx, MAT_ROW_ORIENTED, PETSC_FALSE);
723 MatSetOption(
mtrx, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);
724 MatSetOption(
mtrx, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE);
737 #ifdef __PARALLEL_MODE 747 for (
int i = 1; i <= ndofe; i++ ) {
748 gloc.
at(i) = loc.
at(i) - 1;
759 for (
int b = 0; b < nblocke; ++b ) {
762 if ( gloc[i] != -1 && gloc[i] % this->blocksize != 0 ) {
767 for (
int k = 1; k < this->
blocksize; ++k ) {
768 bool ok_seq = gloc[i] + k == gloc[i+k];
769 bool ok_bc = gloc[i] == -1 && gloc[i+k] == -1;
770 if ( !(ok_bc || ok_seq) ) {
776 bloc[b] = gloc[i] == -1 ? -1 : gloc[i] /
blocksize;
790 #ifdef __PARALLEL_MODE 797 MatSetValues(this->
mtrx, grloc.giveSize(), grloc.givePointer(),
798 gcloc.giveSize(), gcloc.givePointer(), mat.
givePointer(), ADD_VALUES);
802 IntArray grloc(rsize), gcloc(csize);
803 for (
int i = 1; i <= rsize; i++ ) {
804 grloc.at(i) = rloc.
at(i) - 1;
807 for (
int i = 1; i <= csize; i++ ) {
808 gcloc.
at(i) = cloc.
at(i) - 1;
811 MatSetValues(this->
mtrx, rsize, grloc.givePointer(),
814 #ifdef __PARALLEL_MODE 826 return MatAssemblyBegin(this->
mtrx, MAT_FINAL_ASSEMBLY);
833 int val = MatAssemblyEnd(this->
mtrx, MAT_FINAL_ASSEMBLY);
844 MatAssembled(this->
mtrx, & assembled);
846 MatZeroEntries(this->
mtrx);
855 MatNorm(this->
mtrx, NORM_1, & norm);
881 #ifdef __PARALLEL_MODE 884 auto comm = PETSC_COMM_SELF;
886 IntArray prows = rows, pcols = cols;
895 ISCreateGeneral(comm, prows.
giveSize(), prows.
givePointer(), PETSC_USE_POINTER, & is_rows);
896 ISCreateGeneral(comm, pcols.giveSize(), pcols.givePointer(), PETSC_USE_POINTER, & is_cols);
897 MatGetSubMatrix(this->
mtrx, is_rows, is_cols, MAT_INITIAL_MATRIX, & answer->
mtrx);
898 ISDestroy(& is_rows);
899 ISDestroy(& is_cols);
921 PetscViewerSetFormat(PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_ASCII_INFO);
922 MatView(this->
mtrx, PETSC_VIEWER_STDOUT_SELF);
928 PetscViewerSetFormat(PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_ASCII_DENSE);
929 MatView(this->
mtrx, PETSC_VIEWER_STDOUT_SELF);
935 PetscViewerSetFormat(PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_ASCII_MATLAB);
936 MatView(this->
mtrx, PETSC_VIEWER_STDOUT_SELF);
943 #ifdef __PARALLEL_MODE 946 PetscViewerASCIIOpen(PETSC_COMM_SELF, fname, & viewer);
948 MatView(this->
mtrx, viewer);
949 PetscViewerDestroy(& viewer);
956 #ifdef __PARALLEL_MODE 959 VecSetSizes(* answer, this->
leqs, this->
geqs);
960 VecSetFromOptions(* answer);
964 #ifdef __PARALLEL_MODE 975 #ifdef __PARALLEL_MODE 980 VecCreateSeq(PETSC_COMM_SELF, neqs, & locVec);
982 VecScatter n2gvecscat;
984 VecScatterBegin(n2gvecscat, src, locVec, INSERT_VALUES, SCATTER_REVERSE);
985 VecScatterEnd(n2gvecscat, src, locVec, INSERT_VALUES, SCATTER_REVERSE);
986 VecScatterDestroy(& n2gvecscat);
989 VecGetArray(locVec, & ptr);
990 for (
int i = 0; i < neqs; i++ ) {
991 dest.
at(i + 1) = ptr [ i ];
994 VecRestoreArray(locVec, & ptr);
995 VecDestroy(& locVec);
1000 VecGetArray(src, & ptr);
1001 for (
int i = 0; i < neqs; i++ ) {
1002 dest.
at(i + 1) = ptr [ i ];
1004 VecRestoreArray(src, & ptr);
1005 #ifdef __PARALLEL_MODE 1017 #ifdef __PARALLEL_MODE 1024 for (
int i = 0; i < size; i++ ) {
1027 VecSetValue(dest, eqg, ptr [ i ], INSERT_VALUES);
1031 VecAssemblyBegin(dest);
1032 VecAssemblyEnd(dest);
1037 for (
int i = 0; i < size; i++ ) {
1039 VecSetValue(dest, i, ptr [ i ], INSERT_VALUES);
1042 VecAssemblyBegin(dest);
1043 VecAssemblyEnd(dest);
1044 #ifdef __PARALLEL_MODE 1063 return & this->
mtrx;
1073 return MatSetOption(this->
mtrx, op, flag);
virtual void writeToFile(const char *fname) const
Helpful for debugging, writes the matrix to given file.
int giveNumberOfColumns() const
Returns number of columns of receiver.
int nColumns
Number of columns.
void enumerate(int maxVal)
Resizes receiver and enumerates from 1 to the maximum value given.
int giveNumberOfLocalEqs()
Returns number of local eqs; ie.
int setOption(MatOption op, PetscBool flag)
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
Returns number of equations for given domain in active (current time step) time step.
int giveNumberOfDofManagers() const
Returns number of dof managers in domain.
Base class for all matrices stored in sparse format.
virtual void init(EngngModel *, int di, const UnknownNumberingScheme &n)
Initiates the receiver.
int giveDomainIndex() const
double & at(int i)
Coefficient access function.
void createVecGlobal(Vec *answer) const
Creates a global vector that fits the instance of this matrix.
virtual double computeNorm() const
Returns the norm of receiver.
virtual void zero()
Zeroes the receiver.
bool insertSortedOnce(int value, int allocChunk=0)
Inserts given value into a receiver, which is assumed to be sorted.
Ordering from oofem natural ordering (includes all local and shared eqs) to local ordering...
virtual int assembleBegin()
Starts assembling the elements.
virtual void printStatistics() const
Prints the receiver statistics (one-line) to stdout.
bool isParallel() const
Returns true if receiver in parallel mode.
virtual SparseMtrx * GiveCopy() const
Returns a newly allocated copy of receiver.
const int * givePointer() const
Breaks encapsulation.
const double * givePointer() const
Exposes the internal values of the matrix.
KSP ksp
Linear solver context.
#define OOFEM_LOG_DEBUG(...)
Ordering from oofem natural ordering (includes all local and shared eqs) to global ordering...
Class implementing an array of integers.
int & at(int i)
Coefficient access function.
virtual SparseMtrx * giveSubMatrix(const IntArray &rows, const IntArray &cols)
Natural2LocalOrdering * giveN2Lmap()
virtual int giveNewEq(int leq)
Finds the global equation from a local equation.
std::vector< std::unique_ptr< GeneralBoundaryCondition > > & giveBcs()
int giveNumberOfRows() const
Returns number of rows of receiver.
#define OOFEM_LOG_INFO(...)
int giveNumberOfDofs() const
REGISTER_SparseMtrx(CompCol, SMT_CompCol)
SparseMtrxType
Enumerative type used to identify the sparse matrix type.
virtual void init(EngngModel *, int di, const UnknownNumberingScheme &n)
Initiates the receiver.
#define VERBOSEPARALLEL_PRINT(service, str, rank)
virtual void addDiagonal(double x, FloatArray &m)
Adds x * m (treats m as a diagonal matrix, stored as an array)
virtual double & at(int i, int j)
Returns coefficient at position (i,j).
Abstract base class allowing to control the way, how equations are assigned to individual DOFs...
bool giveSymmetryFlag() const
virtual int giveNewEq(int leq)
Finds the global equation from a local equation.
int scatterL2G(const FloatArray &src, Vec dest) const
Scatters and gathers vector in local form to global (parallel) one.
virtual ParallelContext * giveParallelContext(int n)
Returns the parallel context corresponding to given domain (n) and unknown type Default implementatio...
virtual void toFloatMatrix(FloatMatrix &answer) const
Converts receiving sparse matrix to a dense float matrix.
virtual int assemble(const IntArray &loc, const FloatMatrix &mat)
Assembles sparse matrix from contribution of local elements.
Abstract base class for all active boundary conditions.
void resize(int n)
Checks size of receiver towards requested bounds.
MPI_Comm giveParallelComm()
Returns the communication object of reciever.
PETSc library mtrx representation.
SparseMtrxVersionType version
Allows to track if receiver changes.
void add(int val)
Adds given scalar to all values of receiver.
virtual void map2New(IntArray &answer, const IntArray &src, int baseOffset=0)
Class representing vector of real numbers.
Implementation of matrix containing floating point numbers.
int giveNumberOfGlobalEqs()
virtual SparseMtrxType giveType() const
Sparse matrix type identification.
double norm(const FloatArray &x)
bool kspInit
Flag if context initialized.
virtual void timesT(const FloatArray &x, FloatArray &answer) const
Evaluates .
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
IS localIS
Context or scattering/collecting parallel PETSc vectors.
This class provides an communication context for distributed memory parallelism.
virtual void add(double x, SparseMtrx &m)
Adds x * m.
void copyColumn(FloatArray &dest, int c) const
Fetches the values from the specified column.
int giveNumberOfNaturalEqs()
int giveRank() const
Returns domain rank in a group of collaborating processes (0..groupSize-1)
bool newValues
Flag if matrix has changed since last solve.
void beTranspositionOf(const FloatMatrix &src)
Assigns to the receiver the transposition of parameter.
virtual void times(const FloatArray &x, FloatArray &answer) const
Evaluates .
const double * givePointer() const
Gives the pointer to the raw data, breaking encapsulation.
Natural2GlobalOrdering * giveN2Gmap()
virtual bool isAsymmetric() const
Returns true if asymmetric.
int giveNumberOfColumns() const
Returns number of columns of receiver.
std::vector< std::unique_ptr< Element > > & giveElements()
Abstract base class representing the "problem" under consideration.
virtual int buildInternalStructure(EngngModel *eModel, int n, int m, const IntArray &I, const IntArray &J)
Builds internal structure of receiver based on I and J.
int giveSize() const
Returns the size of receiver.
the oofem namespace is to define a context or scope in which all oofem names are defined.
Domain * giveDomain(int n)
Service for accessing particular problem domain.
int scatterG2L(Vec src, FloatArray &dest) const
Scatters global vector to natural one.
DofManager * giveDofManager(int n)
Service for accessing particular domain dof manager.
virtual void printYourself() const
Prints receiver to stdout. Works only for relatively small matrices.
int giveNumberOfRows() const
Returns number of rows of receiver.
This class provides an sparse matrix interface to PETSc Matrices.
virtual int assembleEnd()
Returns when assemble is completed.
#define OOFEM_WARNING(...)
virtual ~PetscSparseMtrx()
virtual void map2New(IntArray &answer, const IntArray &src, int baseOffset=0)
virtual void giveLocationArrays(std::vector< IntArray > &rows, std::vector< IntArray > &cols, CharType type, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s)
Gives a list of location arrays that will be assembled.
void resize(int s)
Resizes receiver towards requested size.