NeoPZ
pzcondensedcompel.cpp
Go to the documentation of this file.
1 
6 #include "pzcondensedcompel.h"
7 #include "pzlog.h"
8 #include "pzstepsolver.h"
9 #include "pzelementgroup.h"
10 #include "pzcmesh.h"
11 
12 #ifdef LOG4CXX
13 static LoggerPtr logger(Logger::getLogger("pz.mesh.tpzcondensedcompel"));
14 #endif
15 
16 #ifdef USING_LAPACK
17 //#define USING_DGER
18 #ifdef USING_MKL
19 #include <mkl.h>
20 #elif MACOSX
21 #include <Accelerate/Accelerate.h>
22 #endif
23 #endif
24 
27 {
28  fKeepMatrix = keepmatrix;
29  if(!ref)
30  {
31  DebugStop();
32  }
33  fReferenceCompEl = ref;
34  fMesh = ref->Mesh();
35  TPZGeoEl *gel = ref->Reference();
36  if (gel) {
37  SetReference(gel->Index());
38  }
39  SetIndex(ref->Index());
40  fMesh->ElementVec()[fIndex] = this;
41  int ncon = ref->NConnects();
42  fIndexes.resize(ncon);
43  for (int ic=0; ic< ncon ; ic++) {
44  fIndexes[ic] = ic;
45  }
46  Resequence();
47 }
48 
49 
51 {
52  if (fMesh->ElementVec()[fIndex] == this) {
54  delete fReferenceCompEl;
55  }
56 }
57 
60 {
61  TPZCompEl *ref = fReferenceCompEl->Clone(mesh);
62  if(!ref)
63  {
64  DebugStop();
65  }
66  fReferenceCompEl = ref;
67  fMesh = ref->Mesh();
68  fKeepMatrix = copy.fKeepMatrix;
69  SetReference(ref->Reference()->Index());
70  SetIndex(ref->Index());
71  fMesh->ElementVec()[fIndex] = this;
72  int ncon = ref->NConnects();
73  fIndexes.resize(ncon);
74  for (int i=0; i<ncon ; ++i) {
75  fIndexes[i] = i;
76  }
77 }
78 
79 
82 {
83  int64_t myindex = fIndex;
84  fMesh->ElementVec()[myindex] = 0;
85  TPZCompEl *ReferenceEl = fReferenceCompEl;
86  int ncon = NConnects();
87  for (int ic=0; ic<ncon ; ic++) {
88  Connect(ic).SetCondensed(false);
89  }
90  TPZGeoEl *gel = Reference();
91  if (gel) {
92  gel->ResetReference();
93  }
94  delete this;
95  ReferenceEl->Mesh()->ElementVec()[myindex] = ReferenceEl;
96  if (gel) {
97  gel->SetReference(ReferenceEl);
98  }
99 }
100 
106 void TPZCondensedCompEl::SetConnectIndex(int inode, int64_t index)
107 {
108  LOGPZ_ERROR(logger,"SetConnectIndex should never be called")
109  DebugStop();
110 }
111 
124  std::map<int64_t,int64_t> & gl2lcConMap,
125  std::map<int64_t,int64_t> & gl2lcElMap) const
126 {
127  TPZCompEl *cel = fReferenceCompEl->ClonePatchEl(mesh,gl2lcConMap,gl2lcElMap);
128  TPZCondensedCompEl *result = new TPZCondensedCompEl(cel);
129  result->fIndexes = fIndexes;
130  result->fCondensed = fCondensed;
131  return result;
132 }
133 
142  TPZSolVec &sol, TPZGradSolVec &dsol,TPZFMatrix<REAL> &axes)
143 {
144  fReferenceCompEl->ComputeSolution(qsi,sol,dsol,axes);
145 }
146 
148 {
150 }
151 
165  TPZVec<REAL> &normal,
166  TPZSolVec &leftsol, TPZGradSolVec &dleftsol,TPZFMatrix<REAL> &leftaxes,
167  TPZSolVec &rightsol, TPZGradSolVec &drightsol,TPZFMatrix<REAL> &rightaxes)
168 {
169  fReferenceCompEl->ComputeSolution(qsi,normal,leftsol,dleftsol,leftaxes,rightsol,drightsol,rightaxes);
170 }
171 
182  const TPZFMatrix<REAL> &axes, TPZSolVec &sol, TPZGradSolVec &dsol)
183 {
184  fReferenceCompEl->ComputeSolution(qsi,phi,dphix,axes,sol,dsol);
185 }
186 
188 {
189  TPZStack<int> condensed;
190  TPZStack<int> notcondensed;
191  int nint=0,next=0;
192  std::set<int64_t> depreceive;
193  int ncon = NConnects();
194  for (int ic=0; ic<ncon; ic++) {
195  TPZConnect &c = Connect(ic);
197  while (dep) {
198  depreceive.insert(dep->fDepConnectIndex);
199  dep = dep->fNext;
200  }
201  }
202  for (int i=0; i<ncon ; ++i) {
203  TPZConnect &c = Connect(i);
204  int64_t cindex = ConnectIndex(i);
205  if(c.NElConnected() == 1 && c.HasDependency() == 0 && depreceive.find(cindex) == depreceive.end())
206  {
207  c.SetCondensed(true);
208  }
209  if (c.IsCondensed()) {
210  if(c.HasDependency())
211  {
212  std::cout << "Not Implemented yet\n";
213  DebugStop();
214  }
215  condensed.Push(i);
216  nint += c.NDof();
217  }
218  else
219  {
220  notcondensed.Push(i);
221  next += c.NDof();
222  }
223  }
224  fNumInternalEqs = nint;
225  fNumTotalEqs = nint+next;
226 
227  int ncond = condensed.size();
228  for (int i=0; i<ncond; ++i) {
229  fIndexes[i] = condensed[i];
230  }
231  for (int i=0; i<notcondensed.size(); ++i) {
232  fIndexes[i+ncond] = notcondensed[i];
233  }
234  //TPZAutoPointer<TPZMatrix<STATE> > k00 = new TPZFMatrix<STATE>(nint,nint,0.);
235  if (fKeepMatrix == false) {
236  nint = 0;
237  }
238  TPZAutoPointer<TPZMatrix<STATE> > k00 = new TPZFMatrix<STATE>(nint, nint, 0.);
239  //TPZStepSolver<STATE> *step = new TPZStepSolver<STATE>(k00);
241  if(0)
242  {
243  TPZAutoPointer<TPZMatrix<STATE> > mat2 = k00->Clone();
244  TPZStepSolver<STATE> *gmrs = new TPZStepSolver<STATE>(mat2);
245  step->SetReferenceMatrix(mat2);
246  gmrs->SetGMRES(20, 20, *step, 1.e-20, 0);
247  TPZAutoPointer<TPZMatrixSolver<STATE> > autogmres = gmrs;
248  }
249  step->SetDirect(ELDLt);
250  TPZAutoPointer<TPZMatrixSolver<STATE> > autostep = step;
251 
252 
253 
254  fCondensed.SetSolver(autostep);
255  if(fKeepMatrix == true)
256  {
257 // fCondensed.Redim(nint+next,nint);
259  }
260 }
261 
263 void TPZCondensedCompEl::TPZCondensedCompEl::Assemble()
264 {
267 
268  fCondensed.Zero();
269  TPZElementMatrix ek,ef;
270 
271  fReferenceCompEl->CalcStiff(ek,ef);
274 
275  int64_t dim = ek.fMat.Rows();
276  for (int64_t i=0; i<dim ; ++i) {
277  for (int64_t j=0; j<dim ; ++j) {
278  fCondensed(i,j) = ek.fMat(i,j);
279  }
280  }
281 
282  fCondensed.SetF(ef.fMat);
284 }
285 
286 
293 {
294  if(fKeepMatrix == false)
295  {
296  fKeepMatrix = true;
299  fKeepMatrix = false;
300  }
301  fReferenceCompEl->CalcStiff(ek,ef);
302 #ifdef LOG4CXX
303  if (logger->isDebugEnabled()) {
304  std::stringstream sout;
305  Print(sout);
306  sout << "Connect indices of element stiffness" << ek.fConnect << std::endl;
307  ek.fMat.Print("EKOrig = ",sout,EMathematicaInput);
308  LOGPZ_DEBUG(logger, sout.str())
309  }
310 #endif
313 #ifdef LOG4CXX
314  if (logger->isDebugEnabled()) {
315  std::stringstream sout;
316  sout << "Element connects\n";
317  int nc = NConnects();
318  for (int ic=0; ic<nc; ic++) {
319  sout << "ic = " << ic << ' ' << " index " << ConnectIndex(ic) << ' ';
320  Connect(ic).Print(*Mesh(),sout);
321  }
322  sout << "Permutations " << fIndexes << std::endl;
323  sout << "Connect indices " << ek.fConnect << std::endl;
324 
325  ek.fMat.Print("EKPermute = ",sout,EMathematicaInput);
326  ef.fMat.Print("EFPermute = ",sout,EMathematicaInput);
327  LOGPZ_DEBUG(logger, sout.str())
328  }
329 #endif
330 
331 // ek.fMat.Print("ek = ",std::cout,EMathematicaInput);
332 // ef.fMat.Print("ef = ",std::cout,EMathematicaInput);
333 
334 #ifdef USING_DGER
335 
336  // Initialization TPZMatRed structure
337  fCondensed.Zero();
338  {
339  TPZFMatrix<STATE> * K00_temp = dynamic_cast<TPZFMatrix<STATE> * >(fCondensed.K00().operator->());
340  K00_temp->InitializePivot();
341  }
342  int64_t dim0 = fCondensed.Dim0();
343  int64_t dim1 = fCondensed.Dim1();
344  int64_t rows = ek.fMat.Rows();
345  int64_t cols = ek.fMat.Cols()+ef.fMat.Cols();
346 
347  TPZFMatrix<STATE> KF(rows,cols); // Local object
348  for(int64_t i=0; i<rows;i++)
349  for (int64_t j=0; j<cols; j++)
350  {
351  if (j<rows)
352  KF(i,j) = ek.fMat(i,j);
353  else
354  KF(i,j) = ef.fMat(i,j-rows);
355  }
356 
357  for (int64_t i=0; i<dim1; i++)
358  for (int64_t j=0; j<dim0; j++)
359  fCondensed.K10().operator()(i,j)=KF(i+dim0,j);
360 
361  for (int64_t i=0; i<dim1; i++)
362  for (int64_t j=0; j<dim1; j++)
363  fCondensed(i+dim0,j+dim0)=KF(i+dim0,j+dim0);
364 
365  fCondensed.SetF(ef.fMat);
366 
367  // LDLt Decomposition
368  for (int64_t i=0; i<rows-dim1; i++)
369  {
370  for(int64_t j=i+1;j<cols;j++)
371  {
372  if (j<rows)
373  {
374  KF(j,i)/=KF(i,i);
375  KF(i,j)/=KF(i,i);
376  }
377  else
378  KF(i,j)/=KF(i,i);
379  }
380 
381 #ifdef STATEdouble
382  cblas_dger (CblasColMajor, rows-i-1, cols-i-1,
383  -KF(i,i), &KF(i+1,i), 1,
384  &KF(i,i+1), rows, &KF(i+1,i+1), rows);
385 #else
386  DebugStop();
387 #endif
388  }
389 
390  for (int64_t i=dim0; i< rows; i++)
391  {
392  ef.fMat(i,0) = KF.GetVal(i,fCondensed.Rows());
393  for (int64_t j=dim0; j< fCondensed.Rows(); j++)
394  {
395  fCondensed(i,j) = KF.GetVal(i,j);
396  }
397  }
398 
400  for (int64_t i=0; i<dim0; i++)
401  for (int64_t j=0; j<dim0; j++)
402  K00->operator()(i, j) = KF(i,j);
403 
405 
406  // if(0){// Implementation not complete
407  //
408  // for (int64_t i=0; i<dim0; i++){ // Substituindo valores obtidos para K01 usando o BLAS
409  // for (int64_t j=0; j<dim1; j++){
410  // fCondensed.K01().operator()(i,j) = KF(i,j+dim0)/KF(i,i);
411  // }
412  //
413  // }
414  //
415  // double scale = 1.0;
416  // cblas_dtrsm(CblasColMajor,CblasLeft,CblasUpper,CblasNoTrans,CblasUnit,dim0,dim1,scale,&fCondensed.K00().operator->()->operator()(0,0),dim0,&fCondensed.K01()(0,0),dim0);
417  // }else{
420  // }
421  fCondensed.SetF(ef.fMat);
424  fCondensed.SetReduced();// Directive that instructs the object to behave as reduced matrix.
425 
426  int64_t dim = ek.fMat.Rows();
427  for (int64_t i = dim0; i<dim; i++) {
428  ef.fMat(i,0) = fCondensed.F1()(i-dim0,0);
429  for (int64_t j = dim0; j<dim; j++) {
430  ek.fMat(i,j) = fCondensed.K11()(i-dim0,j-dim0);
431  }
432  }
433 
434 #ifdef LOG4CXX
435  if(logger->isDebugEnabled())
436  {
437  std::stringstream sout;
438  fCondensed.Print("Condensed = ",sout);
439  ek.fMat.Print("Rigidez",sout);
440  ef.fMat.Print("Carga",sout);
441  LOGPZ_DEBUG(logger, sout.str())
442  }
443 #endif
444 
445 #else
446 
447  fCondensed.Zero();
448 
449  int64_t dim = ek.fMat.Rows();
450  for (int64_t i=0; i<dim ; ++i) {
451  for (int64_t j=0; j<dim ; ++j) {
452  fCondensed(i,j) = ek.fMat(i,j);
453  }
454  }
455 
456  fCondensed.SetF(ef.fMat);
457 
458  int64_t dim1 = fCondensed.Dim1();
459  TPZFNMatrix<200,STATE> K11(dim1,dim1),F1(dim1,ef.fMat.Cols());
460  //const TPZFMatrix<REAL> &k11 = fCondensed.K11Red();
461 
462  fCondensed.K11Reduced(K11, F1);
463 
464 #ifdef PZDEBUG
465  {
466  REAL normk01 = Norm(fCondensed.K01());
467  REAL normf0 = Norm(fCondensed.F0());
468  if(std::isnan(normk01) || std::isnan(normf0))
469  {
470  Print();
471  DebugStop();
472  }
473  }
474 #endif
475 #ifdef LOG4CXX
476  if(logger->isDebugEnabled())
477  {
478  std::stringstream sout;
479  sout << "Index = " << Index() << std::endl;
480  fCondensed.Print("Reduced = ",sout,EMathematicaInput);
481  LOGPZ_DEBUG(logger, sout.str())
482  }
483 #endif
484 
486 
487  //const TPZFMatrix<REAL> &f1 = fCondensed.F1Red();
488  int64_t dim0 = dim-K11.Rows();
489  for (int64_t i=dim0; i<dim; i++) {
490  ef.fMat(i,0) = F1.GetVal(i-dim0,0);
491  for (int64_t j=dim0; j<dim; j++) {
492  ek.fMat(i,j) = K11.GetVal(i-dim0,j-dim0);
493  }
494  }
495 
496 
497 #endif
498 
499 #ifdef USING_DGER
500 #ifdef USING_LAPACK
501  TPZFMatrix<STATE> * K00_temp = dynamic_cast<TPZFMatrix<STATE> * >(fCondensed.K00().operator->());
502  K00_temp->InitializePivot();
503 #endif
504 #endif
505 
506 
507 #ifdef LOG4CXX
508  if(logger->isDebugEnabled())
509  {
510  std::stringstream sout;
511  int nc = NConnects();
512  for (int ic=0; ic<nc; ic++) {
513  sout << "ic = " << ic << ' ';
514  Connect(ic).Print(*Mesh(),sout);
515  }
516  ek.fMat.Print("EK11Reduced",sout,EMathematicaInput);
517  ef.fMat.Print("EF11Reduced",sout,EMathematicaInput);
518  LOGPZ_DEBUG(logger, sout.str())
519  }
520 #endif
521 
522 // fCondensed.Print("Condensed = ",std::cout);
523 // fCondensed.K11().Print("kbar = ",std::cout,EMathematicaInput);
524 // fCondensed.F1().Print("fbar = ",std::cout,EMathematicaInput);
525 // fCondensed.K01().Print("k01 = ",std::cout,EMathematicaInput);
526 // fCondensed.F0().Print("f0 = ",std::cout,EMathematicaInput);
527 
528  if (fKeepMatrix == false) {
529  fCondensed.K00()->Redim(0, 0);
530  fCondensed.K10().Redim(0,0);
531  fCondensed.K11().Redim(0,0);
532  }
533 }
534 
535 
541 {
542  // we need the stiffness matrix computed to compute the residual
543  if (fKeepMatrix == false) {
544  DebugStop();
545  fKeepMatrix = true;
546  fKeepMatrix = false;
547  }
550  fCondensed.SetF(ef.fMat);
551  //const TPZFMatrix<REAL> &f1 = fCondensed.F1Red();
553  fCondensed.F1Red(f1);
554  int64_t dim1 = f1.Rows();
555  int64_t dim = ef.fMat.Rows();
556  int64_t dim0 = dim-dim1;
557  for (int64_t i= dim0; i<dim; i++) {
558  ef.fMat(i,0) = f1.GetVal(i-dim0,0);
559  }
560  if (fKeepMatrix == false) {
561  fCondensed.Redim(0,0);
562  }
563 }
564 
566 bool TPZCondensedCompEl::HasMaterial(const std::set<int> &materialids) const {
567  if(!fReferenceCompEl){
568  return false;
569  }
570  bool has_material_Q = fReferenceCompEl->HasMaterial(materialids);
571  return has_material_Q;
572 }
573 
578 void TPZCondensedCompEl::Print(std::ostream &out) const
579 {
580  out << "Output for a condensed element\n";
581  TPZCompEl::Print(out);
582 
583  TPZElementGroup *eg = dynamic_cast<TPZElementGroup *>(fReferenceCompEl);
584  if(eg)
585  {
586  out << "Index of grouped elements: ";
587  int nel = eg->GetElGroup().size();
588  for(int i=0; i<nel-1; i++){
589  out << eg->GetElGroup()[i]->Index() <<", ";
590  }
591  out << eg->GetElGroup()[nel-1]->Index() <<std::endl;
592  out << "Connect indexes of the contained elements \n";
593  for(int i=0; i<nel; i++){
594  TPZCompEl *cel = eg->GetElGroup()[i];
595  TPZGeoEl *gel = cel->Reference();
596  out << "cel index " << cel->Index() << " cindex ";
597  int nc = cel->NConnects();
598  for (int ic=0; ic<nc; ic++) {
599  out << cel->ConnectIndex(ic) << " ";
600  }
601  if (gel) {
602  out << "\ngelindex " << gel->Index() << " ";
603  out << "matid " << gel->MaterialId();
604  TPZManVector<REAL,3> xi(gel->Dimension()), xco(3);
605  gel->CenterPoint(gel->NSides()-1, xi);
606  gel->X(xi,xco);
607  out << " center x " << xco;
608  }
609  out << std::endl;
610  }
611  }
612  else
613  {
614  if (fReferenceCompEl) {
615  fReferenceCompEl->Print(out);
616  }
617  else
618  {
619  DebugStop();
620  }
621  }
622  out << "Internal index resequencing: " << fIndexes << std::endl;
623 // fCondensed.Print("Condensed matrix",out,EMathematicaInput);
624 }
625 
626 
633 {
634 // if (fKeepMatrix == false) {
635 // fKeepMatrix = true;
636 // fCondensed.K00()->Redim(fNumInternalEqs, fNumInternalEqs);
637 // fCondensed.Redim(fNumTotalEqs, fNumInternalEqs);
638 // TPZElementMatrix ek,ef;
639 // CalcStiff(ek, ef);
640 // fKeepMatrix = false;
641 // }
642  // initialize the solution of the constrained connects
644 
645  // if the matrix has not been condensed then nothing to do
646  if (fCondensed.Rows() != fCondensed.Dim1())
647  {
648  return;
649  }
650  // compute the solution of the internal equations
651  int dim0=0, dim1=0;
652  int nc = NConnects(),nc0 = 0, nc1 = 0;
653  int ic;
654  for (ic=0; ic<nc ; ic++) {
655  int64_t connectindex = ConnectIndex(ic);
656  TPZConnect &con = Connect(ic);
657  int sz = con.NShape()*con.NState();
658  if (con.IsCondensed()) {
659 #ifdef PZDEBUG
660  if (dim1) {
661  DebugStop();
662  }
663 #endif
664  dim0 += sz;
665  nc0++;
666  }
667  else
668  {
669  dim1 += sz;
670  nc1++;
671  }
672  }
673  //TPZBlock<REAL> &bl = Mesh()->Block();
674  TPZBlock<STATE> &bl = Mesh()->Block();
675  int64_t count = 0;
676  //TPZFMatrix<REAL> u1(dim1,1,0.);
677  TPZFMatrix<STATE> u1(dim1,1,0.);
678  //TPZFMatrix<REAL> elsol(dim0+dim1,1,0.);
679  TPZFMatrix<STATE> elsol(dim0+dim1,1,0.);
680  for (ic=nc0; ic<nc ; ic++) {
681  TPZConnect &c = Connect(ic);
682  int64_t seqnum = c.SequenceNumber();
683  int blsize = bl.Size(seqnum);
684  for (int ibl=0; ibl<blsize; ibl++) {
685  u1(count++,0) = bl(seqnum,0,ibl,0);
686  }
687  }
688 #ifdef LOG4CXX
689  if (logger->isDebugEnabled()) {
690  std::stringstream sout;
691  sout << "Computing UGlobal Index " << Index();
692  sout << " Norm fK01 " << Norm(fCondensed.K01()) << std::endl;
693  TPZVec<STATE> u1vec(dim1);
694  for(int i=0; i<u1vec.size(); i++) u1vec[i] = u1(i,0);
695  sout << "u1 " << u1vec;
696  LOGPZ_DEBUG(logger, sout.str())
697  }
698 #endif
699  fCondensed.UGlobal(u1, elsol);
700  count = 0;
701  for (ic=0; ic<nc0 ; ic++) {
702  TPZConnect &c = Connect(ic);
703  int64_t seqnum = c.SequenceNumber();
704  int blsize = bl.Size(seqnum);
705  for (int ibl=0; ibl<blsize; ibl++) {
706  bl(seqnum,0,ibl,0) = elsol(count++,0);
707  }
708  }
709 #ifdef LOG4CXX
710  if (logger->isDebugEnabled()) {
711  std::stringstream sout;
712  sout << "After Computing UGlobal Index " << Index() ;
713  sout << " Norm fK01 " << Norm(fCondensed.K01()) << std::endl;
714  TPZVec<STATE> u1vec(dim1+dim0);
715  for(int i=0; i<u1vec.size(); i++) u1vec[i] = elsol(i,0);
716  sout << "elsol " << u1vec;
717  LOGPZ_DEBUG(logger, sout.str())
718  }
719 #endif
720 // if (fKeepMatrix == false) {
721 // fCondensed.Redim(0,0);
722 // fCondensed.K00()->Redim(0, 0);
723 // }
724 
725 }
726 
728 void TPZCondensedCompEl::BuildCornerConnectList(std::set<int64_t> &connectindexes) const
729 {
730  int nc = NConnects();
731  std::set<int64_t> refconn;
733  for (int ic=0; ic<nc ; ic++) {
734  TPZConnect &c = Connect(ic);
735  if (!c.IsCondensed() || !c.HasDependency()) {
736  int64_t index = ConnectIndex(ic);
737  if (refconn.find(index) != refconn.end()) {
738  connectindexes.insert(index);
739  }
740  }
741  }
742 }
743 
745  return Hash("TPZCondensedCompEl") ^ TPZCompEl::ClassId() << 1;
746 }
virtual void ComputeSolution(TPZVec< REAL > &qsi, TPZSolVec &sol, TPZGradSolVec &dsol, TPZFMatrix< REAL > &axes) override
Computes solution and its derivatives in the local coordinate qsi.
TPZStack< TPZCompEl *, 5 > GetElGroup()
virtual void BuildCornerConnectList(std::set< int64_t > &connectindexes) const =0
adds the connect indexes associated with base shape functions to the set
Represents a set of shape functions associated with a computational element/side. Computational Eleme...
Definition: pzconnect.h:30
Definition: pzmatrix.h:52
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.
TSideMatrix & K10()
Definition: pzmatred.h:120
int MaterialId() const
Returns the material index of the element.
Definition: pzgeoel.h:250
void SetReduced()
changes the declared dimension of the matrix to fDim1
Definition: pzmatred.h:70
int64_t Dim1()
Definition: pzmatred.h:145
virtual void resize(const int64_t newsize)
Definition: pzvec.h:213
void F1Red(TPZFMatrix< TVar > &F1)
Computes the reduced version of the right hand side .
Definition: pzmatred.cpp:171
TPZStack< int64_t > fConnect
Vector of pointers to TPZConnect objects.
Definition: pzelmat.h:39
void UGlobal(const TPZFMatrix< TVar > &U1, TPZFMatrix< TVar > &result)
Computes the complete vector based on the solution u1.
Definition: pzmatred.cpp:342
TPZFMatrix< TVar > & K11()
Definition: pzmatred.h:125
TPZFNMatrix< 1000, STATE > fMat
Pointer to a blocked matrix object.
Definition: pzelmat.h:41
Contains the declaration of the TPZElementGroup class, which implements an computational element whic...
Defines step solvers class. Solver.
Definition: pzmganalysis.h:17
void SetF(const TPZFMatrix< TVar > &F)
Copies the F vector in the internal data structure.
Definition: pzmatred.cpp:140
virtual void CalcStiff(TPZElementMatrix &ek, TPZElementMatrix &ef)
Computes the element stifness matrix and right hand side.
Definition: pzcompel.h:794
void SetIsDecomposed(int val)
Sets current matrix to decomposed state.
Definition: pzmatrix.h:410
virtual int64_t ConnectIndex(int i) const override
Returns the index of the ith connectivity of the element.
int64_t Dim0()
Definition: pzmatred.h:140
Class which groups elements to characterize dense matrices.
void SetF0IsComputed(bool directive)
Sets F0 as computed.
Definition: pzmatred.h:107
int NDof(TPZCompMesh &mesh)
Number of degrees of freedom associated with the object.
Definition: pzconnect.cpp:317
virtual int NSides() const =0
Returns the number of connectivities of the element.
virtual TPZCompEl * ClonePatchEl(TPZCompMesh &mesh, std::map< int64_t, int64_t > &gl2lcConMap, std::map< int64_t, int64_t > &gl2lcElMap) const =0
Method for creating a copy of the element in a patch mesh.
virtual TPZCompEl * ClonePatchEl(TPZCompMesh &mesh, std::map< int64_t, int64_t > &gl2lcConMap, std::map< int64_t, int64_t > &gl2lcElMap) const override
Method for creating a copy of the element in a patch mesh.
virtual void LoadSolution() override
Loads the solution within the internal data structure of the element.
void SetSolver(TPZAutoPointer< TPZMatrixSolver< TVar > > solver)
Definition: pzmatred.cpp:125
virtual int ClassId() const override
Define the class id associated with the class.
Definition: pzcompel.cpp:1129
virtual void BuildCornerConnectList(std::set< int64_t > &connectindexes) const override
adds the connect indexes associated with base shape functions to the set
TPZMatRed< STATE, TPZFMatrix< STATE > > fCondensed
void SetK01IsComputed(bool directive)
Sets K01 as computed.
Definition: pzmatred.h:99
virtual void CenterPoint(int side, TPZVec< REAL > &masscent) const =0
It returns the coordinates from the center of the side of the element in the element coordinate space...
void Print(const char *name=NULL, std::ostream &out=std::cout, const MatrixOutputFormat=EFormatted) const override
Prints the object data structure.
Definition: pzmatred.cpp:454
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
virtual int NConnects() const override
Returns the number of nodes of the element.
void Unwrap()
unwrap the condensed element from the computational element and delete the condensed element ...
virtual TPZCompEl * Clone(TPZCompMesh &mesh) const =0
Method for creating a copy of the element.
TPZFMatrix< TVar > & F0()
Definition: pzmatred.h:135
int64_t SequenceNumber() const
Returns the Sequence number of the connect object.
Definition: pzconnect.h:158
void Push(const T object)
Pushes a copy of the object on the stack.
Definition: pzstack.h:80
int64_t fDepConnectIndex
Definition: pzconnect.h:69
virtual void LoadSolution()
Loads the solution within the internal data structure of the element.
Definition: pzcompel.cpp:200
Structure to reference dependency.
Definition: pzconnect.h:67
TPZCompEl * fReferenceCompEl
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
int64_t Index() const
Returns the index of the element within the element vector of the mesh.
Definition: pzgeoel.h:730
virtual bool HasMaterial(const std::set< int > &materialids) const
Verifies if the material associated with the element is contained in the set.
Definition: pzcompel.cpp:968
TVar Norm(const TPZFMatrix< TVar > &A)
Returns the norm of the matrix A.
virtual int ClassId() const override
Define the class id associated with the class.
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
virtual void SetConnectIndex(int inode, int64_t index) override
Set the index i to node inode.
#define LOGPZ_DEBUG(A, B)
Define log for debug info.
Definition: pzlog.h:87
TPZCompMesh * Mesh() const
Return a pointer to the grid of the element.
Definition: pzcompel.cpp:288
void SetReference(int64_t referenceindex)
Definition: pzcompel.h:164
unsigned int NShape() const
Definition: pzconnect.h:151
int64_t Rows() const
Returns number of rows.
Definition: pzmatrix.h:803
void SetIndex(int64_t index)
Sets element index of the mesh fELementVec list.
Definition: pzcompel.cpp:577
virtual void SetReferenceMatrix(TPZAutoPointer< TPZMatrix< TVar > > matrix)
This method gives a preconditioner to share a matrix with the referring solver object.
Definition: pzsolve.h:132
int HasDependency() const
Returns whether exist dependecy information.
Definition: pzconnect.h:292
const TPZBlock< STATE > & Block() const
Access the block structure of the solution vector.
Definition: pzcmesh.h:213
virtual bool HasMaterial(const std::set< int > &materialids) const override
Verifies if the material associated with the element is contained in the set.
TPZDepend * fNext
Definition: pzconnect.h:71
Class which implements an element which condenses the internal connects.
int64_t Index() const
Returns element index of the mesh fELementVec list.
Definition: pzcompel.h:821
TPZManVector< int64_t, 27 > fIndexes
#define LOGPZ_ERROR(A, B)
Define log for errors (cout)
Definition: pzlog.h:93
virtual int Redim(const int64_t newRows, const int64_t newCols)
Redimensions the matrix reinitializing it with zero.
Definition: pzmatrix.h:289
virtual void ComputeSolution(TPZVec< REAL > &qsi, TPZMaterialData &data)
Definition: pzcompel.h:462
Contains declaration of TPZCompMesh class which is a repository for computational elements...
int Redim(const int64_t newRows, const int64_t newCols) override
Redimension a matrix and ZERO your elements.
Definition: pzfmatrix.h:616
int Redim(const int64_t dim, const int64_t dim00) override
Redim: Set the dimension of the complete matrix and reduced matrix.
Definition: pzmatred.cpp:481
virtual TPZConnect & Connect(int i) const
Returns a pointer to the ith node.
Definition: pzcompel.cpp:298
void SetCondensed(bool flag)
Set the connect as a condensed connect or not.
Definition: pzconnect.h:242
unsigned char NState() const
Number of state variables associated with the connect.
Definition: pzconnect.h:146
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
virtual int64_t ConnectIndex(int i) const =0
Returns the index of the ith connectivity of the element.
This class associates an element matrix with the coeficients of its contribution in the global stiffn...
Definition: pzelmat.h:30
virtual int NConnects() const =0
Returns the number of nodes of the element.
void PermuteGather(TPZVec< int64_t > &permute)
permute the order of the connects
Definition: pzelmat.cpp:502
virtual int Dimension() const =0
Returns the dimension of the element.
virtual void CalcResidual(TPZElementMatrix &ef)
Computes the element right hand side.
Definition: pzcompel.cpp:610
TPZAutoPointer< TPZMatrix< TVar > > K00()
Definition: pzmatred.h:112
TPZGeoEl * Reference() const
Return a pointer to the corresponding geometric element if such exists, return 0 otherwise.
Definition: pzcompel.cpp:1137
virtual void X(TPZVec< REAL > &qsi, TPZVec< REAL > &result) const =0
Return the coordinate in real space of the point coordinate in the master element space...
void ResetReference()
Reset the element referenced by the geometric element to NULL.
Definition: pzgeoel.h:431
virtual void CalcStiff(TPZElementMatrix &ek, TPZElementMatrix &ef) override
Computes the element stifness matrix and right hand side.
TPZDepend * FirstDepend()
Definition: pzconnect.h:296
void K11Reduced(TPZFMatrix< TVar > &K11, TPZFMatrix< TVar > &F1)
Computes the K11 reduced .
Definition: pzmatred.cpp:223
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
int Size(const int block_diagonal) const
Returns block dimension.
Definition: pzblock.h:171
TPZAdmChunkVector< TPZCompEl * > & ElementVec()
Returns a reference to the element pointers vector.
Definition: pzcmesh.h:183
virtual void Print(std::ostream &out=std::cout) const
Prints element data.
Definition: pzcompel.cpp:322
int NElConnected() const
Returns fNElConnected.
Definition: pzconnect.h:264
virtual void Print(std::ostream &out=std::cout) const override
Prints element data.
void SetReference(TPZCompEl *elp)
Make the current element reference to the computational element.
Definition: pzgeoel.h:746
Contains TPZStepSolver class which defines step solvers class.
int64_t Cols() const
Returns number of cols.
Definition: pzmatrix.h:809
virtual void Print(std::ostream &out) const
Definition: pzmatrix.h:253
TPZFMatrix< TVar > & F1()
Definition: pzmatred.h:130
void SetDirect(const DecomposeType decomp)
TPZCondensedCompEl(TPZCompEl *ref, bool keepmatrix=true)
void Print(const TPZCompMesh &mesh, std::ostream &out=std::cout)
Print the information for the connect element.
Definition: pzconnect.cpp:67
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
bool IsCondensed() const
Access method to return the indication whether the connect is condensed or not.
Definition: pzconnect.h:223
Contains the declaration of the TPZCondensedCompEl class, which implements an computational element w...
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
int64_t fIndex
Element index into mesh element vector.
Definition: pzcompel.h:67
virtual void CalcResidual(TPZElementMatrix &ef) override
Computes the element right hand side.
virtual int Zero() override
This method will zero all submatrices associated with this reducable matrix class.
Definition: pzmatred.cpp:505
virtual int SolveDirect(TPZFMatrix< TVar > &F, const DecomposeType dt, std::list< int64_t > &singular)
Solves the linear system using Direct methods.
Definition: pzmatrix.cpp:710
Non abstract class which implements full matrices with preallocated storage with (N+1) entries...
Definition: pzfmatrix.h:716
TPZCompMesh * fMesh
Computational mesh to which the element belongs.
Definition: pzcompel.h:64
void SetGMRES(const int64_t numiterations, const int numvectors, const TPZMatrixSolver< TVar > &pre, const REAL tol, const int64_t FromCurrent)
This class implements a reference counter mechanism to administer a dynamically allocated object...
TSideMatrix & K01()
Definition: pzmatred.h:116