NeoPZ
tpzgensubstruct.cpp
Go to the documentation of this file.
1 
6 #include "tpzgensubstruct.h"
7 #include "pzgengrid.h"
8 #include "pzgmesh.h"
9 #include "pzcmesh.h"
10 #include "pzbndcond.h"
11 #include "TPZGeoElement.h"
12 #include "pzgeopoint.h"
13 #include "pzrefpoint.h"
14 #include "pzsubcmesh.h"
15 
16 #include "tpznodesetcompute.h"
17 #include "pzmetis.h"
18 
19 #include "tpzdohrmatrix.h"
20 #include "tpzdohrsubstruct.h"
22 #include "tpzverysparsematrix.h"
23 
24 #include "pzanalysis.h"
25 
26 #include "pzskylstrmatrix.h"
27 
28 #include "pzsubcmesh.h"
29 
30 #include "tpzpairstructmatrix.h"
31 #include "tpzmatredstructmatrix.h"
32 
33 #include "pzlog.h"
34 
35 #include "TPZfTime.h"
36 #include "TPZTimeTemp.h"
37 
38 #include <sstream>
39 
40 #ifdef LOG4CXX
41 static LoggerPtr logger(Logger::getLogger("substruct.gensubstruct"));
42 #endif
43 
44 TPZGenSubStruct::TPZGenSubStruct(int dimension, int numlevels, int substructlevel) : fMatDist(DistMaterial),fK(8,1.),
45 fDimension(dimension),
46 fNumLevels(numlevels),fSubstructLevel(substructlevel)
47 {
48 }
49 
51 {
52 }
53 
54 #ifndef STATE_COMPLEX
55 REAL co[8][3] = {
57  {0.,0.,0.},
58  {1.,0.,0.},
59  {1.,1.,0.},
60  {0.,1.,0.},
61  {0.,0.,1.},
62  {1.,0.,1.},
63  {1.,1.,1.},
64  {0.,1.,1.}
65 };
66 
67 #include "pzpoisson3d.h"
68 // method which will generate the computational mesh
70 {
71  std::cout << "Generating mesh\n";
72  TPZGeoMesh *gmesh = new TPZGeoMesh;
73  if(fDimension == 2)
74  {
75  TPZVec<int> nx(2,1);
76  TPZVec<REAL> x0(3,0.),x1(3,1.);
77  x1[2] = 0.;
78  TPZGenGrid gen(nx,x0,x1,1,0.);
79  gen.Read(gmesh);
80  } else if(fDimension == 3)
81  {
82  int ic;
83  gmesh->NodeVec().Resize(8);
84  TPZVec<REAL> noco(3);
85  TPZVec<int64_t> nodeindices(8);
86  for(ic=0; ic<8; ic++)
87  {
88  nodeindices[ic] = ic;
89  int i;
90  for(i=0; i<3; i++)
91  {
92  noco[i] = co[ic][i];
93  }
94  gmesh->NodeVec()[ic].Initialize(noco,*gmesh);
95  }
96  int matid = 1;
97  int64_t index;
98  gmesh->CreateGeoElement(ECube,nodeindices,matid,index,0);
99  }
100  this->fCMesh = new TPZCompMesh(gmesh);
101  TPZVec<int64_t> nodeindices(1,0);
103  TPZVec<REAL> convdir(fDimension,1.);
104  TPZMatPoisson3d *matp;
105  int imat;
106  for(imat=0; imat<fK.NElements(); imat++)
107  {
108  matp = new TPZMatPoisson3d(imat+1,fDimension);
109  matp->SetInternalFlux(1.);
110  matp->SetParameters(fK[imat],0.,convdir);
111  {
112  TPZMaterial * mat (matp);
114  }
115  }
116 
117  TPZMaterial * mat (fCMesh->FindMaterial(1));
118  TPZFMatrix<STATE> val1(1,1,0.),val2(1,1,0.);
119  TPZBndCond *bc = new TPZBndCond(mat,-1,0,val1,val2);
120  TPZMaterial * matbc(bc);
122  gmesh->BuildConnectivity();
123  std::cout << "Uniform refine "; std::cout.flush();
124  UniformRefine();
125  std::cout << "AutoBuild "; std::cout.flush();
126  fCMesh->AutoBuild();
127  std::cout << "Number of equations " << fCMesh->NEquations() << std::endl;
128  //tempo.fNumEq = fCMesh->NEquations(); // alimenta timeTemp com o numero de equacoes
129 #ifdef LOG4CXX
130  {
131  std::stringstream str;
132  fCMesh->Print(str);
133  LOGPZ_DEBUG(logger,str.str());
134  }
135 #endif
136  // std::cout << "Identifying corner nodes\n";
137  // IdentifyCornerNodes();
138  return fCMesh;
139 }
140 #endif
141 // divide the geometric elements till num levels is achieved
143 {
144  TPZGeoMesh *gmesh = fCMesh->Reference();
145  int il;
146  for(il=0; il<fNumLevels; il++)
147  {
148  std::cout << il << ' '; std::cout.flush();
149  int nel = gmesh->NElements();
150  int iel;
151  int nk = fK.NElements();
152  for(iel=0; iel<nel; iel++)
153  {
154  TPZGeoEl *gel = gmesh->ElementVec()[iel];
155  if(gel->Level() < fNumLevels && !gel->HasSubElement() && gel->Dimension() > 0)
156  {
157  TPZVec<TPZGeoEl *> subel(4);
158  gel->Divide(subel);
159  if(fMatDist == DistMaterial)
160  {
161  if(il == 0 && gel->MaterialId() == 1)
162  {
163  int nsub = subel.NElements();
164  int is;
165  for(is=0; is<nsub; is++)
166  {
167  subel[is]->SetMaterialId(1+is%nk);
168  }
169  }
170  } else
171  {
172  if(gel->MaterialId() > 0)
173  {
174  int nsub = subel.NElements();
175  int is;
176  for(is=0; is<nsub; is++)
177  {
178  subel[is]->SetMaterialId(1+(rand()%nk));
179  }
180  }
181  }
182  }
183  }
184  }
185 }
186 
187 // divide the elements in substructures
189 {
190  TPZfTime timesubstructuring; // init of timer for substructuring mesh
191 
192  TPZGeoMesh *gmesh = fCMesh->Reference();
193  int nel = gmesh->NElements();
194  int iel;
195  for(iel=0; iel<nel; iel++)
196  {
197  TPZGeoEl *gel = gmesh->ElementVec()[iel];
198  if(gel->Level() == this->fSubstructLevel)
199  {
201  TPZGeoElSide gelside(gel,gel->NSides()-1);
202  gmesh->ResetReference();
204  gelside.HigherLevelCompElementList2(subels,0,0);
205  int nelstack = subels.NElements();
206  if(!nelstack)
207  {
208  TPZCompElSide celside = gelside.Reference();
209  if(celside.Element())
210  {
211  subels.Push(celside);
212  nelstack = 1;
213  }
214  }
215  int64_t index;
216  TPZCompMesh *cmesh = fCMesh.operator->();
217  TPZSubCompMesh *submesh = new TPZSubCompMesh(*cmesh,index);
218  std::cout << '*';
219  std::cout.flush();
220  int sub;
221  for(sub=0; sub<nelstack; sub++)
222  {
223  if(subels[sub].Reference().Dimension() == fDimension)
224  {
225  submesh->TransferElement(cmesh,subels[sub].Element()->Index());
226 #ifdef PZDEBUG
228 #endif
229  }
230  }
231  submesh->ExpandSolution();
232  }
233  }
234 
235  std::cout << timesubstructuring.ReturnTimeString(); // end of timer for substructuring mesh
236 
237  // transfer the point load
238  // find the point elements
239  nel = fCMesh->NElements();
240  TPZCompEl *celpoint = 0;
241  for(iel=0; iel<nel; iel++)
242  {
243  TPZCompEl *cel = fCMesh->ElementVec()[iel];
244  if(!cel) continue;
245  if(cel->Dimension() == 0)
246  {
247  celpoint = cel;
248  break;
249  }
250  }
251  if(celpoint)
252  {
253  int conindex = celpoint->ConnectIndex(0);
254  for(iel=0; iel<nel; iel++)
255  {
256  TPZCompEl *cel = fCMesh->ElementVec()[iel];
257  if(!cel) continue;
258  TPZSubCompMesh *submesh = dynamic_cast<TPZSubCompMesh *>(cel);
259  if(!submesh) continue;
260  int nc = submesh->NConnects();
261  int ic;
262  for(ic=0; ic<nc; ic++)
263  {
264  if(submesh->ConnectIndex(ic) == conindex)
265  {
266  submesh->TransferElement(fCMesh.operator->(),celpoint->Index());
267  celpoint = 0;
268  conindex = -1;
269  break;
270  }
271  }
272  if(!celpoint) break;
273  }
274  }
275 
276  //#define MAKEINTERNAL
277 #ifdef MAKEINTERNAL
278  std::cout << "Making Internal";
279  // make all nodes internal
280  nel = fCMesh->NElements();
281 
282  std::cout << "Make all Internal \n";
283  TPZfTime timeformakeallinternal; // init for timer
285  for(iel=0; iel<nel; iel++)
286  {
287  TPZCompEl *cel = fCMesh->ElementVec()[iel];
288  if(!cel) continue;
289  TPZSubCompMesh *submesh = dynamic_cast<TPZSubCompMesh *>(cel);
290  if(!submesh) continue;
291  std::cout << '-'; std::cout.flush();
292  submesh->MakeAllInternal();
293  }
294  std::cout << " == Finished\n";
296 
297  std::cout << timeformakeallinternal.ReturnTimeString();
298 #endif
299 }
300 
301 // identify cornernodes
303 {
304 #define COMPLETE
305 #ifdef COMPLETE
306  TPZNodesetCompute nodeset;
307  TPZStack<int64_t> elementgraph,elementgraphindex;
308  // fCompMesh->ComputeElGraph(elementgraph,elementgraphindex);
309  int nindep = fCMesh->NIndependentConnects();
310  // int neq = fCMesh->NEquations();
311  fCMesh->ComputeElGraph(elementgraph,elementgraphindex);
312  int nel = elementgraphindex.NElements()-1;
313  TPZMetis renum(nel,nindep);
314  //nodeset.Print(file,elementgraphindex,elementgraph);
315  std::cout << "Convert Graph ";
316  TPZfTime convertgraph;
317  renum.ConvertGraph(elementgraph,elementgraphindex,nodeset.Nodegraph(),nodeset.Nodegraphindex());
318  std::cout << convertgraph.ReturnTimeString();
319  // cout << "nodegraphindex " << nodeset.Nodegraphindex() << endl;
320  // cout << "nodegraph " << nodeset.Nodegraph() << endl;
321  std::cout << "AnalyseGraph ";
322  TPZfTime analysegraph;
323  nodeset.AnalyseGraph();
324  std::cout << analysegraph.ReturnTimeString();
325 #ifdef LOG4CXX
326  {
327  std::stringstream str;
328  nodeset.Print(str);
329  LOGPZ_DEBUG(logger,str.str());
330  }
331 #endif
332  int nnodes = nodeset.Levels().NElements();
333  int maxlev = nodeset.MaxLevel();
334  int in;
335  for(in=0; in<nnodes; in++)
336  {
337  if(nodeset.Levels()[in] == maxlev)
338  {
339  // this->fCornerEqs.insert(in);
340  // int seqnum = fCMesh->ConnectVec()[in].SequenceNumber();
341  int pos = fCMesh->Block().Position(in);
342  int size = fCMesh->Block().Size(in);
343  int ieq;
344  for(ieq=0; ieq<size; ieq++)
345  {
346  this->fCornerEqs.insert(pos+ieq);
347  }
348  }
349  }
350 
351 #else
353  int nindep = fCMesh->NIndependentConnects();
354  int in;
355  for(in=0; in<nindep; in++)
356  {
357  if(fCMesh->ConnectVec()[in].NElConnected() > 2 && fDimension == 2 ||
358  fCMesh->ConnectVec()[in].NElConnected() > 4 && fDimension == 3 )
359  {
360  int seqnum = fCMesh->ConnectVec()[in].SequenceNumber();
361  int pos = fCMesh->Block().Position(seqnum);
362  int size = fCMesh->Block().Size(seqnum);
363  int ieq;
364  for(ieq=0; ieq<size; ieq++)
365  {
366  this->fCornerEqs.insert(pos+ieq);
367  }
368  }
369  }
370 
371 #endif
372 #ifdef LOG4CXX
373  {
374  std::stringstream str;
375  str << "number of corner indices " << fCornerEqs.size() << std::endl;
376  str << " corner connect indices ";
377  std::set<int>::iterator it;
378  for(it=fCornerEqs.begin(); it!=fCornerEqs.end(); it++)
379  {
380  str << *it << " ";
381  }
382  LOGPZ_DEBUG(logger,str.str());
383  }
384 #endif
385 
386 
387 
388 }
389 
390 // initialize the TPZDohrMatrix structure
392 {
393  //Isolate each subcompmesh and put it in the dohrmann matrix
394  TPZDohrMatrix<STATE,TPZDohrSubstruct<STATE> > *dohr = dynamic_cast<TPZDohrMatrix<STATE,TPZDohrSubstruct<STATE> > *>(dohrmatrix.operator->());
396 
397  int neq = fCMesh->NEquations();
398  dohr->Resize(neq,neq);
399  // fCornerEqs was initialized during the mesh generation process
400  dohr->SetNumCornerEqs(this->fCornerEqs.size());
401 
402  int nsub = NSubMesh(fCMesh);
403  assembly->fFineEqs.Resize(nsub);
404  assembly->fCoarseEqs.Resize(nsub);
405  int isub;
406  std::cout << "Computing the system of equations for each substructure\n";
407  for(isub=0; isub<nsub; isub++)
408  {
409  TPZSubCompMesh *submesh = SubMesh(fCMesh, isub);
410  if(!submesh)
411  {
412  continue;
413  }
414  std::cout << '*';
415  // creating the substructure HERE
416  std::cout.flush();
418  // for each subcompmesh, reorder the nodes
419  //TPZAnalysis an(submesh);
420 
421  //keep the original sequence numbers of the connects
422 
423  // int nc = submesh->ConnectVec().NElements();
424  // TPZManVector<int> origseqnum(nc);
425  // int ic;
426  // for(ic=0; ic<nc; ic++) origseqnum[ic] = submesh->ConnectVec()[ic].SequenceNumber();
427 
428  // compute the stiffness matrix
429  int neq = ((TPZCompMesh *)submesh)->NEquations();
430  // int neq = substruct->fStiffness->Rows();
431 
432  substruct->fNEquations = neq;
433 
434  // identify the equation numbers of the submesh
435  std::map<int,int> globinv;
436  // initialize the fGlobalIndex data structure
437  // fGlobalIndex will have -1 entries for internal equations
438  // we need to dismember this (again) in two vectors
439  IdentifyEqNumbers(submesh, substruct->fGlobalEqs,globinv);
440  int next = substruct->fGlobalEqs.NElements();
441  assembly->fFineEqs[isub].Resize(next);
442  int ieq;
443  for(ieq=0; ieq<next; ieq++)
444  {
445  assembly->fFineEqs[isub][ieq] = substruct->fGlobalEqs[ieq].second;
446  }
447 
448 
449  IdentifyEqNumbers(submesh, substruct->fGlobalIndex,globinv);
450 
451  // initialize the fC matrix
452  // associate each column of the fC matrix with a coarse index
453  IdentifySubCornerEqs(globinv,substruct->fCoarseNodes,substruct->fCoarseIndex);
454  substruct->fC.Redim(substruct->fCoarseNodes.NElements(),neq);
455  for(ieq = 0; ieq<substruct->fCoarseNodes.NElements(); ieq++)
456  {
457  substruct->fC(ieq,substruct->fCoarseNodes[ieq]) = 1.;
458  }
459  int ncoarse = substruct->fCoarseIndex.NElements();
460  assembly->fCoarseEqs[isub].Resize(ncoarse);
461  for(ieq=0; ieq<ncoarse; ieq++)
462  {
463  assembly->fCoarseEqs[isub][ieq] = substruct->fCoarseIndex[ieq];
464  }
465 
466  // reorder by internal nodes
467  // the fInternalEqs data structure will not be filled if the connects are made internal
468 
469  // this permutes the nodes of the submesh
470  // This is a lengthy process which should run on the remote processor
471  InitializeMatrices(submesh, substruct, assembly);
472 
473 #ifdef LOG4CXX
474  {
475  std::stringstream sout;
476  /* sout << "Submesh for element " << iel << std::endl;
477  submesh->Print(sout);*/
478  sout << "Substructure for submesh " << isub << std::endl;
479  substruct->Print(sout);
480  LOGPZ_DEBUG(logger,sout.str())
481  }
482 #endif
483  dohr->AddSubstruct(substruct);
484  }
485  std::cout << std::endl;
486 }
487 
488 // initialize the TPZDohrMatrix structure
490 {
491  //Isolate each subcompmesh and put it in the dohrmann matrix
494 
495  int neq = fCMesh->NEquations();
496  dohr->Resize(neq,neq);
497  // fCornerEqs was initialized during the mesh generation process
498  dohr->SetNumCornerEqs(this->fCornerEqs.size());
499 
500  int nsub = NSubMesh(fCMesh);
501  assembly->fFineEqs.Resize(nsub);
502  assembly->fCoarseEqs.Resize(nsub);
503  int isub;
504  std::cout << "Computing the system of equations for each substructure\n";
505  for(isub=0; isub<nsub; isub++)
506  {
507  TPZSubCompMesh *submesh = SubMesh(fCMesh, isub);
508  if(!submesh)
509  {
510  continue;
511  }
512  std::cout << '*';
513  // creating the substructure HERE
514  std::cout.flush();
516 
517  // compute the stiffness matrix
518  int neq = ((TPZCompMesh *)submesh)->NEquations();
519  // int neq = substruct->fStiffness->Rows();
520 
521  substruct->fNEquations = neq;
522 
523  // identify the equation numbers of the submesh
524  std::map<int,int> globinv;
525  // initialize the fGlobalIndex data structure
526  // fGlobalIndex will have -1 entries for internal equations
527  // we need to dismember this (again) in two vectors
528  TPZVec<std::pair<int,int> > globaleqs;
529  IdentifyEqNumbers(submesh, globaleqs ,globinv);
530  int next = globaleqs.NElements();
531  substruct->fNumExternalEquations = next;
532  assembly->fFineEqs[isub].Resize(next);
533  int ieq;
534  for(ieq=0; ieq<next; ieq++)
535  {
536  assembly->fFineEqs[isub][ieq] = globaleqs[ieq].second;
537  }
538 
539  // initialize the permutations from the mesh enumeration to the external enumeration
541  typedef std::pair<ENumbering,ENumbering> Numberingpair;
542  ENumbering tsub,text,tint;
546 
547  TPZVec<int> &toexternal = substruct->fPermutationsScatter[Numberingpair(tsub,text)];
548  TPZVec<int> &fromexternal = substruct->fPermutationsScatter[Numberingpair(text,tsub)];
549  toexternal.Resize(neq,-1);
550  fromexternal.Resize(neq,-1);
551  int nel = globaleqs.NElements();
552 
553  for(ieq=0; ieq<nel; ieq++)
554  {
555  toexternal[globaleqs[ieq].first] = ieq;
556  }
557  int count = nel++;
558  for(ieq=0; ieq<neq; ieq++)
559  {
560  if(toexternal[ieq] == -1) toexternal[ieq] = count++;
561  }
562  for(ieq=0; ieq<neq; ieq++)
563  {
564  fromexternal[toexternal[ieq]] = ieq;
565  }
566 
567  ComputeInternalEquationPermutation(submesh, substruct->fPermutationsScatter[Numberingpair(tsub,tint)], substruct->fPermutationsScatter[Numberingpair(tint,tsub)]);
568  // IdentifyEqNumbers(submesh, substruct->fGlobalIndex,globinv);
569 
570  // initialize the fC matrix
571  // associate each column of the fC matrix with a coarse index
572  IdentifySubCornerEqs(globinv,substruct->fCoarseNodes,assembly->fCoarseEqs[isub]);
573  // int ncoarse = substruct->fCoarseNodes.NElements();
574 
575  // reorder by internal nodes
576  // the fInternalEqs data structure will not be filled if the connects are made internal
577 
578  // this permutes the nodes of the submesh
579  // This is a lengthy process which should run on the remote processor
580  InitializeMatrices(submesh, substruct, assembly);
581 
582 #ifdef LOG4CXX
583  {
584  std::stringstream sout;
585  sout << "The coarse equations are " << assembly->fCoarseEqs[isub] << std::endl;
586  /* sout << "Submesh for element " << iel << std::endl;
587  submesh->Print(sout);*/
588  sout << "Substructure for submesh " << isub << std::endl;
589  substruct->Print(sout);
590  LOGPZ_DEBUG(logger,sout.str())
591  }
592 #endif
593  dohr->AddSubstruct(substruct);
594  }
595  std::cout << std::endl;
596 }
597 
598 // identify the global equations as a pair of local equation and global equation
599 void TPZGenSubStruct::IdentifyEqNumbers(TPZSubCompMesh *sub, TPZVec<std::pair<int,int> > &globaleq, std::map<int,int> &globinv)
600 {
601  int ncon = sub->ConnectVec().NElements();
602  // ncon is the number of connects of the subcompmesh
603  TPZCompMesh *subcomp = (TPZCompMesh *) sub;
604  globaleq.Resize(subcomp->NEquations(),std::pair<int,int>(-1,-1));
605  TPZCompMesh *super = fCMesh.operator->();
606  int count = 0;
607  int ic;
608 #ifdef LOG4CXX_STOP
609  std::stringstream sout;
610  sout << "total submesh connects/glob/loc ";
611 #endif
612  for(ic=0; ic<ncon; ic++)
613  {
614  int glob = sub->NodeIndex(ic,super);
615  // continue is the connect is internal
616  if(glob == -1) continue;
617  int locseq = sub->ConnectVec()[ic].SequenceNumber();
618  int globseq = super->ConnectVec()[glob].SequenceNumber();
619  int locpos = sub->Block().Position(locseq);
620  int globpos = super->Block().Position(globseq);
621  int locsize = sub->Block().Size(locseq);
622  // int globsize = super->Block().Size(globseq);
623  int ieq;
624  for(ieq =0; ieq<locsize; ieq++)
625  {
626 #ifdef LOG4CXX_STOP
627  sout << ic << "/" << globpos+ieq << "/" << locpos+ieq << " ";
628 #endif
629  globaleq[count] = std::pair<int,int>(locpos+ieq,globpos+ieq);
630  count++;
631  globinv[globpos+ieq] = locpos+ieq;
632  }
633  }
634  globaleq.Resize(count);
635 #ifdef LOG4CXX_STOP
636  LOGPZ_DEBUG(logger,sout.str())
637 #endif
638 }
639 
640 
641 // get the global equation numbers of a substructure (and their inverse)
642 void TPZGenSubStruct::IdentifyEqNumbers(TPZSubCompMesh *sub, TPZVec<int> &global, std::map<int,int> &globinv)
643 {
644  int ncon = sub->ConnectVec().NElements();
645  // ncon is the number of connects of the subcompmesh
646  TPZCompMesh *subcomp = (TPZCompMesh *) sub;
647  global.Resize(subcomp->NEquations(),-1);
648  TPZCompMesh *super = fCMesh.operator->();
649  int ic;
650 #ifdef LOG4CXX_STOP
651  std::stringstream sout;
652  sout << "total submesh connects/glob/loc ";
653 #endif
654  for(ic=0; ic<ncon; ic++)
655  {
656  int glob = sub->NodeIndex(ic,super);
657  // continue is the connect is internal
658  if(glob == -1) continue;
659  int locseq = sub->ConnectVec()[ic].SequenceNumber();
660  int globseq = super->ConnectVec()[glob].SequenceNumber();
661  int locpos = sub->Block().Position(locseq);
662  int globpos = super->Block().Position(globseq);
663  int locsize = sub->Block().Size(locseq);
664  // int globsize = super->Block().Size(globseq);
665  int ieq;
666  for(ieq =0; ieq<locsize; ieq++)
667  {
668 #ifdef LOG4CXX_STOP
669  sout << ic << "/" << globpos+ieq << "/" << locpos+ieq << " ";
670 #endif
671  global[locpos+ieq] = globpos+ieq;
672  globinv[globpos+ieq] = locpos+ieq;
673  }
674  }
675 #ifdef LOG4CXX_STOP
676  LOGPZ_DEBUG(logger,sout.str())
677 #endif
678 }
679 
680 
684 void TPZGenSubStruct::ReorderInternalNodes(TPZSubCompMesh *sub,std::map<int,int> &globaltolocal, TPZVec<int> &internalnodes)
685 {
686  TPZVec<int64_t> permute;
687  sub->PermuteInternalFirst(permute);
688  int ninternal = this->NInternalEq(sub);
689  internalnodes.Resize(ninternal);
690  // this datastructure will not be initialized if the connects are made internal
691  internalnodes.Fill(-1);
692  int ncon = sub->ConnectVec().NElements();
693  TPZCompMesh *super = fCMesh.operator->();
694 #ifdef LOG4CXX_STOP
695  std::stringstream sout;
696  sout << "internal submesh connects/glob/loc ";
697 #endif
698  int ic;
699  for(ic=0; ic<ncon; ic++)
700  {
701  int glob = sub->NodeIndex(ic,super);
702  if(glob == -1) continue;
703  int locseq = sub->ConnectVec()[ic].SequenceNumber();
704  int globseq = super->ConnectVec()[glob].SequenceNumber();
705  int locpos = sub->Block().Position(locseq);
706  int globpos = super->Block().Position(globseq);
707  int locsize = sub->Block().Size(locseq);
708  // int globsize = super->Block().Size(globseq);
709  int ieq;
710  for(ieq =0; ieq<locsize; ieq++)
711  {
712  if(locpos+ieq < ninternal)
713  {
714 #ifdef LOG4CXX_STOP
715  sout << ic << "/" << globpos+ieq << "/" << locpos+ieq << " ";
716 #endif
717  internalnodes[locpos+ieq] = globaltolocal[globpos+ieq];
718  }
719  }
720  }
721 #ifdef LOG4CXX_STOP
722  LOGPZ_DEBUG(logger,sout.str())
723 #endif
724 
725 }
726 
727 // computes the permutation vectors from the subcompmesh ordening to the "internal first" ordering
728 // the mesh is modified during this method but is returned to its original state at the end of execution
730  TPZVec<int> &scatterpermute, TPZVec<int> &gatherpermute)
731 {
732  // This permutation vector is with respect to the blocks of the mesh
733  TPZVec<int64_t> scatterpermuteblock;
734  sub->ComputePermutationInternalFirst(scatterpermuteblock);
735  TPZBlock<STATE> destblock = sub->Block();
736  TPZBlock<STATE> &origblock = sub->Block();
737  int nblocks = origblock.NBlocks();
738  if(scatterpermuteblock.NElements() != origblock.NBlocks())
739  {
740  std::cout << __PRETTY_FUNCTION__ << " something seriously wrong!!!\n";
741  }
742  int ib;
743  for(ib=0; ib<nblocks; ib++)
744  {
745  destblock.Set(scatterpermuteblock[ib],origblock.Size(ib));
746  }
747  destblock.Resequence();
748 
749  int neq = ((TPZCompMesh *)sub)->NEquations();
750  scatterpermute.Resize(neq);
751  gatherpermute.Resize(neq);
752  scatterpermute.Fill(-1);
753  gatherpermute.Fill(-1);
754  int ncon = sub->ConnectVec().NElements();
755 #ifdef LOG4CXX_STOP
756  std::stringstream sout;
757  sout << "internal submesh connects/glob/loc ";
758 #endif
759  int ic;
760  for(ic=0; ic<ncon; ic++)
761  {
762  // skip dependent connects
763  TPZConnect &con = sub->ConnectVec()[ic];
764  if(con.HasDependency() || con.IsCondensed()) continue;
765  int locseq = sub->ConnectVec()[ic].SequenceNumber();
766  // skip unused connects
767  if(locseq < 0) continue;
768  int destseq = scatterpermuteblock[locseq];
769  int locpos = origblock.Position(locseq);
770  int destpos = destblock.Position(destseq);
771  int size = origblock.Size(locseq);
772  // int globsize = super->Block().Size(globseq);
773  int ieq;
774  for(ieq =0; ieq<size; ieq++)
775  {
776 #ifdef LOG4CXX_STOP
777  sout << ic << "/" << locpos+ieq << "/" << destpos+ieq << " ";
778 #endif
779  scatterpermute[locpos+ieq] = destpos+ieq;
780  }
781  }
782  int ieq;
783  for(ieq = 0; ieq < neq; ieq++)
784  {
785  gatherpermute[scatterpermute[ieq]] = ieq;
786  }
787 #ifdef LOG4CXX_STOP
788  LOGPZ_DEBUG(logger,sout.str())
789 #endif
790 
791 }
792 
797 {
798  TPZBlock<STATE> prevblock = sub->Block();
799  // This permutation vector is with respect to the blocks of the mesh
800  TPZVec<int64_t> permute;
801  sub->PermuteInternalFirst(permute);
802  blockinvpermute.Resize(permute.NElements());
803  int i;
804  for(i=0; i<permute.NElements(); i++)
805  {
806  blockinvpermute[permute[i]] = i;
807  }
808  int ninternal = NInternalEq(sub);
809  internaleqs.Resize(ninternal);
810  // this datastructure will not be initialized if the connects are made internal
811  internaleqs.Fill(-1);
812  int ncon = sub->ConnectVec().NElements();
813 #ifdef LOG4CXX_STOP
814  std::stringstream sout;
815  sout << "internal submesh connects/glob/loc ";
816 #endif
817  int ic;
818  for(ic=0; ic<ncon; ic++)
819  {
820  int locseq = sub->ConnectVec()[ic].SequenceNumber();
821  int origseq = blockinvpermute[locseq];
822  int locpos = sub->Block().Position(locseq);
823  int origpos = prevblock.Position(origseq);
824  int size = sub->Block().Size(locseq);
825  // int globsize = super->Block().Size(globseq);
826  int ieq;
827  for(ieq =0; ieq<size; ieq++)
828  {
829  if(locpos+ieq < ninternal)
830  {
831 #ifdef LOG4CXX_STOP
832  sout << ic << "/" << globpos+ieq << "/" << locpos+ieq << " ";
833 #endif
834  internaleqs[locpos+ieq] = origpos+ieq;
835  }
836  }
837  }
838 #ifdef LOG4CXX_STOP
839  LOGPZ_DEBUG(logger,sout.str())
840 #endif
841 
842 }
843 
844 // Identify the corner equations associated with a substructure
845 void TPZGenSubStruct::IdentifySubCornerEqs(std::map<int,int> &globaltolocal, TPZVec<int> &cornereqs,
846  TPZVec<int> &coarseindex)
847 {
848 #ifdef LOG4CXX
849  {
850  std::stringstream sout;
851  sout << "Input data for IdentifySubCornerEqs \nglobaltolocal";
852  std::map<int,int>::iterator mapit;
853  for(mapit = globaltolocal.begin(); mapit != globaltolocal.end(); mapit++)
854  {
855  sout << " [" << mapit->first << " , " << mapit->second << "] ";
856  }
857  sout << "\nCorner equations stored in the GenSubStructure data ";
858  std::set<int>::iterator setit;
859  for(setit = fCornerEqs.begin(); setit != fCornerEqs.end(); setit++)
860  {
861  sout << *setit << " , ";
862  }
863  sout << "\ncornereqs " << cornereqs;
864  LOGPZ_DEBUG(logger,sout.str())
865  }
866 #endif
867 
868  // REESCREVER ESTA PARTE
869  std::set<int>::iterator it;
870  std::list<int> subcorn,coarseindexlist;
871  int count = 0;
872  // the corner equations are all corner equations of the supermesh
873  // globaltolocal is the map of one particular submesh
874  for(it = fCornerEqs.begin(); it!= fCornerEqs.end(); it++,count++)
875  {
876  if(globaltolocal.find(*it) != globaltolocal.end())
877  {
878  subcorn.push_back(globaltolocal[*it]);
879  coarseindexlist.push_back(count);
880  }
881  }
882  std::list<int>::iterator lit;
883  cornereqs.Resize(subcorn.size());
884  coarseindex.Resize(coarseindexlist.size());
885  count = 0;
886  for(lit=subcorn.begin(); lit != subcorn.end(); lit++)
887  {
888  cornereqs[count++] = *lit;
889  }
890  count = 0;
891  for(lit = coarseindexlist.begin(); lit != coarseindexlist.end(); lit++)
892  {
893  coarseindex[count++] = *lit;
894  }
895 }
896 
901 {
902  std::list<int64_t> internal;
903  sub->PotentialInternal(internal);
904  std::list<int64_t>::iterator it;
905  int64_t result = 0;
906  for(it=internal.begin(); it != internal.end(); it++)
907  {
908  int64_t seq = sub->ConnectVec()[*it].SequenceNumber();
909  int64_t sz = sub->Block().Size(seq);
910  result += sz;
911  }
912  return result;
913 }
914 
915 // This is a lengthy process which should run on the remote processor
917 {
918  // this should happen in the remote processor
919  TPZSkylineStructMatrix skylstr(submesh);
921 
922  skylstr.EquationFilter().Reset();
923  substruct->fStiffness = TPZAutoPointer<TPZMatrix<STATE> > (skylstr.CreateAssemble(substruct->fLocalLoad,toto));
924 
925  // This should happen in the remote processor
926  substruct->fInvertedStiffness.SetMatrix(substruct->fStiffness->Clone());
927  substruct->fInvertedStiffness.SetDirect(ECholesky);
928 
929  TPZManVector<int64_t> invpermute;
930  TPZGenSubStruct::ReorderInternalNodes2(submesh,substruct->fInternalEqs,invpermute);
931  // compute the stiffness matrix associated with the internal nodes
932  // fInternalEqs indicates the permutation of the global equations to the numbering of the internal equations
933  // this is meaningless if the internal nodes are invisible to the global structure
934  int64_t ninternal = substruct->fInternalEqs.NElements();
935  // THIS SHOULD HAPPEN IN THE REMOTE PROCESSOR
936  skylstr.SetEquationRange(0,ninternal);
937  TPZFMatrix<STATE> rhs;
938  substruct->fInvertedInternalStiffness.SetMatrix(skylstr.CreateAssemble(rhs,toto));
939  substruct->fInvertedInternalStiffness.SetDirect(ECholesky);
940 
941  // put back the original sequence numbers of the connects (otherwise we can't apply a load solution
942 
943  submesh->Permute(invpermute);
944 
945 }
946 
947 // This is a lengthy process which should run on the remote processor
949 {
951  typedef std::pair<ENumbering,ENumbering> pairnumbering;
953  TPZVec<int> &permutescatter = substruct->fPermutationsScatter[fromsub];
954 
955  // create a skyline matrix based on the current numbering of the mesh
956  // put the stiffness matrix in a TPZMatRed object to facilitate the computation of phi and zi
957  TPZSkylineStructMatrix skylstr(submesh);
958  skylstr.EquationFilter().Reset();
959 
960  TPZAutoPointer<TPZMatrix<STATE> > Stiffness = skylstr.Create();
961  TPZMatRed<STATE,TPZFMatrix<STATE> > *matredbig = new TPZMatRed<STATE, TPZFMatrix<STATE> >(Stiffness->Rows()+substruct->fCoarseNodes.NElements(),Stiffness->Rows());
962  matredbig->SetK00(Stiffness);
963  substruct->fMatRedComplete = matredbig;
964 
965 
966  TPZVec<int64_t> permuteconnectscatter;
967 
968  substruct->fNumInternalEquations = submesh->NumInternalEquations();
969  // change the sequencing of the connects of the mesh, putting the internal connects first
970  submesh->PermuteInternalFirst(permuteconnectscatter);
971 
972  // create a "substructure matrix" based on the submesh using a skyline matrix structure as the internal matrix
976 
977  // create a structural matrix which will assemble both stiffnesses simultaneously
978  TPZPairStructMatrix pairstructmatrix(submesh,permutescatter);
979 
980  // reorder the sequence numbering of the connects to reflect the original ordering
981  TPZVec<int64_t> invpermuteconnectscatter(permuteconnectscatter.NElements());
982  int64_t iel;
983  for (iel=0; iel < permuteconnectscatter.NElements(); iel++) {
984  invpermuteconnectscatter[permuteconnectscatter[iel]] = iel;
985  }
986  TPZAutoPointer<TPZMatrix<STATE> > InternalStiffness = matredptr->K00();
987 
988 
989  submesh->Permute(invpermuteconnectscatter);
990 
991  // compute both stiffness matrices simultaneously
992  substruct->fLocalLoad.Redim(Stiffness->Rows(),1);
993  pairstructmatrix.Assemble(Stiffness.operator->(), matredptr, substruct->fLocalLoad);
994  matredbig->SimetrizeMatRed();
995  matredptr->SimetrizeMatRed();
996 
997  substruct->fWeights.Resize(Stiffness->Rows());
998  int64_t i;
999  for(i=0; i<substruct->fWeights.NElements(); i++)
1000  {
1001  substruct->fWeights[i] = Stiffness->GetVal(i,i);
1002  }
1003  // Desingularize the matrix without affecting the solution
1004  int64_t ncoarse = substruct->fCoarseNodes.NElements(), ic;
1005  int64_t neq = Stiffness->Rows();
1006  for(ic=0; ic<ncoarse; ic++)
1007  {
1008  int64_t coarse = substruct->fCoarseNodes[ic];
1009  Stiffness->operator()(coarse,coarse) += 10.;
1010  matredbig->operator()(neq+ic,coarse) = 1.;
1011  matredbig->operator()(coarse,neq+ic) = 1.;
1012  }
1013  //substruct->fStiffness = Stiffness;
1014  TPZStepSolver<STATE> *InvertedStiffness = new TPZStepSolver<STATE>(Stiffness);
1015  InvertedStiffness->SetMatrix(Stiffness);
1016  InvertedStiffness->SetDirect(ECholesky);
1017  matredbig->SetSolver(InvertedStiffness);
1018 
1019 
1020  TPZStepSolver<STATE> *InvertedInternalStiffness = new TPZStepSolver<STATE>(InternalStiffness);
1021  InvertedInternalStiffness->SetMatrix(InternalStiffness);
1022  InvertedInternalStiffness->SetDirect(ECholesky);
1023  matredptr->SetSolver(InvertedInternalStiffness);
1024  matredptr->SetReduced();
1026 
1027  substruct->fMatRed = matred2;
1028 }
1029 
1030 // return the number of submeshes
1032 {
1033  int64_t nel = compmesh->NElements();
1034  TPZCompEl *cel;
1035  int64_t iel, count = 0;
1036  for(iel=0; iel<nel; iel++)
1037  {
1038  cel = compmesh->ElementVec()[iel];
1039  if(!cel) continue;
1040  TPZSubCompMesh *sub = dynamic_cast<TPZSubCompMesh *>(cel);
1041  if(sub) count++;
1042  }
1043  return count;
1044 }
1045 
1046 // return a pointer to the isub submesh
1048 {
1049  int64_t nel = compmesh->NElements();
1050  TPZCompEl *cel;
1051  int64_t iel, count = 0;
1052  for(iel=0; iel<nel; iel++)
1053  {
1054  cel = compmesh->ElementVec()[iel];
1055  if(!cel) continue;
1056  TPZSubCompMesh *sub = dynamic_cast<TPZSubCompMesh *>(cel);
1057  if(sub && isub == count) return sub;
1058  if(sub) count++;
1059  }
1060  return NULL;
1061 }
int64_t NElements() const
Number of computational elements allocated.
Definition: pzcmesh.h:169
void InitializeDohr(TPZAutoPointer< TPZMatrix< STATE > > dohr, TPZAutoPointer< TPZDohrAssembly< STATE > > assembly)
Initialize the TPZDohrMatrix structure.
TPZGeoMesh * Reference() const
Returns a pointer to the geometrical mesh associated.
Definition: pzcmesh.h:209
void AnalyseGraph()
Group the node graph as passed by the parameters.
TPZAutoPointer< TPZCompMesh > fCMesh
computational mesh
Contains the TPZTimeTemp class which takes times.
Contains TPZAnalysis class which implements the sequence of actions to perform a finite element analy...
static void ReorderInternalNodes2(TPZSubCompMesh *sub, TPZVec< int > &internalnodes, TPZVec< int64_t > &invpermute)
int Position(const int block_diagonal) const
Returns the position of first element block dependent on matrix diagonal.
Definition: pzblock.h:177
std::set< int > fCornerEqs
The set of equations which correspond to corner nodes.
Implements computational element and a side. Computational Element.
Definition: pzcompel.h:632
virtual void SetEquationRange(int64_t mineq, int64_t maxeq)
int Set(const int index, const int dim, const int pos=-1)
Modifies existing block dimensions or creates a new block with given index.
Definition: pzblock.cpp:104
Represents a set of shape functions associated with a computational element/side. Computational Eleme...
Definition: pzconnect.h:30
void Print(std::ostream &file) const
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 SetReduced()
changes the declared dimension of the matrix to fDim1
Definition: pzmatred.h:70
int MaterialId() const
Returns the material index of the element.
Definition: pzgeoel.h:250
clarg::argBool bc("-bc", "binary checkpoints", false)
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
int64_t NodeIndex(int64_t nolocal, TPZCompMesh *super)
Gives the id node of one local node in containing mesh.
Definition: pzsubcmesh.cpp:358
Implements renumbering for elements of a mesh. Utility.
Definition: pzmetis.h:17
void ConvertGraph(TPZVec< int64_t > &elgraph, TPZVec< int64_t > &elgraphindex, TPZManVector< int64_t > &nodegraph, TPZManVector< int64_t > &nodegraphindex)
Will convert an element graph defined by elgraph and elgraphindex into a node graph defined by nodegr...
virtual TPZMatrix< STATE > * Create()
Contains the TPZDohrSubstruct class which implements sub structure matrices using Dohrman algorithm...
void SimetrizeMatRed()
If fK00 is simetric, only part of the matrix is accessible to external objects.
Definition: pzmatred.cpp:67
TPZManVector< int64_t > & Nodegraphindex()
clarg::argInt nsub("-nsub", "number of substructs", 4)
int64_t NEquations()
This computes the number of equations associated with non-restrained nodes.
Definition: pzcmesh.cpp:721
int64_t NIndependentConnects()
Number of independent connect objects.
Definition: pzcmesh.cpp:1080
void InitializeMatrices(TPZSubCompMesh *submesh, TPZAutoPointer< TPZDohrSubstruct< STATE > > substruct, TPZDohrAssembly< STATE > &dohrassembly)
This is a lengthy process which should run on the remote processor.
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
TPZGeoElSide Reference() const
Reference to the geometric element.
Definition: pzcompel.cpp:748
Implements a generic geometric element with a uniform refinement pattern. Geometry.
Definition: TPZGeoElement.h:22
Defines step solvers class. Solver.
Definition: pzmganalysis.h:17
Calculate the Times. Utility.
Definition: TPZfTime.h:21
virtual int Dimension() const =0
Dimension of the element.
Implements the generation of a multilayered bi-dimensional geometric grid. Getting Data...
Definition: pzgengrid.h:29
virtual TPZEquationFilter & EquationFilter()
access method for the equation filter
virtual int Resize(const int64_t newRows, const int64_t newCols)
Redimensions a matriz keeping the previous values.
Definition: pzmatrix.h:278
Contains the TPZDohrMatrix class which implements a matrix divided into substructures. Also contains the TPZDohrThreadMultData and TPZDohrThreadMultList structs.
int fSubstructLevel
Level of substructures.
int64_t NElements() const
Number of elements of the mesh.
Definition: pzgmesh.h:129
void InitializeDohrCondense(TPZAutoPointer< TPZMatrix< STATE > > dohr, TPZAutoPointer< TPZDohrAssembly< STATE > > assembly)
Initialize the TPZDohrMatrix structure.
void SetNumCornerEqs(int nc)
int64_t NSubMesh(TPZAutoPointer< TPZCompMesh > compmesh)
Return the number of submeshes.
Utility class which represents an element with its side. The Geometric approximation classes Geometry...
Definition: pzgeoelside.h:83
void AddSubstruct(TPZAutoPointer< TSubStruct > substruct)
It adds a substruct.
int64_t NElements() const
Access method to query the number of elements of the vector.
Implements sub structure matrices using Dohrman algorithm. Sub Structure.
virtual int NSides() const =0
Returns the number of connectivities of the element.
virtual void Divide(TPZVec< TPZGeoEl *> &pv)
Divides the element and puts the resulting elements in the vector.
Definition: pzgeoel.h:742
Implements a matrix divided into substructures. Matrix Sub structure.
Definition: tpzdohrmatrix.h:30
void LoadReferences()
Map this grid in the geometric grid.
Definition: pzcmesh.cpp:482
void SetK00(TPZAutoPointer< TPZMatrix< TVar > > K00)
Sets the matrix pointer of the upper left matrix to K00.
Definition: pzmatred.cpp:134
Contains the TPZMatPoisson3d class.
void ComputeElGraph(TPZStack< int64_t > &elgraph, TPZVec< int64_t > &elgraphindex)
Computes the connectivity graph of the elements, as appropriate for the TPZRenumbering class...
Definition: pzcmesh.cpp:1091
void PermuteInternalFirst(TPZVec< int64_t > &permute)
Permutes the potentially internal connects to the first on the list Respect the previous order of th...
void SetSolver(TPZAutoPointer< TPZMatrixSolver< TVar > > solver)
Definition: pzmatred.cpp:125
Implements SkyLine Structural Matrices. Structural Matrix.
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
Contains the TPZBndCond class which implements a boundary condition for TPZMaterial objects...
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 the TPZGenSubStruct class which is an interface to "feed" the datastructure of the Dohrmann ...
TPZSubCompMesh * SubMesh(TPZAutoPointer< TPZCompMesh > compmesh, int isub)
Return a pointer to the isub submesh.
void Resize(const int newsize)
Increase the size of the chunk vector.
Definition: pzadmchunk.h:280
Contains declaration of TPZMesh class which defines a geometrical mesh and contains a corresponding l...
void ReorderInternalNodes(TPZSubCompMesh *sub, std::map< int, int > &globaltolocal, TPZVec< int > &internalnodes)
static int64_t NInternalEq(TPZSubCompMesh *sub)
virtual TPZMatrix< STATE > * Create()
void Push(const T object)
Pushes a copy of the object on the stack.
Definition: pzstack.h:80
int Level()
Returns the number of ancestors.
Definition: pzgeoel.cpp:319
void SubStructure()
Divide the elements in substructures.
virtual void SetParameters(STATE diff, REAL conv, TPZVec< REAL > &convdir)
Definition: pzpoisson3d.cpp:78
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
void IdentifyEqNumbers(TPZSubCompMesh *sub, TPZVec< std::pair< int, int > > &globaleq, std::map< int, int > &globinv)
Identify the global equations as a pair of local equation and global equation.
To condense matrix divided in sub structures. Sub Structure.
Computes the cardinality of a nodegraph, identifying nodes as vertices, lines, faces or volumes...
Implements a group of computational elements as a mesh and an element. Computational Mesh...
Definition: pzsubcmesh.h:36
virtual void ExpandSolution()
Adapt the solution vector to new block dimensions.
Definition: pzcmesh.cpp:396
#define LOGPZ_DEBUG(A, B)
Define log for debug info.
Definition: pzlog.h:87
This class defines the boundary condition for TPZMaterial objects.
Definition: pzbndcond.h:29
Contains the TPZSkylineStructMatrix class which implements SkyLine Structural Matrices.
virtual void ComputeNodElCon()
Compute the number of elements connected to each connect object.
Definition: pzcmesh.cpp:668
int fDimension
Dimension of the mesh.
TPZVec< int > & Levels()
Returns the level of the nodes.
virtual void SetMatrix(TPZAutoPointer< TPZMatrix< TVar > > Refmat)
Sets a matrix to the current object.
Definition: pzsolve.h:115
TPZMaterial * FindMaterial(int id)
Find the material with identity id.
Definition: pzcmesh.cpp:297
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
Contains the TPZGenGrid class which implements the generation of a multilayered geometric grid (two-d...
Contains the TPZPairStructMatrix class.
const TPZBlock< STATE > & Block() const
Access the block structure of the solution vector.
Definition: pzcmesh.h:213
Contains the TPZDohrSubstructCondense class which condenses matrix divided in sub structures...
REAL co[8][3]
Coordinates of the eight nodes.
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
int fNumLevels
Number of uniform refinements.
void UniformRefine()
Divide the geometric elements till num levels is achieved.
Contains the TPZfTime class which calculates times.
Contains declaration of TPZCompMesh class which is a repository for computational elements...
void PotentialInternal(std::list< int64_t > &connectindices) const
Puts the nodes which can be transferred in an ordered list.
Definition: pzsubcmesh.cpp:789
void IdentifySubCornerEqs(std::map< int, int > &globaltolocal, TPZVec< int > &cornereqs, TPZVec< int > &coarseindex)
Identify the corner equations associated with a substructure.
ETypemesh fMatDist
Contains the TPZGeoPoint class which implements the geometry of a point element or 0-D element...
void IdentifyCornerNodes()
Identify cornernodes.
virtual int64_t ConnectIndex(int i) const =0
Returns the index of the ith connectivity of the element.
virtual int HasSubElement() const =0
Return 1 if the element has subelements.
.. . Sub Structure
int NBlocks() const
Returns number of blocks on diagonal.
Definition: pzblock.h:165
void Assemble(TPZMatrix< STATE > *first, TPZMatrix< STATE > *second, TPZFMatrix< STATE > &rhs)
std::string ReturnTimeString()
When called, returns the time since the creation of the object in a string.
Definition: TPZfTime.cpp:20
void ComputePermutationInternalFirst(TPZVec< int64_t > &permute) const
Computes the permutation vector which puts the internal connects to the first on the list Respect th...
void BuildConnectivity()
Build the connectivity of the grid.
Definition: pzgmesh.cpp:967
TPZAutoPointer< TPZCompMesh > GenerateMesh()
Method which will generate the computational mesh.
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
DESCRIBE PLEASE.
Definition: pzpoisson3d.h:26
virtual void CleanUpUnconnectedNodes()
Delete the nodes which have no elements connected to them.
Definition: pzcmesh.cpp:498
virtual int Dimension() const =0
Returns the dimension of the element.
TPZCompEl * Element() const
Gives a pointer to the reference computational element.
Definition: pzcompel.h:675
This class implements a very simple interface from PZ kernel to GUI. Module: Common.
Implements a simple substructuring of a linear system of equations, composed of 4 submatrices...
Definition: pzmatred.h:34
TPZAutoPointer< TPZMatrix< TVar > > K00()
Definition: pzmatred.h:112
Contains the TPZMatRedStructMatrix class.
This class implements a geometric mesh for the pz environment. Geometry.
Definition: pzgmesh.h:48
virtual int64_t ConnectIndex(int i) const override
Returns the connection index i.
Definition: pzsubcmesh.cpp:346
virtual void MakeAllInternal() override
Makes all mesh connections internal mesh connections.
Definition: pzsubcmesh.cpp:602
int MaxLevel()
Returns the maximum level.
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 int64_t TransferElement(TPZCompMesh *mesh, int64_t elindex) override
Transfer the element elindex belonging to mesh to the current mesh and returns its index...
Definition: pzsubcmesh.cpp:987
Definition: pzeltype.h:61
int64_t NumInternalEquations()
Computes the number of internal equations.
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
TPZGenSubStruct(int dimension, int numlevels, int substructlevel)
Constructor.
Contains TPZVerySparseMatrix class which implements a matrix whose nonzero elements are stored in bin...
TPZVec< REAL > fK
int64_t NElements() const
Returns the number of elements of the vector.
Definition: pzvec.h:190
void SetInternalFlux(STATE flux)
Definition: pzpoisson3d.h:179
void SetDirect(const DecomposeType decomp)
Contains declaration of TPZGeoElement class which implements a generic geometric element with a unifo...
bool VerifyDatastructureConsistency()
Method to verify that the datastructures are consistent.
void ResetReference()
Resets all load references in elements and nodes.
Definition: pzgmesh.cpp:197
Contains the TPZNodesetCompute class which computes the cardinality of a nodegraph.
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
virtual TPZMatrix< STATE > * CreateAssemble(TPZFMatrix< STATE > &rhs, TPZAutoPointer< TPZGuiInterface > guiInterface, unsigned numthreads_assemble, unsigned numthreads_decompose)
Definition: pzstrmatrix.h:100
bool IsCondensed() const
Access method to return the indication whether the connect is condensed or not.
Definition: pzconnect.h:223
virtual short Read(TPZGeoMesh *mesh, int matid=1)
Add nodes and elements to the object mesh.
Definition: pzgengrid.cpp:64
Contains the TPZRefPoint class which implements the uniform refinement of a geometric point element...
Contains TPZMetis class which implements the renumbering for elements of a mesh to minimize the band...
TPZManVector< int64_t > & Nodegraph()
void Reset()
Reset method.
static void ComputeInternalEquationPermutation(TPZSubCompMesh *sub, TPZVec< int > &scatterpermute, TPZVec< int > &gatherpermute)
Computes the permutation vectors from the subcompmesh ordening to the "internal first" ordering...
int Resequence(const int start=0)
Resequences blocks positioning.
Definition: pzblock.cpp:150
virtual void Print(std::ostream &out=std::cout) const
Prints mesh data.
Definition: pzcmesh.cpp:236
TPZAdmChunkVector< TPZGeoEl * > & ElementVec()
Methods for handling pzlists.
Definition: pzgmesh.h:138
void Permute(TPZVec< int64_t > &permute)
Permute the sequence number of the connect objects It is a permute gather operation.
Definition: pzcmesh.cpp:1321
virtual int NConnects() const override
Returns the number of connections.
Definition: pzsubcmesh.cpp:342