NeoPZ
pzespmat.cpp
Go to the documentation of this file.
1 
6 #include <math.h>
7 #include "pzespmat.h"
8 #include "pzfmatrix.h"
9 
10 #include "pzerror.h"
11 #include "stdlib.h"
12 #include "pzlink.h"
13 
14 #include <sstream>
15 #include "pzlog.h"
16 #ifdef LOG4CXX
17 static LoggerPtr logger(Logger::getLogger("pz.matrix.tpzspmatrix"));
18 #endif
19 
20 using namespace std;
21 
22 /*******************/
23 /*** TPZSpMatrix ***/
24 
25 /**************************** PUBLIC ****************************/
26 /**************************/
27 /*** Construtor (int) ***/
28 template<class TVar>
29 TPZSpMatrix<TVar>::TPZSpMatrix(const int64_t rows,const int64_t cols )
30 : TPZRegisterClassId(&TPZSpMatrix::ClassId),
31 TPZMatrix<TVar>( rows, cols )
32 #ifdef WORKPOOL
33 , fWp()
34 #endif
35 {
36  fElem = new TPZLink<TPZNode>[ rows ] ;
37  if ( fElem == NULL )
38  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__, "TPZSpMatrix( dim ) <Error creating Matrix>" );
39 #ifdef WORKPOOL
40  for(int64_t i=0; i<rows; i++) fElem[i].SetWorkPool(&fWp);
41 #endif
42 }
43 
44 
45 /*****************/
46 /*** Destrutor ***/
47 
48 template<class TVar>
50 {
51  Clear();
52 }
53 
54 
55 /***********/
56 /*** Put ***/
57 
58 template<class TVar>
59 int
60 TPZSpMatrix<TVar>::Put(const int64_t row,const int64_t col,const TVar& value )
61 {
62  if ( (row >= this->Rows()) || (col >= this->Cols()) || row <0 || col<0)
63  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__, "Put <indices out of band matrix range>" );
64 
65  return( PutVal( row, col, value ) );
66 }
67 
68 
69 
70 /***********/
71 /*** Get ***/
72 
73 template<class TVar>
74 const TVar &
75 TPZSpMatrix<TVar>::Get(const int64_t row,const int64_t col ) const
76 {
77  if ( (row >= this->Rows()) || (col >= this->Cols()) || row<0 || col<0)
78  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__, "Get <indices out of band matrix range>" );
79 
80  return( GetVal( row, col ) );
81 }
82 
83 
84 
85 /**************/
86 /*** PutVal ***/
87 //
88 // Escreve um elemento na matriz, sem verificar fronteiras.
89 // O valor da linha (row) deve ser maior ou igual ao da coluna
90 // (col).
91 //
92 template<class TVar>
93 int
94 TPZSpMatrix<TVar>::PutVal(const int64_t row,const int64_t col,const TVar & value )
95 {
96  TPZLink<TPZNode> *pRow = &fElem[row];
97  TPZNode node;
98 
99  // Caso a lista esteja vazia, garante que (node.col != col).
100  node.col = -1;
101 
102  // Procura pela posicao que esta ou deveria estar o elemento.
103  pRow->Head();
104  while ( pRow->Get( &node ) && (node.col < col) )
105  pRow->Next();
106 
107  node.elem = value;
108 
109  // Se encontrou a posicao do elemento...
110  if ( node.col == col )
111  {
112  // Se o elemento e' ZERO, remove-o.
113  if ( IsZero( value ) )
114  pRow->Remove();
115 
116  // Se o elemento nao e' ZERO, muda-o.
117  else
118  pRow->Update( node );
119  }
120 
121  // Se nao encontrou a posicao do elemento e ele nao e' ZERO...
122  else if ( ! IsZero( value ) )
123  {
124  node.col = col;
125  pRow->Insert( node );
126  }
127 
128  return( 1 );
129 }
130 
131 
132 
133 /**************/
134 /*** GetVal ***/
135 //
136 // Le um elemento da matriz, sem verificar fronteiras. O valor
137 // da linha (row) deve ser maior ou igual ao da coluna (col).
138 //
139 template<class TVar>
140 const TVar &
141 TPZSpMatrix<TVar>::GetVal(const int64_t row,const int64_t col ) const
142 {
143  TPZLink<TPZNode> *pRow = &fElem[row];
144  TPZNode node;
145 
146  // Caso a lista esteja vazia, garante que (node.col != col).
147  node.col = -1;
148 
149  // Procura pela posicao que esta ou deveria estar o elemento.
150  pRow->Head();
151  while ( pRow->Get( &node ) && (node.col < col) )
152  pRow->Next();
153 
154  // Se encontrou a posicao do elemento...
155  if ( node.col == col )
156  return( pRow->GetNode()->elem );
157  else {
158  static TVar zero;
159  zero = TVar(0);
160  return( zero );
161  }
162 }
163 
164 
165 
166 
167 /******** Operacoes com matrizes FULL ********/
168 
169 /******************/
170 /*** Operator = ***/
171 
172 template<class TVar>
175 {
176  delete [] fElem;
177  fCopy( &A );
178  return( *this );
179 }
180 
181 
182 
183 /******************/
184 /*** Operator + ***/
185 
186 template<class TVar>
189 {
190  TPZSpMatrix<TVar> res( *this );
191  if ( ! res.fAdd( &A ) )
192  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__, "Operator+ (TPZSpMatrix&) <incompatible dimensions>" );
193 
194  return( res );
195 }
196 
197 
198 
199 /******************/
200 /*** Operator - ***/
201 
202 template<class TVar>
205 {
206  TPZSpMatrix<TVar> res( *this );
207  if ( ! res.fSub( &A ) )
208  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__, "Operator+( TPZSpMatrix ) <incompatible dimensions>" );
209 
210  return( res );
211 }
212 
213 /*******************/
214 /*** Operator += ***/
215 
216 template<class TVar>
219 {
220  if ( ! fAdd( &A ) )
221  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__, "operator+( TPZSpMatrix ) <incompatible dimensions>" );
222  return( *this );
223 }
224 
225 
226 
227 /*******************/
228 /*** Operator -= ***/
229 
230 template<class TVar>
233 {
234  if ( ! fSub( &A ) )
235  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__, "operator-( TPZSpMatrix ) <incompatible dimensions>" );
236  return( *this );
237 }
238 
239 /******** Operacoes com valores NUMERICOS ********/
240 /*****************************/
241 /*** Operator * ( REAL ) ***/
242 template<class TVar>
244 TPZSpMatrix<TVar>::operator*(const TVar value ) const
245 {
246  TPZSpMatrix<TVar> res( *this );
247  res.fMult( value );
248  return( res );
249 }
250 
251 /******************************/
252 /*** Operator *= ( REAL ) ***/
253 
254 template<class TVar>
257 {
258  if ( IsZero( value ) )
259  return( Reset() );
260 
261  fMult( value );
262  return( *this );
263 }
264 
265 
266 
267 /******** Metodos genericos ********/
268 
269 /*************/
270 /*** Reset ***/
271 template<class TVar>
274 {
275  for ( int64_t i = 0; i < this->Rows(); i++ )
276  fElem[i].Clear();
277  return( *this );
278 }
279 
280 /**************/
281 /*** Resize ***/
282 //
283 // Muda as dimensoes da matriz, mas matem seus valores antigos. Novas
284 // posicoes sao criadas com ZEROS.
285 //
286 template<class TVar>
287 int
288 TPZSpMatrix<TVar>::Resize(const int64_t newRows,const int64_t newCols )
289 {
290  if ( newRows == this->Rows() )
291  return( 1 );
292 
293  // Cria nova matrix.
294  TPZLink<TPZNode> *newDiag = new TPZLink<TPZNode>[ newRows ] ;
295 
296  // Copia os elementos para a nova matriz.
297  int64_t min = MIN( newRows, this->Rows() );
298  for ( int64_t i = 0; i < min; i++ )
299  newDiag[i] = fElem[i];
300 
301  // Descarta a matriz antiga e valida a nova matriz.
302  delete( fElem );
303  fElem = newDiag;
304  // fRow = fCol = newRows;
305  this->fRow=newRows;this->fCol=newCols;
306  return( 1 );
307 }
308 
309 /*************/
310 /*** Redim ***/
311 //
312 // Muda as dimensoes da matriz e ZERA seus elementos.
313 //
314 template<class TVar>
315 int
316 TPZSpMatrix<TVar>::Redim(const int64_t newRows,const int64_t newCols )
317 {
318  this->fCol = newCols;
319  delete [] fElem;
320  fElem = new TPZLink<TPZNode>[ newRows ];
321  //fRow = fCol = newRows;
322  this->fRow=newRows;
323  return( 1 );
324 }
325 
326 /*************/
327 /*** Zero ***/
328 template<class TVar>
329 int
331 {
332 
333  delete [] fElem;
334  fElem = new TPZLink<TPZNode>[ this->fRow ];
335  this->fDecomposed = 0;
336  return( 1 );
337 }
338 
339 /************************** PROTECTED **************************/
340 
341 /*****************/
342 /*** fAdd (mat) ***/
343 
344 template<class TVar>
345 int
347 {
348  if ( (this->Rows() != A->Rows()) || (this->Cols() != A->Cols()) )
349  return( 0 );
350 
351  TPZNode mNode, aNode;
352  TPZLink<TPZNode> *pm = &fElem[0];
353  TPZLink<TPZNode> *pa = &A->fElem[0];
354 
355  for ( int64_t row = 0; row < this->Rows(); row++, pm++, pa++ )
356  {
357  // Soma uma linha.
358  pm->Head();
359  pa->Head();
360  int64_t mOk = pm->Get( &mNode );
361  int64_t aOk = pa->Get( &aNode );
362 
363  // Enquanto as duas linhas tiverem elementos...
364  while ( mOk && aOk )
365  {
366  // Se as colunas forem iguais, soma os elementos.
367  if ( mNode.col == aNode.col )
368  {
369  mNode.elem += aNode.elem;
370  pm->Update( mNode );
371  pm->Next();
372  pa->Next();
373  mOk = pm->Get( &mNode );
374  aOk = pa->Get( &aNode );
375  }
376 
377  // Se a coluna desta matriz for maior, insere o elemento
378  // da matriz A.
379  else if ( mNode.col > aNode.col )
380  {
381  pm->Insert( aNode );
382 
383  // Volta a lista 'pm' `a posicao anterior.
384  pm->Next();
385 
386  pa->Next();
387  aOk = pa->Get( &aNode );
388  }
389 
390  // Se a coluna da matriz A for maior, apenas avanca a
391  // lista 'pm' para a proxima posicao.
392  else
393  {
394  pm->Next();
395  mOk = pm->Get( &mNode );
396  }
397  }
398 
399  // Se apenas a linha da matriz A tiver elementos,
400  // copia-a para esta matriz.
401  if ( aOk )
402  {
403  while ( pa->Get( &aNode ) )
404  {
405  pm->Append( aNode );
406  pa->Next();
407  }
408  }
409  }
410 
411  return( 1 );
412 }
413 
414 /*****************/
415 /*** fSub (mat) ***/
416 
417 template<class TVar>
418 int
420 {
421  if ( (this->Rows() != A->Rows()) || (this->Cols() != A->Cols()) )
422  return( 0 );
423 
424  TPZNode mNode, aNode;
425  TPZLink<TPZNode> *pm = &fElem[0];
426  TPZLink<TPZNode> *pa = &A->fElem[0];
427 
428  for ( int64_t row = 0; row < this->Rows(); row++, pm++, pa++ )
429  {
430  // Soma uma linha.
431  pm->Head();
432  pa->Head();
433  int64_t mOk = pm->Get( &mNode );
434  int64_t aOk = pa->Get( &aNode );
435 
436  // Enquanto as duas linhas tiverem elementos...
437  while ( mOk && aOk )
438  {
439  // Se as colunas forem iguais, soma os elementos.
440  if ( mNode.col == aNode.col )
441  {
442  mNode.elem -= aNode.elem;
443  pm->Update( mNode );
444  pm->Next();
445  pa->Next();
446  mOk = pm->Get( &mNode );
447  aOk = pa->Get( &aNode );
448  }
449 
450  // Se a coluna desta matriz for maior, insere o elemento
451  // com sinal trocado.
452  else if ( mNode.col > aNode.col )
453  {
454  aNode.elem = -aNode.elem;
455  pm->Insert( aNode );
456 
457  // Volta a lista 'pm' `a posicao anterior.
458  pm->Next();
459 
460  pa->Next();
461  aOk = pa->Get( &aNode );
462  }
463 
464  // Se a coluna da matriz A for maior, apenas avanca a
465  // lista 'pm' para a proxima posicao.
466  else
467  {
468  pm->Next();
469  mOk = pm->Get( &mNode );
470  }
471  }
472 
473  // Se apenas a linha da matriz A tiver elementos,
474  // copia-a para esta matriz (com sinal invertido).
475  if ( aOk )
476  {
477  while ( pa->Get( &aNode ) )
478  {
479  aNode.elem = -aNode.elem;
480  pm->Append( aNode );
481  pa->Next();
482  }
483  }
484  }
485 
486  return( 1 );
487 }
488 
489 /*******************/
490 /*** fCopy (mat) ***/
491 //
492 // Aloca as linhas e copia os elementos.
493 //
494 template<class TVar>
495 int
497 {
498  this->fCol = A->Cols();
499  this->fRow = A->Rows();
500  fElem = new TPZLink<TPZNode>[ this->fRow ];
501 
502  if ( fElem == NULL )
503  return( 0 );
504 
505  TPZLink<TPZNode> *pm = &fElem[0];
506  TPZLink<TPZNode> *pa = &A->fElem[0];
507  for ( int64_t i = 0; i < this->fRow; i++ )
508  *pm++ = *pa++;
509 
510  return( 1 );
511 }
512 
513 /*********************/
514 /*** fMult (value) ***/
515 template<class TVar>
516 int
517 TPZSpMatrix<TVar>::fMult(const TVar value )
518 {
519  TPZNode node;
520  TPZLink<TPZNode> *pm = &fElem[0];
521 
522  for ( int64_t row = 0; row < this->Rows(); row++, pm++ )
523  {
524  pm->Head();
525  while ( pm->Get( &node ) )
526  {
527  node.elem *= value;
528  pm->Update( node );
529  pm->Next();
530  }
531  }
532 
533  return( 1 );
534 }
535 
536 
537 /*************************** PRIVATE ***************************/
538 
539 
540 /****************/
541 /*** Prod Esc ***/
542 template<class TVar>
543 REAL
545  int64_t k )
546 {
547  TVar prod = 0.0;
548 
549  TPZNode node_i, node_j;
550 
551  row_i->Head();
552  row_j->Head();
553 
554  int64_t again_i = row_i->Get( &node_i ) && (node_i.col < k);
555  int64_t again_j = row_j->Get( &node_j ) && (node_j.col < k);
556 
557  while ( again_i && again_j )
558  {
559  if ( node_i.col > node_j.col )
560  {
561  row_j->Next();
562  again_j = row_j->Get( &node_j ) && (node_j.col < k);
563  }
564  else if ( node_i.col < node_j.col )
565  {
566  row_i->Next();
567  again_i = row_i->Get( &node_i ) && (node_i.col < k);
568  }
569  else if ( node_i.col == node_j.col )
570  {
571  prod += node_i.elem * node_j.elem;
572  row_i->Next();
573  row_j->Next();
574  again_i = row_i->Get( &node_i ) && (node_i.col < k);
575  again_j = row_j->Get( &node_j ) && (node_j.col < k);
576  }
577  }
578 
579  // Garante que a lista 'row_i' esteja na coluna 'k'.
580  while ( again_i )
581  {
582  row_i->Next();
583  again_i = row_i->Get( &node_i ) && (node_i.col < k);
584  }
585 
586  // Garante que a lista 'row_j' esteja na coluna 'k'.
587  while ( again_j )
588  {
589  row_j->Next();
590  again_j = row_j->Get( &node_j ) && (node_j.col < k);
591  }
592 
593  return( prod );
594 }
595 
596 template<class TVar>
598  const REAL alpha,const REAL beta,const int opt) const {
599  if ((!opt && this->Cols() != x.Rows()) || this->Rows() != x.Rows())
600  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__, "TPZSpMatrix::MultAdd <matrixs with incompatible dimensions>" );
601  if(x.Cols() != y.Cols() || x.Cols() != z.Cols() || x.Rows() != y.Rows() || x.Rows() != z.Rows()) {
602  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__,"TPZSpMatrix::MultAdd incompatible dimensions\n");
603  }
604  int64_t rows = this->Rows();
605  int64_t xcols = x.Cols();
606  int64_t ic, r;
607  this->PrepareZ(y,z,beta,opt);
608  TVar val;
609  for (ic = 0; ic < xcols; ic++) {
610  if(!opt) {
611  const REAL* firstel = &x.g(0,ic);
612  TPZNode *currentnode;
613  for ( r = 0; r < rows; r++ ) {
614  TPZLink<TPZNode> *runner = &fElem[r];
615  if(!runner->Head()) continue;
616  val = 0.;
617  currentnode = runner->GetNode();
618  while(currentnode) {
619  val += *(firstel+(currentnode->col)) * currentnode->elem;
620  runner->Next();
621  currentnode = runner->GetNode();
622  }
623  z(r,ic) += alpha*val;
624  }
625  } else {
626  TVar * firstelz = &z(0,ic);
627  TPZNode *currentnode;
628  for (r = 0; r<rows; r++) {
629  TVar elx = x.g(r,ic);
630  TPZLink<TPZNode> *runner = &fElem[r];
631  if(!runner->Head()) continue;
632  do {
633  currentnode = runner->GetNode();
634  *(firstelz+(currentnode->col)) += alpha* elx * currentnode->elem;
635  } while (runner->Next());
636  }
637  }
638  }
639 }
640 
641 
642 #ifdef OOPARLIB
643 
644 template<class TVar>
645 int TPZSpMatrix<TVar>::Unpack( TReceiveStorage *buf ){
646  TMatrix::Unpack(buf);
647  int64_t rows;
648  buf->UpkInt(&rows);
649  Redim(rows);
650  int64_t nelem;
651  int64_t col;
652  TVar val;
653  for(int64_t i=0;i<rows;i++) {
654  buf->UpkInt(&nelem);
655  buf->UpkDouble(&val);
656  buf->UpkInt(&col);
657  PutVal(i,col,val);
658  }
659  return 1;
660 }
661 
662 
663 template<class TVar>
664 TSaveable *TPZSpMatrix<TVar>::CreateInstance(TReceiveStorage *buf) {
666  m->Unpack(buf);
667  return m;
668 }
669 
670 template<class TVar>
671 int TPZSpMatrix<TVar>::Pack( TSendStorage *buf ) const {
672  TMatrix::Pack(buf);
673  TPZNode node;
674  TPZLink<TPZNode> *pm = &fElem[0];
675  int64_t rows = Rows();
676  buf->PkInt(&rows);
677  for ( int64_t row = 0; row < rows; row++, pm++ )
678  {
679  int64_t numel = 0;
680  pm->Head();
681  while ( pm->Get( &node ) )
682  {
683  numel++;
684  pm->Next();
685  }
686  buf->PkInt(&numel);
687  pm->Head();
688  while ( pm->Get( &node ) )
689  {
690  // guardar os dados de node
691  buf->PkDouble(&node.elem);
692  buf->PkInt(&node.col);
693  pm->Next();
694  }
695  }
696  return 1;
697 }
698 
699 template<class TVar>
700 int TPZSpMatrix<TVar>::DerivedFrom(const int64_t Classid) const {
701  return TMatrix::DerivedFrom(Classid);
702 }
703 
704 template<class TVar>
705 int TPZSpMatrix<TVar>::DerivedFrom(const char *classname) const {
706 
707  if(!strcmp(ClassName(),classname)) return 1;
708  return TMatrix::DerivedFrom(classname);
709 }
710 
711 #endif
712 
713 template<class TVar>
715  return Hash("TPZSpMatrix") ^ TPZMatrix<TVar>::ClassId() << 1;
716 }
717 template class TPZSpMatrix<REAL>;
TPZSpMatrix & operator-=(const TPZSpMatrix< TVar > &A)
Definition: pzespmat.cpp:232
TPZLink< TPZNode > * fElem
Definition: pzespmat.h:192
int fSub(const TPZSpMatrix< TVar > *const A)
Definition: pzespmat.cpp:419
const int64_t numel
Number of elements to test.
Definition: pzsubcmesh.cpp:47
bool IsZero(long double a)
Returns if the value a is close Zero as the allowable tolerance.
Definition: pzreal.h:668
Contains definitions to LOGPZ_DEBUG, LOGPZ_INFO, LOGPZ_WARN, LOGPZ_ERROR and LOGPZ_FATAL, and the implementation of the inline InitializePZLOG(string) function using log4cxx library or not. It must to be called out of "#ifdef LOG4CXX" scope.
static int Error(const char *msg, const char *msg2=0)
Returns error messages.
Definition: pzmatrix.cpp:1402
int fCopy(const TPZSpMatrix< TVar > *const A)
Definition: pzespmat.cpp:496
const TVar & GetVal(const int64_t row, const int64_t col) const override
Get values without bounds checking This method is faster than "Get" if DEBUG is defined.
Definition: pzespmat.cpp:141
int Clear() override
It clears data structure.
Definition: pzespmat.h:179
int fMult(TVar val)
Definition: pzespmat.cpp:517
Defines PZError.
int64_t fRow
Number of rows in matrix.
Definition: pzmatrix.h:779
void PrepareZ(const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar beta, const int opt) const
Is an auxiliar method used by MultiplyAdd.
Definition: pzmatrix.cpp:135
TPZSpMatrix()
Simple constructor.
Definition: pzespmat.h:57
char fDecomposed
Decomposition type used to decompose the current matrix.
Definition: pzmatrix.h:783
REAL val(STATE &number)
Returns value of the variable.
Definition: pzartdiff.h:23
Defines sparce matrix class. Matrix.
Definition: pzespmat.h:37
virtual void MultAdd(const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const REAL alpha, const REAL beta=0., const int opt=0) const override
Definition: pzespmat.cpp:597
Defines a node.
Definition: pzespmat.h:49
Contains TPZMatrixclass which implements full matrix (using column major representation).
~TPZSpMatrix()
Simple destructor.
Definition: pzespmat.cpp:49
TPZSpMatrix & operator=(const TPZSpMatrix< TVar > &A)
Generic overloaded operator.
Definition: pzespmat.cpp:174
int Zero() override
Zeroes all elements of the matrix.
Definition: pzespmat.cpp:330
Contains TPZSpMatrix class which defines sparse matrix class.
int64_t Rows() const
Returns number of rows.
Definition: pzmatrix.h:803
static TPZSavable * CreateInstance(const int &classId)
Definition: TPZSavable.cpp:120
TPZSpMatrix & Reset()
Desallocate all the elements of the matrix.
Definition: pzespmat.cpp:273
string res
Definition: test.py:151
#define MIN(a, b)
Gets minime value between a and b.
Definition: pzreal.h:85
Full matrix class. Matrix.
Definition: pzfmatrix.h:32
const TVar & Get(const int64_t row, const int64_t col) const override
Get value with bound checking.
Definition: pzespmat.cpp:75
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
TPZSpMatrix & operator+=(const TPZSpMatrix< TVar > &A)
Definition: pzespmat.cpp:218
int64_t fCol
Number of cols in matrix.
Definition: pzmatrix.h:781
int fAdd(const TPZSpMatrix< TVar > *const A)
Definition: pzespmat.cpp:346
TPZSpMatrix operator+(const TPZSpMatrix< TVar > &A) const
Definition: pzespmat.cpp:188
REAL ProdEsc(TPZLink< TPZNode > *row_i, TPZLink< TPZNode > *row_j, int64_t k)
Computes dot product with respect to lines row_i e row_j using only elements that belongs to columns ...
Definition: pzespmat.cpp:544
int Put(const int64_t row, const int64_t col, const TVar &value) override
Put values with bounds checking if DEBUG variable is defined.
Definition: pzespmat.cpp:60
TPZSpMatrix operator-() const
Definition: pzespmat.h:121
int64_t Cols() const
Returns number of cols.
Definition: pzmatrix.h:809
TPZSpMatrix & operator*=(const TVar v)
Definition: pzespmat.cpp:256
int PutVal(const int64_t row, const int64_t col, const TVar &element) override
Put values without bounds checking This method is faster than "Put" if DEBUG is defined.
Definition: pzespmat.cpp:94
int ClassId() const override
Define the class id associated with the class.
Definition: pzespmat.cpp:714
clarg::argString m("-m", "input matrix file name (text format)", "matrix.txt")
TVar & g(const int64_t row, const int64_t col) const
Definition: pzfmatrix.h:594
TPZSpMatrix operator*(const TVar v) const
Numerical values operator.
Definition: pzespmat.cpp:244
int ClassId() const override
Define the class id associated with the class.
Definition: pzmatrix.h:957
int Redim(const int64_t newRows, const int64_t newCols) override
Redimension the matrix to newRows x newCols arrange and zeroes its elements.
Definition: pzespmat.cpp:316
int Resize(const int64_t newRows, const int64_t newCols) override
Redimension the matrix to newRows x newCols arrange.
Definition: pzespmat.cpp:288
Root matrix class (abstract). Matrix.
Definition: pzmatrix.h:60