OOFEM 3.0
Loading...
Searching...
No Matches
petscsparsemtrx.h
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34#ifndef petscsparsemtrx_h
35#define petscsparsemtrx_h
36
37#include "sparsemtrx.h"
38
39#include <petscksp.h>
40
41#define _IFT_PetscSparseMtrx_Name "petsc"
42
43namespace oofem {
48{
49protected:
50 Mat mtrx;
52 MatType mType;
53 int leqs;
54 int geqs;
56 int di;
58
60 KSP ksp;
62 bool kspInit;
65
68
69public:
70 PetscSparseMtrx(int n=0, int m=0);
71
72 virtual ~PetscSparseMtrx();
73
74 std::unique_ptr<SparseMtrx> clone() const override;
75 void times(const FloatArray &x, FloatArray &answer) const override;
76 void timesT(const FloatArray &x, FloatArray &answer) const override;
77 void times(const FloatMatrix &B, FloatMatrix &answer) const override;
78 void timesT(const FloatMatrix &B, FloatMatrix &answer) const override;
79 void times(double x) override;
80 void add(double x, SparseMtrx &m) override;
81 void addDiagonal(double x, FloatArray &m) override;
82 int buildInternalStructure(EngngModel *eModel, int n, int m, const IntArray &I, const IntArray &J) override;
83 int buildInternalStructure(EngngModel *eModel, int di, const UnknownNumberingScheme &s) override;
84 int buildInternalStructure(EngngModel *eModel, int di, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s) override;
85 int assemble(const IntArray &loc, const FloatMatrix &mat) override;
86 int assemble(const IntArray &rloc, const IntArray &cloc, const FloatMatrix &mat) override;
87 int assembleBegin() override;
88 int assembleEnd() override;
89 std::unique_ptr<SparseMtrx> giveSubMatrix(const IntArray &rows, const IntArray &cols) override;
90 bool canBeFactorized() const override { return false; }
91 SparseMtrx *factorized() override { return NULL; }
92 FloatArray *backSubstitutionWith(FloatArray &y) const override { return NULL; }
93 void zero() override;
94 double computeNorm() const override;
95 double &at(int i, int j) override;
96 double at(int i, int j) const override;
97 void toFloatMatrix(FloatMatrix &answer) const override;
98 void printStatistics() const override;
99 void printYourself() const override;
100 void printMatlab() const;
101 SparseMtrxType giveType() const override;
102 bool isAsymmetric() const override;
103 void writeToFile(const char *fname) const override;
104 const char *giveClassName() const override { return "PetscSparseMtrx"; }
105
107 void createVecGlobal(Vec *answer) const;
109 int scatterG2L(Vec src, FloatArray &dest) const;
114 int scatterL2G(const FloatArray &src, Vec dest) const;
115
116 // Internals (should be documented)
117 Mat *giveMtrx();
118 bool giveSymmetryFlag() const;
119 int setOption(MatOption op, PetscBool flag);
120 int giveLeqs();
121 int giveDomainIndex() const;
122
123 friend class PetscSolver;
124};
125} // end namespace oofem
126#endif
std::unique_ptr< SparseMtrx > giveSubMatrix(const IntArray &rows, const IntArray &cols) override
IS localIS
Context or scattering/collecting parallel PETSc vectors.
PetscSparseMtrx(int n=0, int m=0)
void times(const FloatArray &x, FloatArray &answer) const override
bool kspInit
Flag if context initialized.
FloatArray * backSubstitutionWith(FloatArray &y) const override
void add(double x, SparseMtrx &m) override
void addDiagonal(double x, FloatArray &m) override
KSP ksp
Linear solver context.
bool canBeFactorized() const override
Determines, whether receiver can be factorized.
bool newValues
Flag if matrix has changed since last solve.
int assembleBegin() override
Starts assembling the elements.
const char * giveClassName() const override
std::unique_ptr< SparseMtrx > clone() const override
SparseMtrx * factorized() override
int assembleEnd() override
Returns when assemble is completed.
void timesT(const FloatArray &x, FloatArray &answer) const override
int assemble(const IntArray &loc, const FloatMatrix &mat) override
int buildInternalStructure(EngngModel *eModel, int n, int m, const IntArray &I, const IntArray &J) override
SparseMtrx(int n=0, int m=0)
Definition sparsemtrx.h:89
FloatMatrixF< N, M > zero()
Constructs a zero matrix (this is the default behavior when constructing a matrix,...
#define OOFEM_EXPORT
Definition oofemcfg.h:7

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011