OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
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 - 2013 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 
43 namespace oofem {
47 class OOFEM_EXPORT PetscSparseMtrx : public SparseMtrx
48 {
49 protected:
50  Mat mtrx;
51  bool symmFlag;
52  MatType mType;
53  int leqs;
54  int geqs;
55  int blocksize;
56  int di;
58 
60  KSP ksp;
62  bool kspInit;
64  bool newValues;
65 
67  IS localIS, globalIS;
68 
69 public:
70  PetscSparseMtrx(int n, int m);
72 
73  virtual ~PetscSparseMtrx();
74 
75  // Overloaded methods:
76  virtual SparseMtrx *GiveCopy() const;
77  virtual void times(const FloatArray &x, FloatArray &answer) const;
78  virtual void timesT(const FloatArray &x, FloatArray &answer) const;
79  virtual void times(const FloatMatrix &B, FloatMatrix &answer) const;
80  virtual void timesT(const FloatMatrix &B, FloatMatrix &answer) const;
81  virtual void times(double x);
82  virtual void add(double x, SparseMtrx &m);
83  virtual void addDiagonal(double x, FloatArray &m);
84  virtual int buildInternalStructure(EngngModel *eModel, int n, int m, const IntArray &I, const IntArray &J);
85  virtual int buildInternalStructure(EngngModel *eModel, int di, const UnknownNumberingScheme &s);
86  virtual int buildInternalStructure(EngngModel *eModel, int di, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s);
87  virtual int assemble(const IntArray &loc, const FloatMatrix &mat);
88  virtual int assemble(const IntArray &rloc, const IntArray &cloc, const FloatMatrix &mat);
89  virtual int assembleBegin();
90  virtual int assembleEnd();
91  virtual SparseMtrx *giveSubMatrix(const IntArray &rows, const IntArray &cols);
92  virtual bool canBeFactorized() const { return false; }
93  virtual SparseMtrx *factorized() { return NULL; }
94  virtual FloatArray *backSubstitutionWith(FloatArray &y) const { return NULL; }
95  virtual void zero();
96  virtual double computeNorm() const;
97  virtual double &at(int i, int j);
98  virtual double at(int i, int j) const;
99  virtual void toFloatMatrix(FloatMatrix &answer) const;
100  virtual void printStatistics() const;
101  virtual void printYourself() const;
102  void printMatlab() const;
103  virtual SparseMtrxType giveType() const;
104  virtual bool isAsymmetric() const;
105  virtual void writeToFile(const char *fname) const;
106  virtual const char *giveClassName() const { return "PetscSparseMtrx"; }
107 
109  void createVecGlobal(Vec *answer) const;
111  int scatterG2L(Vec src, FloatArray &dest) const;
116  int scatterL2G(const FloatArray &src, Vec dest) const;
117 
118  // Internals (should be documented)
119  Mat *giveMtrx();
120  bool giveSymmetryFlag() const;
121  int setOption(MatOption op, PetscBool flag);
122  int giveLeqs();
123  int giveDomainIndex() const;
124 
125  friend class PetscSolver;
126 };
127 } // end namespace oofem
128 #endif
Base class for all matrices stored in sparse format.
Definition: sparsemtrx.h:60
virtual bool canBeFactorized() const
Determines, whether receiver can be factorized.
KSP ksp
Linear solver context.
Class implementing an array of integers.
Definition: intarray.h:61
SparseMtrxType
Enumerative type used to identify the sparse matrix type.
Abstract base class allowing to control the way, how equations are assigned to individual DOFs...
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
Implements the solution of linear system of equation in the form using solvers from PETSc library...
Definition: petscsolver.h:51
bool kspInit
Flag if context initialized.
virtual SparseMtrx * factorized()
Returns the receiver factorized.
virtual FloatArray * backSubstitutionWith(FloatArray &y) const
Computes the solution of linear system where A is receiver.
IS localIS
Context or scattering/collecting parallel PETSc vectors.
bool newValues
Flag if matrix has changed since last solve.
virtual const char * giveClassName() const
Abstract base class representing the "problem" under consideration.
Definition: engngm.h:181
the oofem namespace is to define a context or scope in which all oofem names are defined.
This class provides an sparse matrix interface to PETSc Matrices.

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:30 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011