OOFEM 3.0
Loading...
Searching...
No Matches
sparsemtrx.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
35#ifndef sparsemtrx_h
36#define sparsemtrx_h
37
38#include "oofemenv.h"
39#include "floatarray.h"
40#include "floatmatrix.h"
41#include "intarray.h"
42#include "error.h"
43#include "sparsemtrxtype.h"
44
45#include <memory>
46
47namespace oofem {
48class EngngModel;
49class TimeStep;
51
63{
64public:
66
67protected:
69 int nRows;
72
83
84public:
89 SparseMtrx(int n=0, int m=0) : nRows(n), nColumns(m), version(0) { }
91 virtual ~SparseMtrx() { }
92
94 SparseMtrxVersionType giveVersion() { return this->version; }
95
102 void checkBounds(int i, int j) const {
103 if ( i <= 0 ) {
104 OOFEM_ERROR("matrix error on rows : %d <= 0", i);
105 } else if ( j <= 0 ) {
106 OOFEM_ERROR("matrix error on columns : %d <= 0", j);
107 } else if ( i > nRows ) {
108 OOFEM_ERROR("matrix error on rows : %d => %d", i, nRows);
109 } else if ( j > nColumns ) {
110 OOFEM_ERROR("matrix error on columns : %d => %d", j, nColumns);
111 }
112 }
113
114 int giveNumberOfRows() const { return nRows; }
116 int giveNumberOfColumns() const { return nColumns; }
118 bool isSquare() const { return nRows == nColumns; }
120 bool isNotEmpty() const { return nRows > 0 && nColumns > 0; }
121
127 virtual std::unique_ptr<SparseMtrx> clone() const
128 {
129 OOFEM_ERROR("Not implemented");
130 }
131
137 virtual void times(const FloatArray &x, FloatArray &answer) const { OOFEM_ERROR("Not implemented"); }
143 virtual void timesT(const FloatArray &x, FloatArray &answer) const { OOFEM_ERROR("Not implemented"); }
149 virtual void times(const FloatMatrix &B, FloatMatrix &answer) const { OOFEM_ERROR("Not implemented"); }
155 virtual void timesT(const FloatMatrix &B, FloatMatrix &answer) const { OOFEM_ERROR("Not implemented"); }
160 virtual void times(double x) { OOFEM_ERROR("Not implemented"); }
166 virtual void add(double x, SparseMtrx &m) { OOFEM_ERROR("Not implemented"); }
172 virtual void addDiagonal(double x, FloatArray &m) {
173 for ( int i = 1; i <= m.giveSize(); ++i ) {
174 this->at(i, i) += x * m.at(i);
175 }
176 }
177
185 virtual int buildInternalStructure(EngngModel *eModel, int n, int m, const IntArray &I, const IntArray &J) { OOFEM_ERROR("Not implemented"); }
198 virtual int buildInternalStructure(EngngModel *eModel, int di, const UnknownNumberingScheme &s) = 0;
208 virtual int buildInternalStructure(EngngModel *eModel, int di, const UnknownNumberingScheme &r_s,
209 const UnknownNumberingScheme &c_s) {
210 OOFEM_ERROR("Not implemented");
211 }
212
220 virtual int assemble(const IntArray &loc, const FloatMatrix &mat) = 0;
231 virtual int assemble(const IntArray &rloc, const IntArray &cloc, const FloatMatrix &mat) = 0;
232
234 virtual int assembleBegin() { return 1; }
236 virtual int assembleEnd() { return 1; }
237
239 virtual bool canBeFactorized() const = 0;
244 virtual SparseMtrx *factorized() { return NULL; }
252 virtual FloatArray *backSubstitutionWith(FloatArray &y) const { return NULL; }
254 virtual void zero() = 0;
255
257 virtual double computeNorm() const {
258 OOFEM_ERROR("Not implemented");
259 }
260
261 virtual std::unique_ptr<SparseMtrx> giveSubMatrix(const IntArray &rows, const IntArray &cols)
262 { OOFEM_ERROR("Not implemented"); }
263
265 virtual double &at(int i, int j) = 0;
267 virtual double at(int i, int j) const = 0;
269 virtual bool isAllocatedAt(int i, int j) const { return false; }
271 virtual void toFloatMatrix(FloatMatrix &answer) const { OOFEM_ERROR("Not implemented"); }
273 virtual void printStatistics() const { OOFEM_LOG_INFO("Not implemented file %s, line %d\n", __FILE__, __LINE__); }
275 virtual void printYourself() const { OOFEM_LOG_INFO("Not implemented file %s, line %d\n", __FILE__, __LINE__); }
277 virtual void writeToFile(const char *fname) const { OOFEM_LOG_INFO("Not implemented file %s, line %d\n", __FILE__, __LINE__); }
279 virtual SparseMtrxType giveType() const = 0;
281 virtual bool isAsymmetric() const = 0;
282
283 virtual const char *giveClassName() const = 0;
285 std :: string errorInfo(const char *func) const { return std :: string(giveClassName()) + func; }
287
288
290 {
291 FloatArray answer;
292 this->times(x, answer);
293 return answer;
294 }
295
297 {
298 FloatArray answer;
299 this->timesT(x, answer);
300 return answer;
301 }
302
303};
304} // end namespace oofem
305#endif // sparsemtrx_h
Vector operator*(const Vector &a, double b)
Definition CSG.h:68
double & at(Index i)
Definition floatarray.h:202
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
virtual int assemble(const IntArray &rloc, const IntArray &cloc, const FloatMatrix &mat)=0
long SparseMtrxVersionType
Definition sparsemtrx.h:65
SparseMtrxVersionType giveVersion()
Return receiver version.
Definition sparsemtrx.h:94
virtual double & at(int i, int j)=0
Returns coefficient at position (i,j).
int nColumns
Number of columns.
Definition sparsemtrx.h:71
virtual SparseMtrxType giveType() const =0
Sparse matrix type identification.
virtual void writeToFile(const char *fname) const
Helpful for debugging, writes the matrix to given file.
Definition sparsemtrx.h:277
SparseMtrxVersionType version
Definition sparsemtrx.h:82
virtual bool canBeFactorized() const =0
Determines, whether receiver can be factorized.
void checkBounds(int i, int j) const
Definition sparsemtrx.h:102
virtual SparseMtrx * factorized()
Definition sparsemtrx.h:244
virtual void add(double x, SparseMtrx &m)
Definition sparsemtrx.h:166
virtual void toFloatMatrix(FloatMatrix &answer) const
Converts receiving sparse matrix to a dense float matrix.
Definition sparsemtrx.h:271
int giveNumberOfRows() const
Returns number of rows of receiver.
Definition sparsemtrx.h:114
virtual FloatArray * backSubstitutionWith(FloatArray &y) const
Definition sparsemtrx.h:252
int giveNumberOfColumns() const
Returns number of columns of receiver.
Definition sparsemtrx.h:116
virtual void addDiagonal(double x, FloatArray &m)
Definition sparsemtrx.h:172
virtual const char * giveClassName() const =0
virtual std::unique_ptr< SparseMtrx > giveSubMatrix(const IntArray &rows, const IntArray &cols)
Definition sparsemtrx.h:261
virtual void timesT(const FloatArray &x, FloatArray &answer) const
Definition sparsemtrx.h:143
virtual bool isAsymmetric() const =0
Returns true if asymmetric.
virtual std::unique_ptr< SparseMtrx > clone() const
Definition sparsemtrx.h:127
virtual int buildInternalStructure(EngngModel *eModel, int n, int m, const IntArray &I, const IntArray &J)
Definition sparsemtrx.h:185
virtual int buildInternalStructure(EngngModel *eModel, int di, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s)
Definition sparsemtrx.h:208
std::string errorInfo(const char *func) const
Error printing helper.
Definition sparsemtrx.h:285
SparseMtrx(int n=0, int m=0)
Definition sparsemtrx.h:89
bool isSquare() const
Returns nonzero if receiver is square matrix.
Definition sparsemtrx.h:118
int nRows
Number of rows.
Definition sparsemtrx.h:69
virtual void timesT(const FloatMatrix &B, FloatMatrix &answer) const
Definition sparsemtrx.h:155
virtual void times(double x)
Definition sparsemtrx.h:160
virtual double at(int i, int j) const =0
Returns coefficient at position (i,j).
virtual int assembleBegin()
Starts assembling the elements.
Definition sparsemtrx.h:234
virtual ~SparseMtrx()
Destructor.
Definition sparsemtrx.h:91
FloatArray trans_mult(const FloatArray &x) const
IML compatibility, .
Definition sparsemtrx.h:296
virtual void printYourself() const
Prints receiver to stdout. Works only for relatively small matrices.
Definition sparsemtrx.h:275
virtual void times(const FloatMatrix &B, FloatMatrix &answer) const
Definition sparsemtrx.h:149
virtual int assemble(const IntArray &loc, const FloatMatrix &mat)=0
bool isNotEmpty() const
Tests for empty matrix.
Definition sparsemtrx.h:120
virtual int assembleEnd()
Returns when assemble is completed.
Definition sparsemtrx.h:236
virtual bool isAllocatedAt(int i, int j) const
Checks whether memory is allocated at position (i,j).
Definition sparsemtrx.h:269
virtual int buildInternalStructure(EngngModel *eModel, int di, const UnknownNumberingScheme &s)=0
virtual double computeNorm() const
Returns the norm of receiver.
Definition sparsemtrx.h:257
virtual void zero()=0
Zeroes the receiver.
virtual void times(const FloatArray &x, FloatArray &answer) const
Definition sparsemtrx.h:137
virtual void printStatistics() const
Prints the receiver statistics (one-line) to stdout.
Definition sparsemtrx.h:273
#define OOFEM_ERROR(...)
Definition error.h:79
#define OOFEM_LOG_INFO(...)
Definition logger.h:143
#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