NeoPZ
TPZFrontStructMatrix.cpp
Go to the documentation of this file.
1 
6 #include "pzstrmatrix.h"
7 #include "pzfstrmatrix.h"
8 #include "TPZFrontStructMatrix.h"
9 
10 #include "pzgmesh.h"
11 #include "pzcmesh.h"
12 #include "pzsubcmesh.h"
13 #include "pzconnect.h"
14 #include "pzadmchunk.h"
15 
16 #include "pzsmfrontalanal.h"
17 
18 
19 #include "pzanalysis.h"
20 #include "pzsolve.h"
21 #include "pzstepsolver.h"
22 #include "TPZFrontMatrix.h"
23 
24 #include "pzdxmesh.h"
25 #include <fstream>
26 
27 using namespace std;
28 
29 #include "pzelmat.h"
30 #include "pzbndcond.h"
31 
32 #include "pzlog.h"
33 
34 #ifdef LOG4CXX
35 static LoggerPtr logger(Logger::getLogger("pz.strmatrix.frontstructmatrix"));
36 static LoggerPtr loggerel(Logger::getLogger("pz.strmatrix.element"));
37 #endif
38 
39 
40 
41 template <class front>
43  int64_t ic;
44 
45  fMesh->ComputeNodElCon();
46 
47  for(ic=0; ic<fMesh->ConnectVec().NElements(); ic++) {
48  TPZConnect &cn = fMesh->ConnectVec()[ic];
49  if(cn.HasDependency() || cn.IsCondensed()) continue;
50  int64_t seqn = cn.SequenceNumber();
51  if(seqn < 0) continue;
52  int64_t firsteq = fMesh->Block().Position(seqn);
53  int64_t lasteq = firsteq+fMesh->Block().Size(seqn);
54  int64_t numactive = fEquationFilter.NumActive(firsteq, lasteq);
55  if(!numactive) continue;
56  if (numactive != lasteq-firsteq) {
57  DebugStop();
58  }
59  TPZManVector<int64_t> firstind(numactive),destindex(numactive);
60  for (int64_t i = 0; i<numactive; i++) {
61  firstind[i] = i;
62  destindex[i] = firsteq+i;
63  }
64  fEquationFilter.Filter(firstind, destindex);
65  for(int64_t ind=0;ind<destindex.size();ind++)
66  numelconnected[destindex[ind] ] = fMesh->ConnectVec()[ic].NElConnected();
67  }
68 }
69 
70 template<class front>
72 }
73 
74 template<class front>
76 }
77 
78 
79 template<class front>
81 
82 
83 
84 template<class front>
86 
87  /* TPZVec <int> numelconnected(fMesh->NEquations(),0);
88  // TPZFrontMatrix<TPZStackEqnStorage, TPZFrontNonSym> *mat = new TPZFrontMatrix<TPZStackEqnStorage, TPZFrontNonSym>(fMesh->NEquations());
89 
90  GetNumElConnected(numelconnected);
91  //mat->SetNumElConnected(numelconnected);
92  return mat;
93  */
94  return 0;
95 }
96 
97 template<class front>
99 
101  return result;
102 }
103 template<class front>
104 void TPZFrontStructMatrix<front>::OrderElement()//TPZVec<int> &elorder)
105 {
106  int64_t numelconnected = 0;
107  int64_t nconnect = fMesh->ConnectVec().NElements();
108  int64_t ic;
109  //firstelconnect contains the first element index in the elconnect vector
110  TPZVec<int64_t> firstelconnect(nconnect+1);
111  firstelconnect[0] = 0;
112  for(ic=0; ic<nconnect; ic++) {
113  numelconnected += fMesh->ConnectVec()[ic].NElConnected();
114  firstelconnect[ic+1] = firstelconnect[ic]+fMesh->ConnectVec()[ic].NElConnected();
115  }
116 
117 #ifdef LOG4CXX
118  if (logger->isDebugEnabled())
119  {
120  std::stringstream sout;
121  sout<<"numelconnected " << numelconnected << endl;
122  sout<< "firstelconnect "<< firstelconnect;
123  LOGPZ_DEBUG(logger,sout.str())
124  }
125 #endif
126  //cout << "numelconnected " << numelconnected << endl;
127  //cout << "firstelconnect ";
128  // for(ic=0; ic<nconnect; ic++) cout << firstelconnect[ic] << ' ';
129  TPZVec<int64_t> elconnect(numelconnected,-1);
130  int64_t el;
131  TPZCompEl *cel;
132  for(el=0; el<fMesh->ElementVec().NElements(); el++) {
133  cel = fMesh->ElementVec()[el];
134  if(!cel) continue;
135  TPZStack<int64_t> connectlist;
136  cel->BuildConnectList(connectlist);
137  int64_t nc = connectlist.NElements();
138  int64_t ic;
139  for(ic=0; ic<nc; ic++) {
140  int64_t cindex = connectlist[ic];
141  elconnect[firstelconnect[cindex]] = el;
142  firstelconnect[cindex]++;
143  }
144  }
145  // for(ic=0; ic<numelconnected; ic++) cout << elconnect[ic] << endl;
146 #ifdef LOG4CXX
147  if (logger->isDebugEnabled())
148  {
149  std::stringstream sout;
150  sout<< "elconnect "<< elconnect;
151  LOGPZ_DEBUG(logger,sout.str())
152  }
153 #endif
154  firstelconnect[0] = 0;
155  for(ic=0; ic<nconnect; ic++) {
156  firstelconnect[ic+1] = firstelconnect[ic]+fMesh->ConnectVec()[ic].NElConnected();
157  }
158  //cout << "elconnect\n";
159  // int no;
160  for(int64_t no=0; no< fMesh->ConnectVec().NElements(); no++) {
161 #ifdef LOG4CXX
162  if (logger->isDebugEnabled())
163  {
164  std::stringstream sout;
165  sout<< "Node index " << no << ' ' << " seq num " << fMesh->ConnectVec()[no].SequenceNumber() << ' ';
166  LOGPZ_DEBUG(logger,sout.str())
167  }
168 #endif
169  //cout << "no numero " << no << ' ' << " seq num " << fMesh->ConnectVec()[no].SequenceNumber() << ' ';
170  // for(ic=firstelconnect[no]; ic<firstelconnect[no+1];ic++) cout << elconnect[ic] << ' ';
171  //cout << endl;
172  }
173 
174 
176  fElementOrder.Fill(-1);
177  TPZVec<int> nodeorder(fMesh->ConnectVec().NElements(),-1);
178  firstelconnect[0] = 0;
179  for(ic=0; ic<nconnect; ic++) {
180  int64_t seqnum = fMesh->ConnectVec()[ic].SequenceNumber();
181  if(seqnum >= 0) nodeorder[seqnum] = ic;
182  }
183  // cout << "nodeorder ";
184  /*for(ic=0; ic<fMesh->ConnectVec().NElements(); ic++) cout << nodeorder[ic] << ' ';
185  cout << endl;
186  cout.flush();*/
187  int64_t seq;
188  int64_t elsequence = 0;
189  TPZVec<int64_t> elorderinv(fMesh->ElementVec().NElements(),-1);
190  for(seq=0; seq<nconnect; seq++) {
191  ic = nodeorder[seq];
192  if(ic == -1) continue;
193  int64_t firstind = firstelconnect[ic];
194  int64_t lastind = firstelconnect[ic+1];
195  int64_t ind;
196  for(ind=firstind; ind<lastind; ind++) {
197  el = elconnect[ind];
198  if(el == -1) {
199  continue;
200  }
201  if(elorderinv[el]==-1) elorderinv[el] = elsequence++;
202  }
203  }
204  // cout << "elorderinv ";
205  // for(seq=0;seq<fMesh->ElementVec().NElements();seq++) cout << elorderinv[seq] << ' ';
206  // cout << endl;
207  elsequence = 0;
208  for(seq=0;seq<fMesh->ElementVec().NElements();seq++) {
209  if(elorderinv[seq] == -1) continue;
210  fElementOrder[elorderinv[seq]] = seq;
211  }
212 #ifdef LOG4CXX
213  if (logger->isDebugEnabled())
214  {
215  std::stringstream sout;
216  sout << "element order " << fElementOrder << std::endl;
217  sout<< "elorderinv " << elorderinv << std::endl;
218  LOGPZ_DEBUG(logger,sout.str())
219  }
220 #endif
221  // cout << "elorder" << endl;
222  // for(ic=0; ic<fMesh->ElementVec().NElements(); ic++) cout << elorder[ic] << endl;
223 
224 }
225 
226 template<class front>
228 
229  int64_t neq = fEquationFilter.NActiveEquations();
230  TPZManVector <int> numelconnected(neq,0);
231  TPZFrontMatrix<STATE,TPZStackEqnStorage<STATE>, front> *mat = new TPZFrontMatrix<STATE,TPZStackEqnStorage<STATE>, front>(neq);//(fMesh->NEquations());
232 
233 // TPZFrontMatrix<STATE,TPZFileEqnStorage<STATE>, front> *mat = new TPZFrontMatrix<STATE,TPZFileEqnStorage<STATE>, front>(neq);
235  // if the frontal matrix is applied to a submesh, we assume there may be rigid body modes
236  TPZSubCompMesh *subcmesh = dynamic_cast<TPZSubCompMesh *> (fMesh);
237  if (subcmesh) {
238  int nrigid = subcmesh->NumberRigidBodyModes();
239  if (nrigid > 0) {
240  mat->GetFront().SetNumRigidBodyModes(nrigid);
241 
242  }
243  }
244 
245  GetNumElConnected(numelconnected);
246  mat->SetNumElConnected(numelconnected);
247 
248  OrderElement();
249 
250  Assemble(*mat,rhs,guiInterface);
251 
252 #ifdef LOG4CXX
253  if (logger->isDebugEnabled())
254  {
255  std::stringstream sout;
256  mat->FinishWriting();
257  mat->ReOpen();
258 // mat->Print("Frontal matrix", sout);
259  LOGPZ_DEBUG(logger,sout.str())
260  }
261 #endif
262 
263  mat->FinishWriting();
264  mat->ReOpen();
265  return mat;
266 }
267 
268 template<class front>
270 
271  int64_t iel;
272  int64_t numel = 0, nelem = fMesh->NElements();
274  TPZManVector<int64_t> destinationindex(0);
275  TPZManVector<int64_t> sourceindex(0);
276 
278 
279 
281  TPZVec<int> elorder(fMesh->NEquations(),0);
282 
283  OrderElement();
284 
285 
286  for(iel=0; iel < nelem; iel++) {
287 
288  if(guiInterface) if(guiInterface->AmIKilled()){
289  break;
290  }
291 
292  if(fElementOrder[iel] < 0) continue;
293  TPZCompEl *el = elementvec[fElementOrder[iel]];
294  if(!el) continue;
295  // int dim = el->NumNodes();
296 
297  //Builds elements stiffness matrix
298  el->CalcStiff(ek,ef);
299  ek.ComputeDestinationIndices();
300  FilterEquations(ek.fSourceIndex,ek.fDestinationIndex);
301  //ek.fMat->Print(out);
302  //ef.fMat->Print();
303  if(!f_quiet)
304  {
305  if(!(numel%20)) cout << endl << numel;
306  // if(!(numel%20)) cout << endl;
307  cout << '*';
308  cout.flush();
309  }
310  numel++;
311 
312  if(!el->HasDependency()) {
313  stiffness.AddKel(ek.fMat,ek.fSourceIndex,ek.fDestinationIndex);
314  rhs.AddFel(ef.fMat,ek.fSourceIndex, ek.fDestinationIndex);
315  }
316  else {
317  //ek.Print(*fMesh,cout);
318  ek.ApplyConstraints();
319  ef.ApplyConstraints();
320  stiffness.AddKel(ek.fConstrMat,ek.fSourceIndex,ek.fDestinationIndex);
321  rhs.AddFel(ef.fConstrMat,ek.fSourceIndex,ek.fDestinationIndex);
322  /*
323  if(ek.fConstrMat->Decompose_LU() != -1) {
324  el->ApplyConstraints(ek,ef);
325  ek.Print(*this,check);
326  check.flush();
327  }
328  */
329  }
330 
331  }//fim for iel
332  if(!f_quiet)
333  {
334  cout << endl;
335  }
336 }
337 
338 
339 template<class front>
341 
342  int64_t iel;
343  int64_t numel = 0, nelem = fMesh->NElements();
345 
347 
348 
350  TPZVec<int> elorder(fMesh->NEquations(),0);
351 
352  OrderElement();
353 
354  for(iel=0; iel < nelem; iel++) {
355 
356  int64_t elindex = fElementOrder[iel];
357  if(elindex < 0) continue;
358  TPZCompEl *el = elementvec[elindex];
359  if(!el) continue;
360  TPZMaterial * mat = el->Material();
361  if(mat)
362  {
363  int matid = mat->Id();
364  if(this->ShouldCompute(matid) == false)
365  {
366  continue;
367  }
368  }
369  else
370  {
371  }
372 
373  // int dim = el->NumNodes();
374 
375  //Builds elements stiffness matrix
376  el->CalcStiff(ek,ef);
377 
378 
379  std::cout<< " assemblando elemento frontal " << iel <<std::endl;
380 
381 #ifdef LOG4CXX
382  if (logger->isDebugEnabled())
383  {
384  std::stringstream sout;
385  ek.fMat.Print("Element stiffness Frontal",sout);
386  LOGPZ_DEBUG(logger,sout.str())
387  }
388 #endif
389 
390  AssembleElement(el, ek, ef, stiffness, rhs);
391 #ifdef LOG4CXX
392  if(loggerel->isDebugEnabled())
393  {
394  std::stringstream sout;
395  ek.fMat.Print("Element stiffness depois de assemblada ",sout);
396  LOGPZ_DEBUG(loggerel,sout.str())
397  }
398 #endif
399 
400  if(!f_quiet)
401  {
402  cout << '*';
403  if(!(numel%20)) {
404  cout << " " << (100*iel/nelem) << "% Elements assembled " << endl;
405  cout.flush();
406  }
407  }
408  numel++;
409 
410  }//fim for iel
411 
412 }
413 
414 //Verificar declaracao dos parametros !!!!!
415 template<class front>
417 
418 
419  if(!el->HasDependency()) {
420  //ek.fMat->Print("stiff has no constraint",test);
421  //ef.fMat->Print("rhs has no constraint",test);
422  //test.flush();
425 #ifdef LOG4CXX
426  if (logger->isDebugEnabled())
427  {
428  std::stringstream sout;
429  sout << "Element index " << el->Index() << " Unconstrained destination index " << ek.fDestinationIndex;
430  LOGPZ_DEBUG(logger,sout.str())
431  }
432 #endif
433  //ek.Print(*fMesh,cout);
434  stiffness.AddKel(ek.fMat,ek.fSourceIndex,ek.fDestinationIndex);
435  rhs.AddFel(ef.fMat,ek.fSourceIndex, ek.fDestinationIndex); // ??????????? Erro
436  }
437  else
438  {
439  ek.ApplyConstraints();
440  ef.ApplyConstraints();
443 #ifdef LOG4CXX
444  if (logger->isDebugEnabled())
445  {
446  std::stringstream sout;
447  sout << "Element index " << el->Index() << " Constrained destination index " << ek.fDestinationIndex;
448  LOGPZ_DEBUG(logger,sout.str())
449  }
450 #endif
451  stiffness.AddKel(ek.fConstrMat,ek.fSourceIndex,ek.fDestinationIndex);
453  }
454 }
455 
456 
457 
461 template<class front>
463 {
464  this->f_quiet = quiet;
465 }
466 
467 #ifndef STATE_COMPLEX
468 #include "pzmat2dlin.h"
469 
470 template<class front>
472  int refine = 5;
473  int order = 2;
474 
475  TPZGeoMesh gmesh;
476  TPZCompMesh cmesh(&gmesh);
477  double coordstore[4][3] = {{0.,0.,0.},{1.,0.,0.},{1.,1.,0.},
478  {0.,1.,0.}};
479 
480  int i,j;
481  TPZVec<REAL> coord(3,0.);
482  for(i=0; i<4; i++) {
483  // initializar as coordenadas do no em um vetor
484  for (j=0; j<3; j++) coord[j] = coordstore[i][j];
485 
486  // identificar um espa�o no vetor onde podemos armazenar
487  // este vetor
488 
489  // initializar os dados do n�
490  gmesh.NodeVec ()[i].Initialize (i,coord,gmesh);
491  }
492  int el;
493  //TPZGeoEl *gel;
494  for(el=0; el<1; el++) {
495 
496  // initializar os indices dos n�s
497  TPZVec<int64_t> indices(4);
498  for(i=0; i<4; i++) indices[i] = i;
499  // O proprio construtor vai inserir o elemento na malha
500  int64_t index;
501  /*gel = */gmesh.CreateGeoElement(EQuadrilateral, indices, 1, index);
502  }
503  gmesh.BuildConnectivity ();
504 
505  TPZVec<TPZGeoEl *> subel;
506  //gel->Divide(subel);
507 
508 
509 
510  cout << "Refinement ";
511  cin >> refine;
512  cout << endl;
513  DebugStop();
514 // UniformRefine(refine,gmesh);
515 
516 
517 
518  TPZMat2dLin *meumat = new TPZMat2dLin(1);
519  TPZFMatrix<STATE> xk(1,1,1.),xc(1,2,0.),xf(1,1,1.);
520  meumat->SetMaterial (xk,xc,xf);
521  TPZMaterial * meumatptr(meumat);
522  cmesh.InsertMaterialObject(meumatptr);
523 
524  TPZFMatrix<STATE> val1(1,1,0.),val2(1,1,0.);
525  TPZMaterial * bnd = meumat->CreateBC (meumatptr,-4,0,val1,val2);
526  cmesh.InsertMaterialObject(bnd);
527 
528 
529 
530  cout << "Interpolation order ";
531  cin >> order;
532  cout << endl;
533 
534  // TPZCompEl::gOrder = order;
535  cmesh.SetDefaultOrder(order);
536 
537  cmesh.AutoBuild();
538  // cmesh.AdjustBoundaryElements();
539  cmesh.InitializeBlock();
540 
541  ofstream output("outputNon.dat");
542  // ofstream output2("outputNon.dat");
543  cmesh.Print(output);
544  TPZAnalysis an(&cmesh,true,output);
545  // TPZAnalysis an2(&cmesh,output);
546 
547  TPZVec<int> numelconnected(cmesh.NEquations(),0);
548  int64_t ic;
549  cout << "Numero de Equacoes -> " << cmesh.NEquations() << endl;
550  cout.flush();
551 
552  ofstream out("cmeshBlock_out.txt");
553  // cmesh.Print(out);
554  // cmesh.Block().Print("Block",out);
555  for(ic=0; ic<cmesh.ConnectVec().NElements(); ic++) {
556  TPZConnect &cn = cmesh.ConnectVec()[ic];
557  if(cn.HasDependency()) continue;
558  int64_t seqn = cn.SequenceNumber();
559  if(seqn < 0) continue;
560  int64_t firsteq = cmesh.Block().Position(seqn);
561  int64_t lasteq = firsteq+cmesh.Block().Size(seqn);
562  int64_t ind;
563  int temp = cmesh.ConnectVec()[ic].NElConnected();
564  for(ind=firsteq;ind<lasteq;ind++) {
565  numelconnected[ind] = temp;//cmesh.ConnectVec()[ic].NElConnected();
566  }
567  }
568  // << "nequations " << numelconnected.NElements();
569  // for(ic=0;ic<numelconnected.NElements(); ic++) cout << numelconnected[ic] <<' ';
570  // cout << endl;
571  // cout.flush();
572 
573  // TPZFrontMatrix<TPZFileEqnStorage, TPZFrontNonSym> *mat = new TPZFrontMatrix<TPZFileEqnStorage, TPZFrontNonSym>(cmesh.NEquations());
574  //TPZFrontMatrix<TPZStackEqnStorage, TPZFrontNonSym> *mat = new TPZFrontMatrix<TPZStackEqnStorage, TPZFrontNonSym>(cmesh.NEquations());
575  //TPZFrontMatrix<TPZStackEqnStorage> *mat = new TPZFrontMatrix<TPZStackEqnStorage>(cmesh.NEquations());
576 
578 
579  // TPZFStructMatrix mat2(&cmesh);
580  // mat->SetNumElConnected(numelconnected);
581  //mat = CreateAssemble();
582 
583 
584  an.SetStructuralMatrix(mat);
585  // an2.SetStructuralMatrix(mat2);
586 
588  sol.SetDirect(ECholesky);
589  // TPZStepSolver sol2;
590  // sol2.SetDirect(ECholesky);
591  // sol.SetDirect(ELU);
592 
593 
594  an.SetSolver(sol);
595  // an2.SetSolver(sol2);
596  // mat->SetNumElConnected(numelconnected);
597  // mat->SetFileName("longhin.bin");
598  // an.Solver().SetDirect(ELU);
599  // mat->FinishWriting();
600  // mat->SetFileName('r',"longhin.bin");
601  // cout << "******************************************************************************************************AQUI 1" << endl;
602  an.Run(output);
603  an.Print("solution of frontal solver", output);
604  // cout << "******************************************************************************************************AQUI 2" << endl;
605  // an2.Run(output2);
606  // an2.Print("solution of frontal solver", output2);
607  /*
608  TPZVec<char *> scalnames(1);
609  scalnames[0] = "state";
610 
611  TPZVec<char *> vecnames(0);
612 
613  TPZDXGraphMesh graph(&cmesh,2,meumat,vecnames,scalnames);
614  ofstream *dxout = new ofstream("poisson.dx");
615  graph.SetOutFile(*dxout);
616  graph.SetResolution(0);
617 
618  //an.DefineGraphMesh(2, scalnames, vecnames, plotfile);
619  //an.Print("FEM SOLUTION ",output);
620  //an.PostProcess(1);
621  int istep = 0,numstep=1;
622 
623  graph.DrawMesh(numstep+1);
624  graph.DrawSolution(0,0);
625 
626  TPZAnalysis an2(&cmesh,output);
627  TPZFMatrix<STATE> *full = new TPZFMatrix(cmesh.NEquations(),cmesh.NEquations(),0.);
628  an2.SetMatrix(full);
629  an2.Solver().SetDirect(ELU);
630  an2.Run(output);
631  an2.Print("solution of full matrix", output);
632 
633  // full->Print("full decomposed matrix");
634  */
635  output.flush();
636  cout.flush();
637  return 0;
638 
639 }
640 #endif
641 
644 template<class front>
646 {
647  int64_t nconnect = this->fMesh->ConnectVec().NElements();
648  TPZManVector<int64_t> permute(nconnect);
650  int64_t i;
651  for(i=0; i<nconnect; i++)
652  {
653  permute[i] = i;
654  }
655  TPZCompEl *cel;
656  int64_t nelem = fElementOrder.NElements();
657  int64_t el;
658  int64_t connectcount = 0;
659  for(i=0; i<nelem; i++)
660  {
661  el = fElementOrder[i];
662  if(el<0) continue;
663  cel = fMesh->ElementVec()[el];
664  if(!cel) continue;
665  std::set<int64_t> indepconnects, depconnects;
666  cel->BuildConnectList(indepconnects,depconnects);
667  std::set<int64_t>::iterator it;
668  for(it=indepconnects.begin(); it != indepconnects.end(); it++)
669  {
670  TPZConnect &nod = fMesh->ConnectVec()[*it];
671  int nelcon = nod.NElConnected()-1;
672  int64_t seqnum = nod.SequenceNumber();
673  if(nelcon == 0) permute[seqnum]= connectcount++;
674  nod.DecrementElConnected();
675  }
676  }
677  for(i=0; i<nconnect; i++)
678  {
679  TPZConnect &nod = fMesh->ConnectVec()[i];
680  if(nod.SequenceNumber() < 0 || nod.NElConnected() <= 0) continue;
681  if(permute[nod.SequenceNumber()] < connectcount)
682  {
683  std::cout << __PRETTY_FUNCTION__ << " very fishy\n";
684  DebugStop();
685  }
686  }
687 #ifdef LOG4CXX
688  if (logger->isDebugEnabled())
689  {
690  std::stringstream sout;
691  sout << __PRETTY_FUNCTION__ << " permutation " << permute;
692  LOGPZ_DEBUG(logger,sout.str())
693  }
694 #endif
695  fMesh->Permute(permute);
696 }
697 
698 template<class TVar>
699 class TPZFrontSym;
700 template<class TVar>
701 class TPZFrontNonSym;
702 
705 
void GetNumElConnected(TPZVec< int > &numelconnected)
Returns a vector containing all elements connected to a degree of freedom.
int64_t NElements() const
Number of computational elements allocated.
Definition: pzcmesh.h:169
void OrderElement()
It is applied over fElementOrder putting it in the correct order.
Contains TPZAnalysis class which implements the sequence of actions to perform a finite element analy...
const int64_t numel
Number of elements to test.
Definition: pzsubcmesh.cpp:47
int Position(const int block_diagonal) const
Returns the position of first element block dependent on matrix diagonal.
Definition: pzblock.h:177
Represents a set of shape functions associated with a computational element/side. Computational Eleme...
Definition: pzconnect.h:30
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.
void SetStructuralMatrix(TPZAutoPointer< TPZStructMatrix > strmatrix)
Set structural matrix as auto pointer for analysis.
Definition: pzanalysis.cpp:82
void Print(const std::string &name, std::ostream &out)
Print connect and solution information.
Definition: pzanalysis.cpp:447
Responsible for the frontal method as a whole. Frontal.
void Assemble(TPZMatrix< STATE > &stiffness, TPZFMatrix< STATE > &rhs, TPZAutoPointer< TPZGuiInterface > guiInterface)
Assemble a stiffness matrix.
void DecrementElConnected()
Decrement fNElConnected.
Definition: pzconnect.h:262
virtual TPZGeoEl * CreateGeoElement(MElementType type, TPZVec< int64_t > &cornerindexes, int matid, int64_t &index, int reftype=1)
Generic method for creating a geometric element. Putting this method centrally facilitates the modifi...
Definition: pzgmesh.cpp:1296
virtual void BuildConnectList(std::set< int64_t > &indepconnectlist, std::set< int64_t > &depconnectlist)
Builds the list of all connectivities related to the element including the connects pointed to by dep...
Definition: pzcompel.cpp:435
Contains the TPZFrontStructMatrix class which responsible for a interface among Finite Element Packag...
TPZFront< TVar > & GetFront() override
returns a pointer to the front matrix
int64_t NEquations()
This computes the number of equations associated with non-restrained nodes.
Definition: pzcmesh.cpp:721
void AddFel(TPZFMatrix< TVar > &rhs, TPZVec< int64_t > &destination)
Performs a right hand side assemblage.
Definition: pzfmatrix.cpp:209
TPZFNMatrix< 1000, STATE > fMat
Pointer to a blocked matrix object.
Definition: pzelmat.h:41
Defines step solvers class. Solver.
Definition: pzmganalysis.h:17
void SetSolver(TPZMatrixSolver< STATE > &solver)
Set solver matrix.
virtual void CalcStiff(TPZElementMatrix &ek, TPZElementMatrix &ef)
Computes the element stifness matrix and right hand side.
Definition: pzcompel.h:794
TPZEquationFilter fEquationFilter
Object which will determine which equations will be assembled.
Abstract class implements storage and decomposition process of the frontal matrix involving simmetry ...
Declarates the TPZBlock<REAL>class which implements block matrices.
int NumberRigidBodyModes()
Return the number of rigid body modes associated with the internal degrees of freedom.
int64_t NElements() const
Access method to query the number of elements of the vector.
Implements a chunk vector with free store administration. Utility.
Definition: TPZStream.h:39
Contains the TPZStructMatrixOR class which responsible for a interface among Matrix and Finite Elemen...
TPZManVector< int64_t > fDestinationIndex
Definition: pzelmat.h:53
void ComputeDestinationIndices()
Definition: pzelmat.cpp:126
virtual ~TPZFrontStructMatrix()
Class destructor.
Contains the TPZFStructMatrix class which implements Full Structural Matrices.
void SetDefaultOrder(int order)
Definition: pzcmesh.h:403
Refines geometrical mesh (all the elements) num times.
Definition: pzstrmatrix.h:35
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
Contains declaration of TPZElementMatrix struct which associates an element matrix with the coeficien...
TPZCompMesh * fMesh
Pointer to the computational mesh from which the matrix will be generated.
DecomposeType fDecomposeType
Used Decomposition method.
void AdjustSequenceNumbering()
Resequence the connects according to the element order.
Contains the TPZBndCond class which implements a boundary condition for TPZMaterial objects...
virtual TPZMaterial * Material() const
Identify the material object associated with the element.
Definition: pzcompel.cpp:959
TPZVec< int > fElementOrder
This vector contains an ordered list.
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
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
Contains declaration of TPZMesh class which defines a geometrical mesh and contains a corresponding l...
int64_t SequenceNumber() const
Returns the Sequence number of the connect object.
Definition: pzconnect.h:158
void AssembleElement(TPZCompEl *el, TPZElementMatrix &ek, TPZElementMatrix &ef, TPZMatrix< STATE > &stiffness, TPZFMatrix< STATE > &rhs)
Computes element matrices.
Implements the sequence of actions to perform a finite element analysis. Analysis.
Definition: pzanalysis.h:32
int64_t NActiveEquations() const
Retorna o numero de equacoes ativas do sistema.
Implements a group of computational elements as a mesh and an element. Computational Mesh...
Definition: pzsubcmesh.h:36
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
void FinishWriting()
Finishes writing of a binary file and closes it.
#define LOGPZ_DEBUG(A, B)
Define log for debug info.
Definition: pzlog.h:87
void SetNumElConnected(TPZVec< int > &numelconnected)
Initializes the number of elements connected to each equation.
Contains declaration of TPZConnect class which represents a set of shape functions associated with a ...
virtual void ComputeNodElCon()
Compute the number of elements connected to each connect object.
Definition: pzcmesh.cpp:668
void AssembleNew(TPZMatrix< STATE > &stiffness, TPZFMatrix< STATE > &rhs, TPZAutoPointer< TPZGuiInterface > guiInterface)
Assemble a stiffness matrix according to rhs.
TPZAdmChunkVector< TPZGeoNode > & NodeVec()
Definition: pzgmesh.h:140
virtual void AutoBuild(const std::set< int > *MaterialIDs)
Creates the computational elements, and the degree of freedom nodes.
Definition: pzcmesh.cpp:308
int HasDependency() const
Returns whether exist dependecy information.
Definition: pzconnect.h:292
virtual void AddKel(TPZFMatrix< TVar > &elmat, TPZVec< int64_t > &destinationindex)
Add a contribution of a stiffness matrix.
Definition: pzmatrix.cpp:506
const TPZBlock< STATE > & Block() const
Access the block structure of the solution vector.
Definition: pzcmesh.h:213
Responsible for a interface among Finite Element Package and Matrices package to frontal method...
void SetNumRigidBodyModes(int nrigid)
Indicate the first equation dedicated to rigid body modes.
Definition: TPZFront.h:63
Contains the TPZMat2dLin class which implements a bi-dimensional linear problem.
TPZAdmChunkVector< TPZConnect > & ConnectVec()
Return a reference to the connect pointers vector.
Definition: pzcmesh.h:198
int64_t Index() const
Returns element index of the mesh fELementVec list.
Definition: pzcompel.h:821
void ReOpen()
Reopen the binary file.
virtual bool ShouldCompute(int matid) const
Establish whether the element should be computed.
TPZManVector< int64_t > fSourceIndex
Definition: pzelmat.h:53
Contains declaration of TPZCompMesh class which is a repository for computational elements...
TPZStructMatrix * Clone()
Clones a TPZFrontStructMatrix.
void SetMaterial(TPZFMatrix< STATE > &xkin, TPZFMatrix< STATE > &xcin, TPZFMatrix< STATE > &xfin)
Definition: pzmat2dlin.h:44
TPZMatrix< STATE > * Create()
Returns a pointer to TPZMatrix<STATE>
virtual void Run(std::ostream &out=std::cout)
Calls the appropriate sequence of methods to build a solution or a time stepping sequence.
Definition: pzanalysis.cpp:920
This class associates an element matrix with the coeficients of its contribution in the global stiffn...
Definition: pzelmat.h:30
TPZFNMatrix< 1000, STATE > fConstrMat
Pointer to the constrained matrix object.
Definition: pzelmat.h:48
virtual void SetDecomposeType(DecomposeType dectype)=0
Set the decomposition type.
void ApplyConstraints()
Apply the constraints applied to the nodes by transforming the tangent matrix and right hand side...
Definition: pzelmat.cpp:186
void BuildConnectivity()
Build the connectivity of the grid.
Definition: pzgmesh.cpp:967
virtual TPZBndCond * CreateBC(TPZMaterial *reference, int id, int typ, TPZFMatrix< STATE > &val1, TPZFMatrix< STATE > &val2)
Creates an object TPZBndCond derived of TPZMaterial.
Contains declaration of TPZSubCompMesh class which implements a group of computational elements as a ...
int InsertMaterialObject(TPZMaterial *mat)
Insert a material object in the datastructure.
Definition: pzcmesh.cpp:287
This class implements a geometric mesh for the pz environment. Geometry.
Definition: pzgmesh.h:48
Implements a bi-dimensional linear problem.
Definition: pzmat2dlin.h:22
Contains TPZSubMeshFrontalAnalysis class which implements the analysis for substructuring.
Implements computational mesh. Computational Mesh.
Definition: pzcmesh.h:47
int Size(const int block_diagonal) const
Returns block dimension.
Definition: pzblock.h:171
Contains the TPZAbstractFrontMatrix class which implements a matrix stored in a frontal decomposition...
TPZAdmChunkVector< TPZCompEl * > & ElementVec()
Returns a reference to the element pointers vector.
Definition: pzcmesh.h:183
Contains TPZSolver class which defines a abstract class of solvers which will be used by matrix class...
int NElConnected() const
Returns fNElConnected.
Definition: pzconnect.h:264
int Id() const
Definition: TPZMaterial.h:170
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
Contains TPZStepSolver class which defines step solvers class.
virtual void FilterEquations(TPZVec< int64_t > &origindex, TPZVec< int64_t > &destindex) const
Filter out the equations which are out of the range.
int64_t NElements() const
Returns the number of elements of the vector.
Definition: pzvec.h:190
void SetDirect(const DecomposeType decomp)
Contains the TPZDXGraphMesh class which implements the interface of the graphmesh to the OpenDX graph...
void InitializeBlock()
Resequence the block object, remove unconnected connect objects and reset the dimension of the soluti...
Definition: pzcmesh.cpp:391
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
virtual int HasDependency()
Returns 1 if the element has at least one dependent node. Returns 0 otherwise.
Definition: pzcompel.cpp:565
Abstract class implements storage and decomposition process of the frontal matrix involving non-simet...
virtual void Print(std::ostream &out=std::cout) const
Prints mesh data.
Definition: pzcmesh.cpp:236
void Permute(TPZVec< int64_t > &permute)
Permute the sequence number of the connect objects It is a permute gather operation.
Definition: pzcmesh.cpp:1321