NeoPZ
pzblockdiag.cpp
Go to the documentation of this file.
1 
6 #include "pzfmatrix.h"
7 #include "pzblockdiag.h"
8 
9 
10 #include <math.h>
11 #include <stdlib.h>
12 #include <stdio.h>
13 
14 #include "pzlog.h"
15 #include <sstream>
16 #include "pzstack.h"
17 
18 #ifdef LOG4CXX
19 static LoggerPtr logger(Logger::getLogger("pz.StrMatrix"));
20 #endif
21 
22 using namespace std;
23 
24 template<class TVar>
26 
27  int64_t firstpos = fBlockPos[i];
28  int64_t bsize = fBlockSize[i];
29 
30  int64_t r,c;
31  for(r=0; r<bsize; r++) {
32  for(c=0; c<bsize; c++) {
33  fStorage[firstpos+r+bsize*c] += block(r,c);
34  }
35  }
36 }
37 
38 template<class TVar>
40  int64_t firstpos = fBlockPos[i];
41  int64_t bsize = fBlockSize[i];
42 
43  int64_t r,c;
44  for(r=0; r<bsize; r++) {
45  for(c=0; c<bsize; c++) {
46  fStorage[firstpos+r+bsize*c] = block(r,c);
47  }
48  }
49 }
50 
51 template<class TVar>
53  int64_t firstpos = fBlockPos[i];
54  int64_t bsize = fBlockSize[i];
55  block.Redim(bsize,bsize);
56  int64_t r,c;
57  for(r=0; r<bsize; r++) {
58  for(c=0; c<bsize; c++) {
59  block(r,c) = fStorage[firstpos+r+bsize*c];
60  }
61  }
62 }
63 template<class TVar>
65  int64_t nblock = blocksize.NElements();
66 #ifdef LOG4CXX
67  if (logger->isDebugEnabled())
68  {
69  std::stringstream sout;
70  sout << "Number of blocks \t" << nblock;
71  LOGPZ_DEBUG(logger,sout.str());
72  }
73 #endif
74  fBlockSize = blocksize;
75  fBlockPos.Resize(nblock+1,0);
76  int64_t b;
77  int64_t ndata = 0;
78  int64_t neq = 0;
79  int bsize;
80  for(b=0; b<nblock; b++) {
81  bsize = blocksize[b];
82  fBlockPos[b+1] = fBlockPos[b]+bsize*bsize;
83  ndata += bsize*bsize;
84  neq += bsize;
85  }
86 #ifdef LOG4CXX
87  if(ndata > 10000000)
88  {
89  std::stringstream sout;
90  sout << "Calling fStorage.Resize(ndata,0.) with ndata = " << ndata;
91  LOGPZ_DEBUG(logger,sout.str());
92  }
93 #endif
94 
95  fStorage.Fill(0.,0);
96  fStorage.Resize(ndata,0.);
97  this->fDecomposed = 0;
98  this->fRow = neq;
99  this->fCol = neq;
100 }
101 
102 template<class TVar>
104  if(mat.Rows() != this->Rows()) {
105  cout << "TPZBlockDiagonal::BuildFromMatrix wrong data structure\n";
106  return;
107  }
108  int64_t nblock = fBlockSize.NElements();
109  int64_t b,eq=0;
110  for(b=0; b<nblock; b++) {
111  int r,c,bsize = fBlockSize[b];
112  int64_t pos = fBlockPos[b];
113  for(r=0; r<bsize; r++){
114  for(c=0; c<bsize; c++) {
115  fStorage[pos+r+c*bsize] = mat.GetVal(eq+r,eq+c);
116  }
117  }
118  eq += bsize;
119  }
120 }
121 
122 /*******************/
123 /*** Constructor ***/
124 
125 template<class TVar>
128 TPZMatrix<TVar>(), fStorage(), fBlockPos(1,0),fBlockSize()
129 {
130 
131 }
132 
133 template<class TVar>
136 TPZMatrix<TVar>(), fStorage(), fBlockPos(1,0),fBlockSize()
137 {
138  Initialize(blocksize);
139 }
140 
141 /********************/
142 /*** Constructors ***/
143 
144 template<class TVar>
147 TPZMatrix<TVar>(), fBlockSize(blocksizes)
148 {
149  int64_t nblock = blocksizes.NElements();
150  fBlockPos.Resize(nblock+1,0);
151  int64_t b;
152  int64_t ndata = 0;
153  int64_t neq = 0;
154  int bsize;
155  for(b=0; b<nblock; b++) {
156  bsize = blocksizes[b];
157  fBlockPos[b+1] = fBlockPos[b]+bsize*bsize;
158  ndata += bsize*bsize;
159  neq += bsize;
160  }
161  fStorage.Resize(ndata,0.);
162  this->fRow = neq;
163  this->fCol = neq;
164  int64_t pos;
165  int64_t eq = 0, r, c;
166  for(b=0; b<nblock; b++) {
167  bsize = fBlockSize[b];
168  pos = fBlockPos[b];
169  for(r=0; r<bsize; r++) {
170  for(c=0; c<bsize; c++) {
171  fStorage[pos+r+bsize*c]= glob.GetVal(eq+r,eq+c);
172  }
173  }
174  eq += bsize;
175  }
176 
177 }
178 
179 /*********************************/
180 /*** Constructor( TPZBlockDiagonal& ) ***/
181 template<class TVar>
184 TPZMatrix<TVar>( A.Dim(), A.Dim() ), fStorage(A.fStorage),
186 {
187 }
188 
189 /******************/
190 /*** Destructor ***/
191 template<class TVar>
193 {
194 }
195 
196 /***********/
197 /*** Put ***/
198 template<class TVar>
199 int TPZBlockDiagonal<TVar>::Put(const int64_t row,const int64_t col,const TVar& value )
200 {
201  // cout << "TPZBlockDiagonal.Put should not be called\n";
202  if ( (row >= Dim()) || (col >= Dim()) || row<0 || col<0)
203  {
204  cout << "TPZBlockDiagonal::Put: " << row << "," << col << "," << Dim();
205  cout << "\n";
206  return( 0 );
207  }
208 
209  return( PutVal( row, col, value ) );
210 }
211 
212 /***********/
213 /*** PutVal ***/
214 template<class TVar>
215 int TPZBlockDiagonal<TVar>::PutVal(const int64_t row,const int64_t col,const TVar& value )
216 {
217  int64_t b = 0;
218  int64_t nb = fBlockSize.NElements();
219  if(nb==0) {
220  cout << "TPZBlockDiagonal::PutVal called with parameters out of range\n";
221  return -1;
222  }
223  int64_t eq=0;
224  int64_t bsize = fBlockSize[b];
225  while(eq+bsize <= row && b < nb) {
226  eq+=bsize;
227  b++;
228  bsize = fBlockSize[b];
229  }
230  if(b==nb) {
231  cout << "TPZBlockDiagonal::PutVal wrong data structure\n";
232  return -1;
233  }
234  if(col < eq || col >= eq+bsize) {
235  if(value != TVar(0.)) {
236  cout << "TPZBlockDiagonal::PutVal, indices row col out of range\n";
237  return -1;
238  } else {
239  return 0;
240  }
241  }
242  fStorage[fBlockPos[b]+row-eq+bsize*(col-eq)] = value;
243  return 0;
244 }
245 
246 
247 
248 /***********/
249 /*** Get ***/
250 
251 template<class TVar>
252 const TVar&
253 TPZBlockDiagonal<TVar>::Get(const int64_t row,const int64_t col ) const
254 {
255  if ( (row >= Dim()) || (col >= Dim()) )
256  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__, "TPZBlockDiagonal::Get <indices out of band matrix range>" );
257 
258  return( GetVal( row, col ) );
259 }
260 
261 template<class TVar>
262 TVar &
263 TPZBlockDiagonal<TVar>::operator()(const int64_t row, const int64_t col) {
264 
265  int64_t b = 0;
266  int64_t nb = fBlockSize.NElements();
267  if(nb==0) {
268  cout << "TPZBlockDiagonal::operator() called with parameters out of range\n";
269  static TVar zero = 0.;
270  return zero;
271  }
272  int64_t eq=0;
273  int64_t bsize = fBlockSize[b];
274  while(eq+bsize <= row && b < nb) {
275  eq+=bsize;
276  b++;
277  bsize = fBlockSize[b];
278  }
279  if(b==nb) {
280  cout << "TPZBlockDiagonal::operator() wrong data structure\n";
281  static TVar zero = 0.;
282  return zero;
283  }
284  if(col < eq || col >= eq+bsize) {
285  cout << "TPZBlockDiagonal::operator(), indices row col out of range\n";
286  static TVar zero = 0.;
287  return zero;
288  }
289  return fStorage[fBlockPos[b]+row-eq+bsize*(col-eq)];
290 }
291 
292 /***********/
293 /*** GetVal ***/
294 template<class TVar>
295 const TVar &TPZBlockDiagonal<TVar>::GetVal(const int64_t row,const int64_t col ) const
296 {
297  int64_t b = 0;
298  int64_t nb = fBlockSize.NElements();
299  if(nb==0) {
300  cout << "TPZBlockDiagonal::GetVal called with parameters out of range\n";
301  }
302  int64_t eq=0;
303  int64_t bsize = fBlockSize[b];
304  while(eq+bsize <= row && b < nb) {
305  eq+=bsize;
306  b++;
307  bsize = fBlockSize[b];
308  }
309  if(b==nb) {
310  cout << "TPZBlockDiagonal::GetVal wrong data structure\n";
311  }
312  if(col < eq || col >= eq+bsize) {
313  //cout << "TPZBlockDiagonal::GetVal, indices row col out of range\n";
314  static TVar zero = 0.;
315  return zero;
316  }
317  return fStorage[fBlockPos[b]+row-eq+bsize*(col-eq)];
318 }
319 
320 
321 /******** Operacoes com MATRIZES GENERICAS ********/
322 
323 /*******************/
324 /*** MultiplyAdd ***/
325 //
326 // perform a multiply add operation to be used by iterative solvers
327 //
328 
329 template<class TVar>
331  const TVar alpha,const TVar beta ,const int opt) const {
332  // Computes z = beta * y + alpha * opt(this)*x
333  // z and x cannot overlap in memory
334 
335  if ((!opt && this->Cols() != x.Rows()) || this->Rows() != x.Rows())
336  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__, "TPZBlockDiagonal::MultAdd <matrixs with incompatible dimensions>" );
337  if(x.Cols() != y.Cols() || x.Cols() != z.Cols() || x.Rows() != y.Rows() || x.Rows() != z.Rows()) {
338  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__,"TPZBlockDiagonal::MultAdd incompatible dimensions\n");
339  }
340 
341  this->PrepareZ(y,z,beta,opt);
342  int64_t xcols = x.Cols();
343  int64_t nb= fBlockSize.NElements();
344  int64_t ic, b, eq=0;
345  int bsize, r, c;
346  if(opt == 0) {
347  for (ic = 0; ic < xcols; ic++) {
348  eq=0;
349  for(b=0; b<nb; b++) {
350  bsize = fBlockSize[b];
351  int64_t pos = fBlockPos[b];
352  for(r=0; r<bsize; r++) {
353  for(c=0; c<bsize; c++) {
354  z(eq+r,ic) += alpha*fStorage[pos+r+bsize*c]*x.GetVal((eq+c),ic);
355  }
356  }
357  eq += bsize;
358  }
359  }
360  } else {
361  cout << "xcols \t" << xcols << "\n";
362  for (ic = 0; ic < xcols; ic++) {
363  eq=0;
364  for(b=0; b<nb; b++) {
365  bsize = fBlockSize[b];
366  int64_t pos = fBlockPos[b];
367  for(r=0; r<bsize; r++) {
368  for(c=0; c<bsize; c++) {
369  z(eq+r,ic) += alpha*fStorage[pos+r+bsize*c]*x.GetVal((eq+c),ic);
370  }
371  }
372  eq+=bsize;
373  }
374  }
375  }
376 }
377 
378 
379 /***************/
380 /**** Zero ****/
381 template<class TVar>
383 {
384 
385  fStorage.Fill(0.,0);
386  this->fDecomposed = 0;
387 
388  return( 1 );
389 }
390 
391 
392 
393 /********************/
394 /*** Transpose () ***/
395 template<class TVar>
397 {
398  T->Resize( Dim(), Dim() );
399 
400  int64_t b, eq = 0, pos;
401  int bsize, r, c;
402  int64_t nb = fBlockSize.NElements();
403  for ( b=0; b<nb; b++) {
404  pos= fBlockPos[b];
405  bsize = fBlockSize[b];
406  for(r=0; r<bsize; r++) {
407  for(c=0; c<bsize; c++) {
408  T->PutVal(eq+r,eq+c,fStorage[pos+c+r*bsize]);
409  }
410  }
411  eq += bsize;
412  }
413 }
414 
415 
416 /*****************/
417 /*** Decompose_LU ***/
418 //fElem[ fBand * (2*row + 1) + col ]
419 template<class TVar>
420 int TPZBlockDiagonal<TVar>::Decompose_LU(std::list<int64_t> &singular)
421 {
422  return Decompose_LU();
423 }
424 
425 template<class TVar>
427 {
428 
429  LOGPZ_DEBUG(logger, "TPZBlockDiagonal::Decompose_LU");
430 
431  if ( this->fDecomposed && this->fDecomposed == ELU) {
432  return ELU;
433  } else if(this->fDecomposed) {
434  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__,"TPZBlockDiagonal::Decompose_LU is already decomposed with other scheme");
435  }
436 
437  int64_t b,nb,pos;
438  int bsize;
439  nb = fBlockSize.NElements();
440  for(b=0;b<nb; b++) {
441  pos = fBlockPos[b];
442  bsize = fBlockSize[b];
443  if(!bsize) continue;
444 
445 #ifdef LOG4CXX
446  if (logger->isDebugEnabled())
447  {
448  std::stringstream mess;
449  mess << "TPZBlockDiagonal::Decompose_LU() - bsize = " << bsize << ", bsize*bsize = " << bsize*bsize;
450  LOGPZ_DEBUG(logger,mess.str());
451  }
452 #endif
453 
454  TPZFMatrix<TVar> temp(bsize,bsize,&fStorage[pos],bsize*bsize);
455  std::list<int64_t> singular;
456  temp.Decompose_LU(singular);
457  }
458  this->fDecomposed = ELU;
459  return 1;
460 }
461 
462 template<class TVar>
463 int
465 {
466  if(this->fDecomposed != ELU) {
467  TPZMatrix<TVar>::Error(__PRETTY_FUNCTION__,"TPZBlockDiagonal::Decompose_LU is decomposed with other scheme");
468  }
469 
470  int64_t b,nb,pos,bsize,eq=0;
471  nb = fBlockSize.NElements();
472  int64_t c, nc = B->Cols();
473  for(c=0; c<nc; c++) {
474  eq = 0;
475  for(b=0;b<nb; b++) {
476  pos = fBlockPos[b];
477  bsize = fBlockSize[b];
478  if(!bsize) continue;
479  TPZFMatrix<TVar> BTemp(bsize,1,&(B->operator()(eq,c)),bsize);
480  TVar *ptr = fStorage.begin()+pos;
481  TPZFMatrix<TVar>::Substitution(ptr,bsize,&BTemp);
482  eq+= bsize;
483  }
484  }
485  return 1;
486 }
487 
488 /************************** Private **************************/
489 
490 /*************/
491 /*** Clear ***/
492 template<class TVar>
494 {
495  fStorage.Resize(0);
496  fBlockPos.Resize(0);
497  fBlockSize.Resize(0);
498  this->fRow = 0;
499  this->fCol = 0;
500  this->fDecomposed = 0;
501  return( 1 );
502 }
503 
504 template<class TVar>
506 
507  cout << "Entering the main program\n";
508  cout.flush();
509  TPZFMatrix<TVar> ref(7,7,0.);
510  int r,c;
511  for(r=0; r<7; r++) {
512  for(c=0; c<7; c++) {
513  ref(r,c) = ((TVar)(float)(5+r*c+3*r));
514  }
515  ref(r,r) += (TVar)1000;
516  }
517  TPZVec<int> blocksize(3);
518  blocksize[0] = 2;
519  blocksize[1] = 4;
520  blocksize[2] = 1;
521  TPZBlockDiagonal bd1(blocksize,ref);
522  TPZBlockDiagonal bd2(bd1);
523  ref.Print("original matrix",std::cout);
524  bd1.Solve_LU(&ref);
525  bd1.Solve_LU(&ref);
526  ref.Print("after inverting the diagonal",std::cout);
527  TPZFMatrix<TVar> ref2;
528  bd2.Multiply(ref,ref2);
529  bd2.Multiply(ref2,ref);
530  ref.Print("restoring the original matrix",std::cout);
531  return 1;
532 
533 }
534 
535 template<class TVar>
536 void TPZBlockDiagonal<TVar>::Print(const char *msg, std::ostream &out, const MatrixOutputFormat format) const {
537 
538  if(format != EFormatted)
539  {
540  TPZMatrix<TVar>::Print(msg,out,format);
541  return;
542  }
543  out << "TPZBlockDiagonal matrix ";
544  if(msg) out << msg;
545  out << std::endl;
546 
547  int64_t nblock = fBlockSize.NElements();
548  out << "Number of blocks " << nblock << std::endl;
549  int64_t b,bsize,pos;
550  for(b=0; b<nblock; b++) {
551  bsize = fBlockSize[b];
552  out << "block number " << b << " size : " << bsize << std::endl;
553  int64_t r,c;
554  pos = fBlockPos[b];
555  for(c=0; c<bsize; c++) {
556  for(r=0; r<bsize ; r++) {
557  out << fStorage[pos+r+bsize*c] << ' ';
558  }
559  out << std::endl;
560  }
561  }
562 }
563 
567 template<class TVar>
569 {
570  if(!mat)
571  {
572  cout << "TPZBlockDiagonal::UpdateFrom" << " called with zero argument\n";
573  return;
574  }
575  this->fDecomposed = ENoDecompose;
576  int64_t nblock = fBlockSize.NElements();
577  int64_t b,bsize,pos,firsteq = 0;
578  for(b=0; b<nblock; b++) {
579  bsize = fBlockSize[b];
580  // int r,c;
581  pos = fBlockPos[b];
582  TPZFMatrix<TVar> block(bsize,bsize,&fStorage[pos],bsize*bsize);
583  mat->GetSub(firsteq,firsteq,bsize,bsize,block);
584  firsteq += bsize;
585  }
586 }
587 
589 template<class TVar>
590 void TPZBlockDiagonal<TVar>::AutoFill(int64_t neq, int64_t jeq, int symmetric) {
591 
592  if (neq != jeq) {
593  DebugStop();
594  }
595  TPZStack<int> blsizes;
596  int64_t totalsize = 0;
597  while (totalsize < neq) {
598  int64_t blsize = (neq*rand())/RAND_MAX;
599  blsize = blsize < neq-totalsize ? blsize : neq-totalsize;
600  blsizes.Push(blsize);
601  totalsize += blsize;
602  }
603  Initialize(blsizes);
604  // Initialize the blocksizes!!
605  int64_t b, bsize, eq = 0, pos;
606  int64_t nb = fBlockSize.NElements(), r, c;
607  for ( b=0; b<nb; b++) {
608  pos= fBlockPos[b];
609  bsize = fBlockSize[b];
610  for(c=0; c<bsize; c++) {
611  float sum = 0.;
612  r=0;
613  if (symmetric == 1) {
614  for (r=0; r<c; r++) {
615  fStorage[pos+c+r*bsize]=fStorage[pos+r+c*bsize];
616  sum += fabs(fStorage[pos+r+c*bsize]);
617  }
618  }
619  for(; r<bsize; r++) {
620  float val = ((float)rand())/RAND_MAX;
621  fStorage[pos+c+r*bsize] = (TVar)(val);
622  if(c!= r) sum += fabs(val);
623  }
624  //totototo
625 // if(r==4)
626 // {
627 // std::cout << "sum " << sum << std::endl;
628 // }
629  if (fabs(fStorage[pos+c+c*bsize]) < sum) {
630  fStorage[pos+c+c*bsize] = (TVar)(sum + (float)1.);
631  }
632  }
633  eq += bsize;
634  }
635 }
636 
637 template<class TVar>
639  return Hash("TPZBlockDiagonal") ^ TPZMatrix<TVar>::ClassId() << 1;
640 }
641 
642 template class TPZBlockDiagonal<float>;
643 template class TPZBlockDiagonal<double>;
644 template class TPZBlockDiagonal<long double>;
645 
646 template class TPZBlockDiagonal<std::complex<float> >;
virtual void UpdateFrom(TPZAutoPointer< TPZMatrix< TVar > > mat) override
Updates the values of the matrix based on the values of the matrix.
void MultAdd(const TPZFMatrix< TVar > &x, const TPZFMatrix< TVar > &y, TPZFMatrix< TVar > &z, const TVar alpha=1., const TVar beta=0., const int opt=0) const override
Computes z = alpha * opt(this)*x + beta * y.
expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ fabs
Definition: tfadfunc.h:140
int Solve_LU(TPZFMatrix< TVar > *B, std::list< int64_t > &singular)
Solves the linear system using LU method .
Definition: pzmatrix.h:900
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.
virtual int Substitution(TPZFMatrix< TVar > *B) const override
Makes the backward and forward substitutions whether the matrix was LU decomposed.
static int Error(const char *msg, const char *msg2=0)
Returns error messages.
Definition: pzmatrix.cpp:1402
void Transpose(TPZMatrix< TVar > *const T) const override
It makes *T the transpose of current matrix.
MatrixOutputFormat
Defines output format.
Definition: pzmatrix.h:55
virtual int Resize(const int64_t newRows, const int64_t newCols)
Redimensions a matriz keeping the previous values.
Definition: pzmatrix.h:278
void BuildFromMatrix(TPZMatrix< TVar > &matrix)
Builds a block from matrix.
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
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
Contains TPZBlockDiagonal class which defines block diagonal matrices.
int Clear() override
Clean data matrix. Zeroes number of columns and rows.
virtual const TVar & GetVal(const int64_t, const int64_t) const
Get values without bounds checking This method is faster than "Get" if DEBUG is defined.
Definition: pzmatrix.cpp:56
~TPZBlockDiagonal()
Simple destructor.
virtual void Print(const char *message, std::ostream &out=std::cout, const MatrixOutputFormat format=EFormatted) const override
Prints current matrix data.
virtual void Resize(const int64_t newsize, const T &object)
Resizes the vector object reallocating the necessary storage, copying the existing objects to the new...
Definition: pzvec.h:373
int Put(const int64_t row, const int64_t col, const TVar &value) override
Put values with bounds checking if DEBUG variable is defined.
void Push(const T object)
Pushes a copy of the object on the stack.
Definition: pzstack.h:80
Contains TPZMatrixclass which implements full matrix (using column major representation).
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
#define LOGPZ_DEBUG(A, B)
Define log for debug info.
Definition: pzlog.h:87
Definition: pzmatrix.h:52
int64_t Rows() const
Returns number of rows.
Definition: pzmatrix.h:803
TPZVec< int > fBlockSize
Stores block sizes data.
Definition: pzblockdiag.h:142
void SetBlock(int64_t i, TPZFMatrix< TVar > &block)
Sets a block in the current matrix.
Definition: pzblockdiag.cpp:39
void AutoFill(int64_t dim, int64_t dimj, int symmetric)
static int Substitution(const TVar *ptr, int64_t rows, TPZFMatrix< TVar > *B)
Definition: pzfmatrix.cpp:1323
T * begin() const
Casting operator. Returns The fStore pointer.
Definition: pzvec.h:450
virtual int Decompose_LU() override
const TVar & GetVal(const int64_t row, const int64_t col) const override
This method don&#39;t make verification if the element exist. It is fast than Get.
void AddBlock(int64_t i, TPZFMatrix< TVar > &block)
Adds a block to current matrix.
Definition: pzblockdiag.cpp:25
Full matrix class. Matrix.
Definition: pzfmatrix.h:32
int Redim(const int64_t newRows, const int64_t newCols) override
Redimension a matrix and ZERO your elements.
Definition: pzfmatrix.h:616
A simple stack.
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
const TVar & Get(const int64_t row, const int64_t col) const override
Get value with bound checking.
int64_t fCol
Number of cols in matrix.
Definition: pzmatrix.h:781
void Initialize(const TPZVec< int > &blocksize)
Initializes current matrix based on blocksize.
Definition: pzblockdiag.cpp:64
Defines block diagonal matrices. Matrix.
Definition: pzblockdiag.h:21
int Zero() override
Zeroes all the elements of the matrix.
virtual int PutVal(const int64_t, const int64_t, const TVar &val)
Put values without bounds checking This method is faster than "Put" if DEBUG is defined.
Definition: pzmatrix.h:154
TPZVec< TVar > fStorage
Stores matrix data.
Definition: pzblockdiag.h:138
int PutVal(const int64_t row, const int64_t col, const TVar &value) override
This method don&#39;t make verification if the element exist. It is fast than Put.
TPZBlockDiagonal()
Simple constructor.
void GetBlock(int64_t i, TPZFMatrix< TVar > &block)
Gets a block from current matrix.
Definition: pzblockdiag.cpp:52
TVar & operator()(const int64_t row, const int64_t col)
static int main()
This method checks the working of the class.
virtual void Multiply(const TPZFMatrix< TVar > &A, TPZFMatrix< TVar > &res, int opt=0) const
It mutiplies itself by TPZMatrix<TVar>A putting the result in res.
Definition: pzmatrix.cpp:1916
TPZVec< int64_t > fBlockPos
Stores blocks data.
Definition: pzblockdiag.h:140
void Fill(const T &copy, const int64_t from=0, const int64_t numelem=-1)
Will fill the elements of the vector with a copy object.
Definition: pzvec.h:460
int64_t Cols() const
Returns number of cols.
Definition: pzmatrix.h:809
virtual void Print(std::ostream &out) const
Definition: pzmatrix.h:253
int64_t NElements() const
Returns the number of elements of the vector.
Definition: pzvec.h:190
int64_t Dim() const override
Returns the dimension of the matrix if the matrix is square.
Definition: pzblockdiag.h:61
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: pzfmatrix.h:566
int ClassId() const override
Define the class id associated with the class.
Definition: pzmatrix.h:957
Root matrix class (abstract). Matrix.
Definition: pzmatrix.h:60
This class implements a reference counter mechanism to administer a dynamically allocated object...
int ClassId() const override
Define the class id associated with the class.