NeoPZ
pzsubcmesh.cpp
Go to the documentation of this file.
1 
6 #include "pzsubcmesh.h"
7 #include "pzgmesh.h"
8 #include "pzcompel.h"
9 #include "TPZInterfaceEl.h"
11 #include "pzcmesh.h"
12 #include "pzelmat.h"
13 #include "pznonlinanalysis.h"
14 #include "pzskylmat.h"
15 #include "pzsolve.h"
16 #include "pzstepsolver.h"
17 #include "pzskylstrmatrix.h"
19 #include "pzfstrmatrix.h"
20 #include "TPZFrontStructMatrix.h"
22 #include "pzsmfrontalanal.h"
23 #include "pzsmanal.h"
24 #include "pzbndcond.h"
25 #include "pzvisualmatrix.h"
26 
27 #include "pzcondensedcompel.h"
28 #include "pzelementgroup.h"
29 
30 
31 #ifndef STATE_COMPLEX
32 //#include "pzmathyperelastic.h"
33 #endif
34 
35 #include <stdio.h>
36 
37 #include <sstream>
38 #include <iterator>
39 #include "pzlog.h"
40 
41 #ifdef LOG4CXX
42 static LoggerPtr logger(Logger::getLogger("pz.mesh.subcmesh"));
43 static LoggerPtr logger2(Logger::getLogger("pz.mesh.tpzcompmesh"));
44 #endif
45 
47 const int64_t numel=1;
48 
49 // Construction/Destruction
50 #ifndef STATE_COMPLEX
51 
53 static REAL angle = 0.2;
54 
55 
57 static void Forcing(const TPZVec<REAL> &x, TPZVec<STATE> &disp) {
58  disp[0] = -(x[1]-0.5)*sin(angle)+(x[0]-0.5)*cos(angle)-(x[0]-0.5);
59  disp[1] = (x[1]-0.5)*cos(angle)+(x[0]-0.5)*sin(angle)-(x[1]-0.5);
60  disp[2] = 0.;
61 }
62 
64  // int index;
65 
66  //Create the Geometric Mesh
67  TPZGeoMesh geo;
68 
69  //Define the output file name
70  std::ofstream output("output.dat");\
71 
72  //Set the basic coordinate nodes
73  double coordstore[4][2] = {{0.,0.},{1.,0.},{1.,1.},{0.,1.}};
74 
75  TPZVec<REAL> coord(3,0.);
76  int i,j;
77 
78  //Set the node coordinates
79  for(i=0; i<4*(numel+1); i++) {
80  for (j=0; j<2; j++) {
81  coord[j] = coordstore[i%4][j];
82  coord[2] = i/4;
83  }
84  // int nodeindex = geo.NodeVec().AllocateNewElement();
85  geo.NodeVec()[i].Initialize(i,coord,geo);
86  }
87 
88  // create the elements
89  TPZGeoEl *gel[numel];
90  TPZVec<int64_t> indices(8);
91 
92  // Set the connectivities
93  for(i=0; i<numel; i++) {
94  // initialize node indexes
95  for(j=0; j<8; j++) indices[j] = 4*i+j;
96  int64_t index;
97  gel[i] = geo.CreateGeoElement(ECube,indices,1, index);
98  }
99  // TPZGeoElBC t3(gel[0],20,-1,geo);
100  // TPZGeoElBC t4(gel[numel-1],25,-2,geo);
101  geo.BuildConnectivity();
102 
103  //Create the computacional mesh
104  TPZCompMesh mesh(&geo);
105 
106  // Insert the materials
107  TPZMaterial * meumat = 0;
108  DebugStop();
109  //new TPZMatHyperElastic(1,1.e5,0.25);
110  mesh.InsertMaterialObject(meumat);
111 
112  //int numeq;
113  TPZVec<int> skyline;
114 
115  // Insert the boundary conditions
116  TPZFMatrix<STATE> val1(3,3,0.),val2(3,1,0.);
117  TPZMaterial * bnd = meumat->CreateBC (meumat,-1,0,val1,val2);
118  mesh.InsertMaterialObject(bnd);
119  bnd = meumat->CreateBC (meumat,-2,0,val1,val2);
120 
122  mesh.InsertMaterialObject(bnd);
123 
124  mesh.AutoBuild();
125  mesh.InitializeBlock();
126 
127 
128  //numeq = mesh.NEquations();
129 
130  // Teste 1 colocar os elementos, inclusive intermedi�rios, como sub elementos
131  TPZSubCompMesh *sub[numel];
132 
133  int64_t index = -1;
134  for (i=0;i<numel;i++){
135  sub[i] = new TPZSubCompMesh(mesh,index);
136  }
137 
138  //Teste 2 - Passar todos os sub elementos para os subelementos
139 
140  for (i=0;i<numel;i++){
141  sub[i]->TransferElement(&mesh,i);
142  }
143 
144  for (i=0;i<numel;i++){
145  sub[i]->MakeAllInternal();
146  // sub[i]->Prints(output);
147  }
148 
149  // mesh.ComputeNodElCon();
150  // mesh.Print(output);
151  TPZNonLinearAnalysis an(&mesh,output);
152 
153 
154  // mesh.Print(output);
155  output.flush();
156  // TPZFMatrix<REAL> *rhs = new TPZFMatrix(skyline);
157  TPZSkylineStructMatrix strskyl(&mesh);
158  an.SetStructuralMatrix(strskyl);
159  an.Solution().Zero();
161  // sol.ShareMatrix(an.Solver());
162  sol.SetDirect(ELDLt);
163  an.SetSolver(sol);
164  // an.Solver().SetDirect(ELDLt);
165  an.IterativeProcess(output,0.00001,5);
166 
167  //mesh.Print(output);
168  sub[0]->LoadSolution();
169  sub[0]->SetName("sub[0]");
170  sub[0]->Print(output);
171  output.flush();
172 
173  return 0;
174 }
175 #else
176 int TPZSubCompMesh::main() {
177  return 0;
178 }
179 #endif
180 
181 
183 fSingularConnect(-1) {
184  SetDimModel(mesh.Dimension());
185  fAnalysis = NULL;
186 
187 }
188 
190 
191  fAnalysis = NULL;
192 }
193 
195 
196  // THIS ROUTINE NEEDS TO INCLUDE THE DELETION OF THE LIST POINTERS
198  if (ref){
199  ref->ResetReference();
200  this->LoadReferences();
201  }
202 #ifdef PZDEBUG
203  ComputeNodElCon();
204 #endif
205  int64_t i, nelem = this->NElements();
206 
207  //deleting subcompmesh
208  for(i=0; i<nelem; i++){
209  TPZCompEl *el = fElementVec[i];
210  TPZSubCompMesh * subc = dynamic_cast<TPZSubCompMesh*>(el);
211  if(subc){
212  delete subc;
213  }
214  }
215 
216 
217  //unwrapping condensed compel
218  for(i=0; i<nelem; i++){
219  TPZCompEl *el = fElementVec[i];
220  TPZCondensedCompEl * cond = dynamic_cast<TPZCondensedCompEl*>(el);
221  if(cond){
222  cond->Unwrap();
223  }
224  }
225 
226  //unwrapping element groups
227  for(i=0; i<nelem; i++){
228  TPZCompEl *el = fElementVec[i];
229  TPZElementGroup * group = dynamic_cast<TPZElementGroup*>(el);
230  if(group){
231  group->Unwrap();
232  }
233  }
234 
235  //deleting interface elements
236  for(i=0; i<nelem; i++){
237  TPZCompEl *el = fElementVec[i];
238  TPZInterfaceElement * face = dynamic_cast<TPZInterfaceElement*>(el);
239  if(face){
240  delete el;
241  }
242  }
243 
244  //deleting other elements
245  for(i=0; i<nelem; i++) {
246  TPZCompEl *el = fElementVec[i];
247  if(!el) continue;
248  delete el;
249  }
250 
251  fElementVec.Resize(0);
253  fConnectVec.Resize(0);
255 
256  MaterialVec().clear();
257 }
258 
259 
261  return Mesh();
262 }
263 
264 
266 
268  int64_t pos1=0, pos2, comind;
269  TPZCompMesh *father = FatherMesh();
270  s1.Push(this);
271  while (father){
272  s1.Push((father));
273  pos1++;
274  father = s1[pos1]->FatherMesh();
275  }
276  pos2 = 0;
277  s2.Push(mesh);
278  father = mesh->FatherMesh();
279  while (father){
280  s2.Push(father);
281  pos2++;
282  father = s2[pos2]->FatherMesh();
283  }
284  if (s1[pos1] != s2[pos2]) return 0;
285  comind=0; //The first mesh is common for all submeshes
286  for (; pos1>=0 && pos2>=0; pos1--, pos2--) {
287  if((s1[pos1])!=(s2[pos2])) {
288  comind=pos1+1;
289  return (s1[comind]);
290  }
291  }
292  return (pos1 >=0 ) ? (s1[pos1+1]) : s2[pos2+1];
293 }
294 
299 {
301  std::map<int64_t,int64_t>::iterator it;
302  for (it=fFatherToLocal.begin(); it!=fFatherToLocal.end(); it++) {
303  fConnectVec[it->second].IncrementElConnected();
304  }
305  if(fSingularConnect >= 0)
306  {
307  fConnectVec[fSingularConnect].IncrementElConnected();
308  }
309  /*
310  int ic;
311  for(ic = 0; ic< fConnectVec.NElements(); ic++)
312  {
313  if(fExternalLocIndex[ic] != -1)
314  {
315  fConnectVec[ic].IncrementElConnected();
316  }
317  }
318  */
319 }
320 
325 {
326  TPZCompMesh::ComputeNodElCon(nelconnected);
327  std::map<int64_t,int64_t>::const_iterator it;
328  for (it=fFatherToLocal.begin(); it!=fFatherToLocal.end(); it++) {
329  nelconnected[it->second]++;
330  }
331  /* int ic;
332  for(ic = 0; ic< fConnectVec.NElements(); ic++)
333  {
334  if(fExternalLocIndex[ic] != -1)
335  {
336  nelconnected[ic]++;
337  }
338  }
339  */
340 }
341 
343  return fConnectIndex.NElements();
344 }
345 
346 int64_t TPZSubCompMesh::ConnectIndex(int i) const{
347  return fConnectIndex[i];
348 }
349 
351  return -1;
352 }
353 
354 
355 //void TPZSubCompMesh::SetMaterial(TPZMaterial * mat){
356 //}
357 
358 int64_t TPZSubCompMesh::NodeIndex(int64_t nolocal, TPZCompMesh *super)
359 {
360  if(super == this) return nolocal;
361  TPZCompMesh *root = CommonMesh(super);
362  if(!root || fExternalLocIndex[nolocal] == -1) return -1;
363  int64_t result = fConnectIndex[fExternalLocIndex[nolocal]];
364 
365  if(root == FatherMesh())
366  {
367  return result;
368  }
369  else
370  {
371  TPZCompMesh *father = FatherMesh();
372  TPZSubCompMesh *sub = dynamic_cast<TPZSubCompMesh *> (father);
373  if(sub)
374  {
375  return sub->NodeIndex(result,super);
376  }
377  else
378  {
379  return -1;
380  }
381  }
382  /* int rootindex = PutinSuperMesh(nolocal,root);
383  return neighbour->GetFromSuperMesh(rootindex,root);*/
384 }
385 
386 int64_t TPZSubCompMesh::AllocateNewConnect(int nshape, int nstate, int order){
387 
388  int64_t connectindex = TPZCompMesh::AllocateNewConnect(nshape, nstate, order);
389  int64_t seqnum = fConnectVec[connectindex].SequenceNumber();
390  int blocksize = nshape*nstate;
391  fBlock.Set(seqnum,blocksize);
392  fConnectVec[connectindex].SetOrder(order,connectindex);
393  int64_t i,oldsize = fExternalLocIndex.NElements();
394 
395  if(oldsize <= connectindex) {
396  fExternalLocIndex.Resize(connectindex+1);
397  for(i=oldsize; i<=connectindex;i++) fExternalLocIndex[i] = -1;
398  } else {
399  fExternalLocIndex[connectindex] = -1;
400  }
401  return connectindex;
402 }
403 
405 
406  int64_t connectindex = TPZCompMesh::AllocateNewConnect(connect);
407  int64_t seqnum = fConnectVec[connectindex].SequenceNumber();
408  int nshape = connect.NShape();
409  int nstate = connect.NState();
410  int blocksize = nshape*nstate;
411  fBlock.Set(seqnum,blocksize);
412  int64_t i,oldsize = fExternalLocIndex.NElements();
413 
414  if(oldsize <= connectindex) {
415  fExternalLocIndex.Resize(connectindex+1);
416  for(i=oldsize; i<=connectindex;i++) fExternalLocIndex[i] = -1;
417  } else {
418  fExternalLocIndex[connectindex] = -1;
419  }
420  return connectindex;
421 }
422 
423 
424 void TPZSubCompMesh::MakeExternal(int64_t local){
425  if(fExternalLocIndex[local] == -1) {
426  //Allocate the dependent nodes of the selected local node in father mesh
427  int64_t extconnect;
428  int64_t lastext = fConnectIndex.NElements();
429  fConnectIndex.Resize(lastext+1);
430  //Allocate the selected local node in father mesh
431  TPZConnect &c = fConnectVec[local];
432  extconnect = FatherMesh()->AllocateNewConnect(c);
433 
434  fConnectIndex[lastext] = extconnect;
435  fExternalLocIndex[local] = lastext;
436  fFatherToLocal[extconnect] = local;
437  TPZConnect::TPZDepend *listdepend = fConnectVec[local].FirstDepend();
438  while(listdepend) {
439  int64_t depindex = listdepend->fDepConnectIndex;
440  MakeExternal(listdepend->fDepConnectIndex);
441  int64_t depextind = fConnectIndex[fExternalLocIndex[depindex]];
442  int64_t r = listdepend->fDepMatrix.Rows();
443  int64_t c = listdepend->fDepMatrix.Cols();
444  FatherMesh()->ConnectVec()[extconnect].AddDependency(extconnect,depextind,listdepend->fDepMatrix,0,0,r,c);
445  fConnectVec[local].RemoveDepend(local,depindex);
446  listdepend = fConnectVec[local].FirstDepend();
447  }
448  } else {
449  if(fConnectVec[local].FirstDepend() ) {
450  std::cout << "TPZSubCompMesh iconsistent data structure !";
451  }
452  }
453 }
454 
455 int64_t TPZSubCompMesh::PutinSuperMesh(int64_t local, TPZCompMesh *super){
456  if(super == this) return local;
457  if(fExternalLocIndex[local] == -1) MakeExternal(local);
458  return FatherMesh()->PutinSuperMesh(fConnectIndex[fExternalLocIndex[local]],super);
459 }
460 
461 int64_t TPZSubCompMesh::GetFromSuperMesh(int64_t superind, TPZCompMesh *super){
462  if(super == this) return superind;
463  if(super != FatherMesh())
464  {
465  superind = FatherMesh()->GetFromSuperMesh(superind,super);
466  super = FatherMesh();
467  }
468  std::map<int64_t,int64_t>::iterator it = fFatherToLocal.find(superind);
469  // int i,nc = fConnectIndex.NElements();
470  // for(i=0; i<nc; i++) if(fConnectIndex[i] == superind) break;
471  // if(i== nc) {
472  if(it == fFatherToLocal.end())
473  {
474  TPZConnect &c = super->ConnectVec()[superind];
475  int64_t gl = AllocateNewConnect(c);
477  fConnectIndex[fConnectIndex.NElements()-1] = superind;
479  fFatherToLocal[superind] = gl;
480  return gl;
481  } else {
482  int64_t j;
483  j = it->second;
484 
485 #ifdef LOG4CXX2
486  if(logger->isDebugEnabled())
487  {
488  std::stringstream sout;
489  sout << "Connect in fathermesh " << superind << " existing connect found : corresponds to connect " << j << " in subcompmesh";
490  LOGPZ_DEBUG(logger,sout.str())
491  }
492 #endif
493 
494  return j;
495  }
496 }
497 
498 void TPZSubCompMesh::Print(std::ostream &out) const {
499 
500  out << "Sub Mesh" << (void *) this;
501  TPZCompEl::Print(out);
502  TPZCompMesh::Print(out);
503  out.flush();
504  int64_t i;
505  for (i=0; i<fConnectVec.NElements(); i++){
506  out << "Node[" << i <<"]\t" << fExternalLocIndex[i];
507  if (fExternalLocIndex[i] != -1) out << " Index in father mesh:\t" << fConnectIndex[fExternalLocIndex[i]];
508  out << std::endl;
509  }
510  std::map<int64_t,int64_t>::const_iterator it;
511  out << "Global to Local map ";
512  for (it=fFatherToLocal.begin(); it!=fFatherToLocal.end(); it++) {
513  out << it->first << '|' << it->second << ' ';
514  }
515  out << std::endl;
516 }
517 
524 {
525  if (fExternalLocIndex[local] == -1) return;
526  TPZCompMesh *father = FatherMesh();
527  int64_t superind = fConnectIndex[fExternalLocIndex[local]];
528 #ifdef PZDEBUG
529  if(father->ConnectVec()[superind].NElConnected() != 1)
530  {
531  std::cout << __PRETTY_FUNCTION__ << " number of elements connected to connect " << superind <<
532  " = " << father->ConnectVec()[superind].NElConnected() << std::endl;
533  }
534  if(father && RootMesh(local) != father) {
535  std::cout << "ERROR";
536  }
537 #endif
538  TPZConnect::TPZDepend *listdepend = father->ConnectVec()[superind].FirstDepend();
539  while(listdepend) {
540  int64_t depfatherindex = listdepend->fDepConnectIndex;
541  int64_t depindexlocal = GetFromSuperMesh(depfatherindex,father);
542  int64_t r = listdepend->fDepMatrix.Rows();
543  int64_t c = listdepend->fDepMatrix.Cols();
544  ConnectVec()[local].AddDependency(local,depindexlocal,listdepend->fDepMatrix,0,0,r,c);
545  father->ConnectVec()[superind].RemoveDepend(superind,depfatherindex);
546  listdepend = father->ConnectVec()[superind].FirstDepend();
547  }
548 }
549 
550 void TPZSubCompMesh::MakeInternal(int64_t local){
551  TransferDependencies(local);
552  int64_t i;
553  int64_t localindex = fExternalLocIndex[local];
554  int64_t fatherindex = fConnectIndex[localindex];
555  for (i=fExternalLocIndex[local]; i<fConnectIndex.NElements()-1; i++){
556  fConnectIndex[i]= fConnectIndex[i+1];
557  }
558  for(i=0; i<fConnectVec.NElements(); i++) {
559  if(fExternalLocIndex[i] != -1 && fExternalLocIndex[i] > localindex) fExternalLocIndex[i]--;
560  }
562  fFatherToLocal.erase(fatherindex);
563  fExternalLocIndex[local]= -1;
564 }
565 
567  TransferDependencies(local);
568  int64_t localindex = fExternalLocIndex[local];
569  int64_t fatherindex = fConnectIndex[localindex];
570  Mesh()->ConnectVec()[fatherindex].RemoveDepend();
571  fConnectIndex[localindex] = -1;
572  fFatherToLocal.erase(fatherindex);
573  fExternalLocIndex[local]= -1;
574 }
575 
576 static void GatherDependency(TPZCompMesh &cmesh, TPZConnect &start, std::set<int64_t> &dependency)
577 {
578  if (!start.HasDependency()) {
579  return;
580  }
581  TPZConnect::TPZDepend *depend = start.FirstDepend();
582  while (depend) {
583  dependency.insert(depend->fDepConnectIndex);
584  TPZConnect &c = cmesh.ConnectVec()[depend->fDepConnectIndex];
585  GatherDependency(cmesh, c, dependency);
586  depend = depend->fNext;
587  }
588 }
589 
591  if (fExternalLocIndex[local] == -1) return this;
592  else return (FatherMesh()->RootMesh(fConnectIndex[fExternalLocIndex[local]]));
593  //return NULL;
594 }
595 
603  // TPZStack<int> stack;
604  //TPZVec<int> nelcon;
605  TPZCompMesh *father = FatherMesh();
606  TPZAdmChunkVector<TPZConnect> &connectvec = father->ConnectVec();
607 #ifdef PZDEBUG
608  //father->ComputeNodElCon();
609 #endif
610 #ifdef LOG4CXX
611  if (logger->isDebugEnabled())
612  {
613  std::stringstream sout;
614  sout << "Connect indexes " << fConnectIndex;
615  LOGPZ_DEBUG(logger,sout.str())
616  }
617 #endif
618 
619  //father->ComputeNodElCon(nelcon);
620  //#ifdef PZDEBUG
621  // int in;
622  // int nn = nelcon.NElements();
623  // for (in=0; in<nn; in++) {
624  // if(father->ConnectVec()[in].NElConnected() != nelcon[in])
625  // {
626  // std::cout << "NelConnected " << in << " " << father->ConnectVec()[in].NElConnected() << " != " << nelcon[in] << std::endl;
627  // }
628  // }
629  //#endif
630  //TPZCompMesh::Print();
631  //father->Print();
632 
633  // cantransfer contains indices in the local mesh
634  std::set<int64_t> cantransfer;
635  std::set<int64_t> delaytransfer;
636  std::map<int64_t,int64_t>::iterator it;
637  for (it=fFatherToLocal.begin(); it!=fFatherToLocal.end(); it++) {
638  // put the candidate nodes in the stack
639  if (father->ConnectVec()[it->first].NElConnected() == 1)
640  {
641  cantransfer.insert(it->second);
642  }
643  }
644  // look for dependent nodes
645  while (cantransfer.size() || delaytransfer.size())
646  {
647  std::set<int64_t>::iterator itset;
648  for (itset = cantransfer.begin(); itset != cantransfer.end(); itset++) {
649  // this doesnt make sense : connectvec is a vector of the father mesh
650  int64_t csubmeshindex = *itset;
651  int64_t elindex = this->fExternalLocIndex[csubmeshindex];
652  int64_t fatherindex = fConnectIndex[elindex];
653 #ifdef PZDEBUG
654  if (fFatherToLocal[fatherindex] != csubmeshindex) {
655  DebugStop();
656  }
657 #endif
658  TPZConnect &con = connectvec[fatherindex];
659 
660  std::set<int64_t> dependset;
661  GatherDependency(*father, con, dependset);
662  for (std::set<int64_t>::iterator it = dependset.begin(); it != dependset.end(); it++) {
663  int64_t connectindex = *it;
664  if (fFatherToLocal.find(connectindex) != fFatherToLocal.end()) {
665  int64_t submeshconnectindex = fFatherToLocal[connectindex];
666  if (cantransfer.find(submeshconnectindex) != cantransfer.end())
667  {
668  delaytransfer.insert(submeshconnectindex);
669  cantransfer.erase(submeshconnectindex);
670  }
671  }
672  }
673  }
674 #ifdef LOG4CXX
675  if (logger->isDebugEnabled())
676  {
677  std::stringstream sout;
678  sout << " connect indexes that will be transferred ";
679  std::copy(cantransfer.begin(), cantransfer.end(), std::ostream_iterator<int64_t>(sout, " "));
680  sout << std::endl;
681  sout << " connect indexes that are delayed ";
682  std::copy(delaytransfer.begin(), delaytransfer.end(), std::ostream_iterator<int64_t>(sout, " "));
683  LOGPZ_DEBUG(logger, sout.str())
684  }
685 #endif
686  for (itset=cantransfer.begin(); itset!=cantransfer.end(); itset++)
687  {
688 #ifdef LOG4CXX
689  if (logger->isDebugEnabled())
690  {
691  std::stringstream sout;
692  int64_t localindex = fExternalLocIndex[*itset];
693  int64_t fatherindex = fConnectIndex[localindex];
694  father->ConnectVec()[fatherindex].Print(*father,sout);
695  sout << "Making the connect index " << *itset << " internal " << " index in the father mesh " << fatherindex << std::endl;
696  sout << "Connect indexes " << fConnectIndex;
697  LOGPZ_DEBUG(logger,sout.str())
698  }
699 #endif
700  MakeInternalFast(*itset);
701  }
702  cantransfer = delaytransfer;
703  delaytransfer.clear();
704  }
705  /*
706  for (i=0;i<fConnectVec.NElements();i++){
707  if (fExternalLocIndex[i]==-1) continue;
708  // put the candidate nodes in the stack
709  if (father->ConnectVec()[fConnectIndex[fExternalLocIndex[i]]].NElConnected() == 1) stack.Push(i);
710  }
711  */
712  // put the independent connects first
713 
714  /*
715  while(stack.NElements()) {
716  int locind = stack.Pop();
717  int can = 0;
718  for(j=0; j<stack.NElements();j++)
719  {
720  int jlocind = stack[j];
721  if (jlocind == locind) continue;
722  TPZConnect &conj = father->ConnectVec()[fConnectIndex[fExternalLocIndex[stack[j]]]];
723  // special procedure when the node has dependencies
724  if (conj.FirstDepend())
725  {
726  TPZConnect::TPZDepend *listdepend = conj.FirstDepend();
727  // if the node upon which locind is dependent is already on the stack, no further analysis required
728  if (listdepend->HasDepend(fConnectIndex[fExternalLocIndex[locind]]))
729  {
730  #ifdef LOG4CXX
731  {
732  std::stringstream sout;
733  sout << "Connect " << locind << " cannot be made internal because of " << jlocind;
734  LOGPZ_DEBUG(logger,sout.str())
735  }
736  #endif
737  break;
738  }
739  }
740  }
741  // no element on the stack is listed as dependent from the current node
742  if (j == stack.NElements())
743  {
744  can=1;
745  }
746  // we found an element in the dependency list. Let s check it first
747  else
748  {
749  // put the node upon which the current node depends in the current position and the dependent node at the end
750  int jlocind = stack[j];
751  stack[j] = locind;
752  stack.Push(jlocind);
753  }
754  // if the node is not internal to the fathermesh, don't put it on the stack
755  if(can && RootMesh(locind) != FatherMesh()) can = 0;
756  if (can)
757  {
758  #ifdef LOG4CXX
759  {
760  std::stringstream sout;
761  sout << "Making the connect index " << locind << " internal";
762  LOGPZ_DEBUG(logger,sout.str())
763  }
764  #endif
765  MakeInternal(locind);
766  }
767  }
768  */
769 
771  int64_t count = 0;
772  for (it=fFatherToLocal.begin(); it!=fFatherToLocal.end(); it++)
773  {
774  fConnectIndex[count] = it->first;
775  fExternalLocIndex[it->second] = count;
776  count++;
777  }
778 #ifdef PZDEBUG
779  if (count != (int64_t)fFatherToLocal.size()) {
780  DebugStop();
781  }
782 #endif
784  //TPZCompMesh::Print();
785  //father->Print();
786  //std::cout.flush();
787 }
788 
789 void TPZSubCompMesh::PotentialInternal(std::list<int64_t> &connectindices) const {
790  int64_t i;
791  TPZCompMesh *father = FatherMesh();
792  TPZVec<int> nelconnected;
793  father->ComputeNodElCon(nelconnected);
794  //TPZCompMesh::Print();
795  //father->Print();
796  for (i=0;i<fConnectVec.NElements();i++){
797  if (fExternalLocIndex[i]==-1)
798  {
799  connectindices.push_back(i);
800  }
801  else
802  {
803  int64_t extcon = this->fConnectIndex[fExternalLocIndex[i]];
804  if(father->ConnectVec()[extcon].NElConnected() == 1)
805  {
806  connectindices.push_back(i);
807  }
808  }
809  }
810 }
811 
812 
813 void TPZSubCompMesh::SetConnectIndex(int inode, int64_t index){
814  fConnectIndex[inode] = index;
815 }
816 
817 int64_t TPZSubCompMesh::TransferElementFrom(TPZCompMesh *mesh, int64_t elindex){
818  if(mesh == this) return elindex;
819 #ifdef PZDEBUG
820  if (! IsAllowedElement(mesh,elindex)) {
821  std::cout <<"TPZSubCompMesh::TransferElementFrom ERROR: trying to transfer an element not allowed" << std::endl;
822  DebugStop();
823  return -1;
824  }
825 #endif
826  if (mesh != FatherMesh()){
827  elindex = FatherMesh()->TransferElementFrom(mesh,elindex);
828  mesh = FatherMesh();
829  }
830 #ifdef PZDEBUG
831  if (CommonMesh(mesh) != mesh){
832  std::cout <<"TPZSubCompMesh::TransferElementFrom ERROR: mesh is not supermesh" << std::endl;
833  DebugStop();
834  return -1;
835  }
836 #endif
837  TPZCompMesh *father = FatherMesh();
838  TPZCompEl *cel = father->ElementVec()[elindex];
839  if (!cel) {
840  std::cout <<"TPZSubCompMesh::TransferElementFrom ERROR: element not existing" << std::endl;
841  DebugStop();
842  return -1;
843  }
844  TPZInterfaceElement *interf = dynamic_cast<TPZInterfaceElement *> (cel);
845  TPZMultiphysicsInterfaceElement *multinterf = dynamic_cast<TPZMultiphysicsInterfaceElement *>(cel);
846  TPZCompEl *left = 0;
847  TPZCompEl *right = 0;
848  if (interf) {
849  left = interf->LeftElement();
850  right = interf->RightElement();
851  }
852  if (multinterf) {
853  left = multinterf->LeftElement();
854  right = multinterf->RightElement();
855  }
856 
857 // if(!interf && !multinterf)
858  {
859  int ncon = cel->NConnects();
860  for (int i=0; i<ncon; i++){
861  int64_t superind = cel->ConnectIndex(i);
862  int64_t subindex = GetFromSuperMesh(superind,father);
863  cel->SetConnectIndex(i,subindex);
864  }
865  }
866 
867 // else
868 // {
869 // int nleftcon = left->NConnects();
870 // {
871 // TPZCompMesh *comm = CommonMesh(left->Mesh());
872 // int ncon = nleftcon;
873 // for (int ic=0; ic<ncon ; ic++) {
874 // int64_t superind = left->ConnectIndex(ic);
875 // int64_t commind = left->Mesh()->PutinSuperMesh(superind, comm);
876 // int64_t subindex = GetFromSuperMesh(commind, comm);
877 // if (multinterf) {
878 // cel->SetConnectIndex(ic, subindex);
879 // }
880 // }
881 // }
882 // {
883 // TPZCompMesh *comm = CommonMesh(right->Mesh());
884 // int ncon = right->NConnects();
885 // for (int ic=0; ic<ncon ; ic++) {
886 // int64_t superind = right->ConnectIndex(ic);
887 // int64_t commind = right->Mesh()->PutinSuperMesh(superind, comm);
888 // int64_t subindex = GetFromSuperMesh(commind, comm);
889 // if (multinterf) {
890 // cel->SetConnectIndex(ic+nleftcon, subindex);
891 // }
892 // }
893 // }
894 // }
895 
896  if(cel->Reference())
897  {
898  TPZMaterial * matfather;
899  matfather = cel->Material();
900  if(!matfather)
901  {
902  // I don't know what to do...
903  DebugStop();
904  }
905  int matid = matfather->Id();
906  TPZMaterial * matthis = FindMaterial(matid);
907 
908  // perform a "shallow copy" of the material
909  if (!matthis) {
910  MaterialVec()[matfather->Id()] = matfather;
911  }
912 
913  // for boundary conditions we need to copy the referred material too
914  TPZBndCond *bnd = dynamic_cast<TPZBndCond *>(matfather);
915  if (bnd) {
916  TPZMaterial *ref = bnd->Material();
917  int refid = ref->Id();
918  TPZMaterial *matthis = FindMaterial(refid);
919  if(!matthis)
920  {
921  MaterialVec()[refid] = ref;
922  }
923  }
924  }
925  cel->SetMesh(this);
926  /*
927  if(cel->Reference())
928  {
929  TPZMaterial * mat = cel->Material();
930  if(!mat)
931  {
932  father->CopyMaterials(*this);
933  }
934  }
935  */
936  // int blocksize=mesh->ConnectVec()[elindex].NDof((TPZCompMesh *)mesh);
937  int64_t newelind = fElementVec.AllocateNewElement();
938  fElementVec[newelind] = cel;
939  cel->SetIndex(newelind);
940  father->ElementVec()[elindex] = 0;
941  father->ElementVec().SetFree(elindex);
942  return newelind;
943 }
944 
945 int64_t TPZSubCompMesh::TransferElementTo(TPZCompMesh *mesh, int64_t elindex){
946 #ifdef PZDEBUG
947  TPZCompMesh *common = CommonMesh(mesh);
948  if ( common!= mesh){
949  std::cout <<"TPZSubCompMesh::TransferElementTo ERROR: mesh is not supermesh" << std::endl;
950  return -1;
951  }
952 #endif
953  if(mesh == this) return elindex;
954 
955  if (mesh != FatherMesh()){
956  FatherMesh()->TransferElementTo(mesh,elindex);
957  }
958 
959 
960  TPZCompMesh *father = FatherMesh();
961  if(elindex >= ElementVec().NElements()){
962  std::cout <<"TPZSubCompMesh::TransferElementTo ERROR: not possible transfer non existing element" << std::endl;
963  return -1;
964  }
965  TPZCompEl *cel = ElementVec()[elindex];
966  if (!cel) {
967  std::cout <<"TPZSubCompMesh::TransferElementTo ERROR: not possible transfer null element" << std::endl;
968  return -1;
969  }
970  int i,ncon = cel->NConnects();
971  for (i=0; i<ncon; i++){
972  int64_t subindex = cel->ConnectIndex(i);
973  MakeExternal(subindex);
974  int64_t superind = fConnectIndex[fExternalLocIndex[subindex]];
975  cel->SetConnectIndex(i,superind);
976  }
977  // int blocksize=father->ConnectVec()[elind].NDof(father);
978  int64_t newelind = father->ElementVec().AllocateNewElement();
979  father->ElementVec()[newelind] = cel;
980  cel->SetMesh(father);
981  cel->SetIndex(newelind);
982  ElementVec()[elindex] = 0;
983  ElementVec().SetFree(elindex);
984  return newelind;
985 }
986 
987 int64_t TPZSubCompMesh::TransferElement(TPZCompMesh *mesh, int64_t elindex){
988  TPZCompMesh *comm = CommonMesh(mesh);
989  int64_t newelind = mesh->TransferElementTo(comm,elindex);
990  int64_t ell=TransferElementFrom(comm,newelind);
991  return ell;
992 }
993 
994 int TPZSubCompMesh::IsAllowedElement(TPZCompMesh *mesh, int64_t elindex){
995  if (CommonMesh(mesh) == mesh){
996  TPZCompMesh *father = this;
997  while(father->FatherMesh() != mesh) {
998  father = father->FatherMesh();
999  }
1000  TPZSubCompMesh *sub = (TPZSubCompMesh *) father;
1001  int64_t index = sub->Index();
1002  return (elindex != index);
1003  }
1004  return 1;
1005 }
1006 
1009 {
1010  if(fAnalysis)
1011  {
1012  }
1013  else
1014  {
1015  std::cout << "The SubCompMesh needs a configured analysis\n";
1016  DebugStop();//this->SetAnalysis();
1017  }
1018  std::set<int> matids = fAnalysis->StructMatrix()->MaterialIds();
1019  if(!NeedsComputing(matids))
1020  {
1021  return;
1022  }
1023  int i=0;
1026  fAnalysis->Assemble();
1027 
1028  //Trying to get a derived Analysis which is a SubMeshAnalysis.
1029  //It could be better done with an abstract class SubMeshAnalysis which defines CondensedSolution method
1030  TPZSubMeshAnalysis * castedAnal = dynamic_cast<TPZSubMeshAnalysis *>(fAnalysis.operator->());
1031  if(!castedAnal)
1032  {
1033  DebugStop();
1034  }
1035 
1036  TPZAutoPointer<TPZMatrix<STATE> > ReducableStiff = castedAnal->Matrix();
1037  if (!ReducableStiff) {
1038  DebugStop();
1039  }
1040  TPZMatRed<STATE, TPZFMatrix<STATE> > *matred = dynamic_cast<TPZMatRed<STATE, TPZFMatrix<STATE> > *> (ReducableStiff.operator->());
1041  if(!matred) DebugStop();
1042 
1043  matred->SetF(fAnalysis->Rhs());
1044 
1045 }
1046 
1049 {
1050  TPZBlock<STATE> &block = Mesh()->Block();
1051  TPZMaterial * mat = MaterialVec().begin()->second;
1052  int nstate = mat->NStateVariables();
1053  int numloadcases = mat->NumLoadCases();
1054  ef.fMesh = Mesh();
1056 
1057  int nelemnodes = NConnects();
1058  ef.fBlock.SetNBlocks(nelemnodes);
1059  for (int i = 0; i < nelemnodes ; i++) {
1060  //int nodeindex = ConnectIndex(i);
1061  int64_t seqnum = Connect(i).SequenceNumber();
1062  ef.fBlock.Set(i,block.Size(seqnum));
1063  }
1064  ef.fConnect.Resize(nelemnodes);
1065 
1066  for(int i=0; i<nelemnodes; ++i){
1067  (ef.fConnect)[i] = ConnectIndex(i);
1068  }
1069 }
1070 
1071 
1072 
1074  if(fAnalysis)
1075  {
1076  }
1077  else
1078  {
1079  std::cout << "The SubCompMesh needs a configured analysis\n";
1080  DebugStop();//this->SetAnalysis();
1081  }
1082  std::set<int> matids = fAnalysis->StructMatrix()->MaterialIds();
1083  if(!NeedsComputing(matids))
1084  {
1085  ek.Reset();
1086  ef.Reset();
1087  return;
1088  }
1089  int i=0;
1092 
1093 
1094 
1095  TPZBlock<STATE> &block = Mesh()->Block();
1096  // TPZFMatrix<REAL> &MeshSol = Mesh()->Solution();
1097  // clean ek and ef
1098 
1099  // int nmeshnodes = fConnectVec.NElements();
1100  int64_t numeq=0, numeq2=0;
1101  //??
1102  int64_t ic;
1103  for (ic=0; ic<fConnectIndex.NElements(); ic++) {
1104  int64_t conindex = fConnectIndex[ic];
1105  TPZConnect &cn = Mesh()->ConnectVec()[conindex];
1106  if (cn.SequenceNumber()<0 || cn.HasDependency()) {
1107  DebugStop();
1108  }
1109  int64_t seqnum = cn.SequenceNumber();
1110  int blsize = Mesh()->Block().Size(seqnum);
1111  numeq2 += blsize;
1112  }
1113  int64_t numeq3 = Mesh()->NEquations();
1114  {
1115  int ftlsize = fFatherToLocal.size();
1116  int ncon = fConnectIndex.NElements();
1117 
1118  if(ftlsize != ncon)
1119  {
1120  std::cout << "Number of connects of the submesh out of sink\n";
1121  DebugStop();
1122  }
1123  }
1124  std::map<int64_t,int64_t>::iterator it;
1125  for (it=fFatherToLocal.begin(); it!=fFatherToLocal.end(); it++)
1126  {
1127  i = it->second;
1128  // for (i=0; i< nmeshnodes; i++){
1129  // if(fExternalLocIndex[i] == -1) {
1130  TPZConnect &df = fConnectVec[i];
1131  if (fExternalLocIndex[i] == -1) {
1132  LOGPZ_ERROR(logger,"fExternalLocIndex and fFatherToLocal are out of sink")
1133  DebugStop();
1134  }
1135  if(df.HasDependency() || !df.NElConnected() || df.SequenceNumber() == -1) continue;
1136  int64_t seqnum = df.SequenceNumber();
1137  numeq += Block().Size(seqnum);
1138  // }
1139  }
1140  int64_t nconstrconnects = 0;
1141  int globeq2 = 0;
1142  for (ic=0; ic<fConnectVec.NElements(); ic++) {
1143  TPZConnect &cn = fConnectVec[ic];
1144  if (cn.HasDependency() || cn.IsCondensed()) {
1145  nconstrconnects++;
1146  }
1147  else if(cn.SequenceNumber() >= 0) {
1148  globeq2 += Block().Size(cn.SequenceNumber());
1149  }
1150 
1151  }
1152 
1153  if (numeq != numeq2 || numeq > numeq3) {
1154  DebugStop();
1155  }
1156 
1157  // check whether the connects are properly enumerated
1158 #ifdef PZDEBUG
1159  int64_t numextconnects = fConnectIndex.NElements();
1160  int64_t nconnects = fConnectVec.NElements();
1161  int64_t numintconnects = nconnects-numextconnects-nconstrconnects;
1162  {
1163  int64_t globeq = TPZCompMesh::NEquations();
1164  if (globeq2 != globeq) {
1165  DebugStop();
1166  }
1167  int64_t numinteq = globeq - numeq;
1168  int64_t in;
1169  // verify whether the block structure is resequenced...
1170  for (in=0; in<nconnects-1; in++) {
1171  int blsize = Block().Size(in);
1172  int64_t pos1 = Block().Position(in);
1173  int64_t pos2 = Block().Position(in+1);
1174  if (pos2-pos1 != blsize) {
1175  DebugStop();
1176  }
1177  }
1178 
1179  int64_t numinteq2 = 0;
1180  if(numintconnects != 0) numinteq2 = Block().Position(numintconnects-1)+Block().Size(numintconnects-1);
1181  if (numinteq != numinteq2) {
1182  DebugStop();
1183  }
1184  //int globeq = TPZCompMesh::NEquations();
1185 
1186  for(in=0; in<nconnects; in++)
1187  {
1188  TPZConnect &df = ConnectVec()[in];
1189  if( ! df.NElConnected() || df.SequenceNumber() == -1) continue;
1190  int64_t seqnum = df.SequenceNumber();
1191  int64_t eq = Block().Position(seqnum);
1192  int blsize = Block().Size(seqnum);
1193  if((eq < numinteq || seqnum < numintconnects) && fExternalLocIndex[in] != -1 )
1194  {
1195  std::stringstream sout;
1196  sout << "Connect " << in << " has equation " << eq << " but is external";
1197  LOGPZ_ERROR(logger,sout.str())
1198  DebugStop();
1199  }
1200  if ((eq >= numinteq || seqnum >= numintconnects) && blsize && !df.HasDependency() && !df.IsCondensed() && fExternalLocIndex[in] == -1) {
1201  std::stringstream sout;
1202  sout << "Connect " << in << " has equation " << eq << " but is internal and has no dependencies ";
1203  df.Print(*this,sout);
1204  LOGPZ_ERROR(logger,sout.str())
1205  DebugStop();
1206  }
1207  if((eq < globeq || seqnum < nconnects-nconstrconnects) && (df.HasDependency() || df.IsCondensed()) )
1208  {
1209  std::stringstream sout;
1210  sout << "Connect " << in << " with dependency was not put at the end of the stack equation " << eq << " global equations " << globeq;
1211  LOGPZ_ERROR(logger,sout.str())
1212  DebugStop();
1213  }
1214  }
1215  }
1216 #endif
1217  //??
1218 
1219  TPZMaterial * mat = MaterialVec().begin()->second;
1220  int nstate = mat->NStateVariables();
1221  int numloadcases = mat->NumLoadCases();
1222  ek.fMesh = Mesh();
1223  ef.fMesh = Mesh();
1224 
1225  ek.fMat.Redim(numeq,numeq);
1226  ef.fMat.Redim(numeq,numloadcases);
1229 
1230  int nelemnodes = NConnects();
1231 
1232  ek.fBlock.SetNBlocks(nelemnodes);
1233  ef.fBlock.SetNBlocks(nelemnodes);
1234  for (i = 0; i < nelemnodes ; i++) {
1235  //int nodeindex = ConnectIndex(i);
1236  int64_t seqnum = Connect(i).SequenceNumber();
1237  ek.fBlock.Set(i,block.Size(seqnum));
1238  ef.fBlock.Set(i,block.Size(seqnum));
1239  }
1240  ek.fConnect.Resize(nelemnodes);
1241  ef.fConnect.Resize(nelemnodes);
1242 
1243  for(i=0; i<nelemnodes; ++i){
1244  (ef.fConnect)[i] = ConnectIndex(i);
1245  (ek.fConnect)[i] = ConnectIndex(i);
1246  }
1247  if (! fAnalysis){
1248  TPZFStructMatrix local(this);
1249  TPZAutoPointer<TPZMatrix<STATE> > stiff = local.CreateAssemble(ef.fMat,NULL);
1250  ek.fMat = *(stiff.operator->());
1251  // TPZStructMatrix::Assemble(ek.fMat,ef.fMat,*this,-1,-1);
1252  }
1253  else{
1254  //if(!fAnalysis->Solver().Matrix())
1255  {
1256  fAnalysis->Run(std::cout);
1257  if(fAnalysis->AmIKilled()){
1258  return;
1259  }
1260  }
1261 
1262  TPZSubMeshFrontalAnalysis *sman = dynamic_cast<TPZSubMeshFrontalAnalysis *> (fAnalysis.operator->());
1263  if(sman)
1264  {
1265  TPZAbstractFrontMatrix<STATE> *frontmat = dynamic_cast<TPZAbstractFrontMatrix<STATE> *> (fAnalysis->Solver().Matrix().operator->());
1266  if(frontmat)
1267  {
1268  sman->SetFront(frontmat->GetFront());
1269  }
1270  }
1271 
1272  //Trying to get a derived Analysis which is a SubMeshAnalysis.
1273  //It could be better done with an abstract class SubMeshAnalysis which defines CondensedSolution method
1274  TPZSubMeshAnalysis * castedAnal = dynamic_cast<TPZSubMeshAnalysis *>(fAnalysis.operator->());
1275  if(castedAnal){
1276  castedAnal->CondensedSolution(ek.fMat,ef.fMat);
1277  }
1278 
1279  TPZSubMeshFrontalAnalysis * castedAnalFrontal = dynamic_cast<TPZSubMeshFrontalAnalysis *>(fAnalysis.operator->());
1280  if(castedAnalFrontal){
1281  castedAnalFrontal->CondensedSolution(ek.fMat,ef.fMat);
1282  }
1283 
1284  if(!castedAnal && !castedAnalFrontal){
1285  DebugStop();
1286  }
1287 
1288  }
1289 #ifdef LOG4CXX
1290  if(logger->isDebugEnabled())
1291  {
1292  std::stringstream sout;
1293  sout << "Substructure stiffness matrix\n";
1294  ek.Print(sout);
1295  sout << "Substructure right hand side\n";
1296  ef.Print(sout);
1297  LOGPZ_DEBUG(logger,sout.str())
1298  }
1299 #endif
1300 
1301  //ek.fMat->Print();
1302 }
1303 
1309 {
1310  TPZFMatrix<STATE> rhs;
1312  TPZSubMeshAnalysis * castedAnal = dynamic_cast<TPZSubMeshAnalysis *>(fAnalysis.operator->());
1313 
1314  if (!castedAnal) {
1315  DebugStop();
1316  }
1317  InitializeEF(ef);
1318  castedAnal->ReducedRightHandSide(ef.fMat);
1319 // TPZCompMesh::CalcResidual(ef);
1320 // ef.PermuteGather(fIndexes);
1321 // fCondensed.SetF(ef.fMat);
1322 // //const TPZFMatrix<REAL> &f1 = fCondensed.F1Red();
1323 // TPZFNMatrix<100,STATE> f1(fCondensed.Dim1(),ef.fMat.Cols());
1324 // fCondensed.F1Red(f1);
1325 // int64_t dim1 = f1.Rows();
1326 // int64_t dim = ef.fMat.Rows();
1327 // int64_t dim0 = dim-dim1;
1328 // for (int64_t i= dim0; i<dim; i++) {
1329 // ef.fMat(i,0) = f1.GetVal(i-dim0,0);
1330 // }
1331 }
1332 
1333 
1334 
1335 void TPZSubCompMesh::SetAnalysisSkyline(int numThreads, int preconditioned, TPZAutoPointer<TPZGuiInterface> guiInterface){
1336  fAnalysis = new TPZSubMeshAnalysis(this);
1337  fAnalysis->SetGuiInterface(guiInterface);
1339 
1340  if(numThreads > 0){
1341  str = new TPZSkylineStructMatrix(this);
1342  str->SetNumThreads(numThreads);
1343  }
1344  else{
1345  str = new TPZSkylineStructMatrix(this);
1346  }
1347 
1348  SaddlePermute();
1349 #ifdef LOG4CXX
1350  if (logger->isDebugEnabled())
1351  {
1352  std::stringstream sout;
1353  Print(sout);
1354  LOGPZ_DEBUG(logger, sout.str())
1355  }
1356 #endif
1358 
1359 
1360 
1361  str->SetNumThreads(numThreads);
1362  int64_t numinternal = NumInternalEquations();
1363  str->EquationFilter().SetMinMaxEq(0, numinternal);
1364  TPZAutoPointer<TPZMatrix<STATE> > mat = str->Create();
1365  str->EquationFilter().Reset();
1366  TPZAutoPointer<TPZMatrix<STATE> > mat2 = mat->Clone();
1367 
1369  TPZStepSolver<STATE> *step = new TPZStepSolver<STATE>(mat);
1370  TPZStepSolver<STATE> *gmrs = new TPZStepSolver<STATE>(mat2);
1371  step->SetReferenceMatrix(mat2);
1372  step->SetDirect(ELDLt);
1373  gmrs->SetGMRES(20, 20, *step, 1.e-20, 0);
1374  TPZAutoPointer<TPZMatrixSolver<STATE> > autostep = step;
1375  TPZAutoPointer<TPZMatrixSolver<STATE> > autogmres = gmrs;
1376  if(preconditioned)
1377  {
1378  fAnalysis->SetSolver(autogmres);
1379  }
1380  else
1381  {
1382  fAnalysis->SetSolver(autostep);
1383  }
1384 
1385 
1386 #ifdef PZDEBUG2
1387  {
1388  TPZFMatrix<REAL> fillin;
1389  int resolution = 100;
1390  ComputeFillIn(resolution,fillin);
1391 #ifdef USING_BOOST
1392  std::string out("matrix_boost.vtk");
1393 #else
1394  std::string out("matrix_native.vtk");
1395 #endif
1396  VisualMatrix(fillin,out);
1397  }
1398 #endif
1399 
1400 }
1401 
1402 void TPZSubCompMesh::SetAnalysisSkyline(int numThreads, int preconditioned, TPZAutoPointer<TPZRenumbering> renumber){
1404  fAnalysis->SetRenumber(renumber);
1405  fAnalysis->SetCompMesh(this, true);
1407 
1408  if(numThreads > 0){
1409  str = new TPZSkylineStructMatrix(this);
1410  str->SetNumThreads(numThreads);
1411  }
1412  else{
1413  str = new TPZSkylineStructMatrix(this);
1414  }
1415 
1416  SaddlePermute();
1417 #ifdef LOG4CXX
1418  if (logger->isDebugEnabled())
1419  {
1420  std::stringstream sout;
1421  Print(sout);
1422  LOGPZ_DEBUG(logger, sout.str())
1423  }
1424 #endif
1426 
1427 
1428 
1429  str->SetNumThreads(numThreads);
1430  int64_t numinternal = NumInternalEquations();
1431  str->EquationFilter().SetMinMaxEq(0, numinternal);
1432  TPZAutoPointer<TPZMatrix<STATE> > mat = str->Create();
1433  str->EquationFilter().Reset();
1434  TPZAutoPointer<TPZMatrix<STATE> > mat2 = mat->Clone();
1435 
1437  TPZStepSolver<STATE> *step = new TPZStepSolver<STATE>(mat);
1438  TPZStepSolver<STATE> *gmrs = new TPZStepSolver<STATE>(mat2);
1439  step->SetReferenceMatrix(mat2);
1440  step->SetDirect(ELDLt);
1441  gmrs->SetGMRES(20, 20, *step, 1.e-20, 0);
1442  TPZAutoPointer<TPZMatrixSolver<STATE> > autostep = step;
1443  TPZAutoPointer<TPZMatrixSolver<STATE> > autogmres = gmrs;
1444  if(preconditioned)
1445  {
1446  fAnalysis->SetSolver(autogmres);
1447  }
1448  else
1449  {
1450  fAnalysis->SetSolver(autostep);
1451  }
1452 
1453 
1454 #ifdef PZDEBUG
1455  {
1456  TPZFMatrix<REAL> fillin;
1457  int resolution = 100;
1458  ComputeFillIn(resolution,fillin);
1459 #ifdef USING_BOOST
1460  std::string out("matrix_boost.vtk");
1461 #else
1462  std::string out("matrix_native.vtk");
1463 #endif
1464  VisualMatrix(fillin,out);
1465  }
1466 #endif
1467 
1468 }
1469 
1471 
1473  fAnalysis->SetGuiInterface(guiInterface);
1474 
1475 #ifdef PZDEBUG2
1476  {
1477  TPZFMatrix<REAL> fillin;
1478  int resolution = 100;
1479  ComputeFillIn(resolution,fillin);
1480 #ifdef USING_BOOST
1481  std::string out("matrix_boost.vtk");
1482 #else
1483  std::string out("matrix_native.vtk");
1484 #endif
1485  VisualMatrix(fillin,out);
1486  }
1487 #endif
1488 
1489  TPZAutoPointer<TPZStructMatrix> fstr = NULL;
1490  if(numThreads){
1492  static_cast<TPZParFrontStructMatrix<TPZFrontNonSym<STATE> > *>(fstr.operator->())
1493  ->SetNumThreads(numThreads+2); //o frontal tem dois threads auxiliares
1494  }
1495  else{
1496  fstr = new TPZFrontStructMatrix<TPZFrontNonSym<STATE> >(this);
1497  }
1498 
1499  fstr->SetNumThreads(numThreads);
1501 
1502  TPZStepSolver<STATE> solver;
1503  solver.SetDirect(ELU);
1504  fAnalysis->SetSolver(solver);
1505 
1506  LOGPZ_DEBUG(logger2, __PRETTY_FUNCTION__)
1508 }
1509 
1515 {
1516  // map from sequence number of the pontentially internal nodes to the node indices
1517  // first the independent nodes, then the dependent nodes
1518  std::map<int64_t,int64_t> independent;
1519  std::list<int64_t> internal;
1520  this->PotentialInternal(internal);
1521 #ifdef LOG4CXX
1522  {
1523  std::stringstream sout;
1524  sout << "Index = " << Index() << " Internal connects ic/seqnum";
1525  std::list<int64_t>::iterator it;
1526  for(it=internal.begin(); it!= internal.end(); it++)
1527  {
1528  sout << *it << "/" << ConnectVec()[*it].SequenceNumber() << " ";
1529  }
1530  if (logger->isDebugEnabled())
1531  {
1532  LOGPZ_DEBUG(logger, sout.str())
1533  }
1534  }
1535 #endif
1536  TPZCompMesh *father = this->FatherMesh();
1537  std::list<int64_t>::iterator it;
1538  for(it=internal.begin(); it!= internal.end(); it++)
1539  {
1540  int64_t locind = *it;
1541  int64_t externallocindex = this->fExternalLocIndex[locind];
1542  if(externallocindex > 0)
1543  {
1544  int64_t superind = fConnectIndex[externallocindex];
1545  if(father->ConnectVec()[superind].FirstDepend())
1546  {
1547  }
1548  else
1549  {
1550  independent[ConnectVec()[locind].SequenceNumber()] = locind;
1551  }
1552  }
1553  else if (!ConnectVec()[locind].FirstDepend())
1554  {
1555  independent[ConnectVec()[locind].SequenceNumber()] = locind;
1556  }
1557  }
1558 #ifdef LOG4CXX
1559  {
1560  std::stringstream sout;
1561  sout << "Mesh Address " << (void *) this << " Index = " << Index() << " \nIndependent connect sequence numbers and indices ";
1562  std::map<int64_t,int64_t>::iterator mapit;
1563  for(mapit=independent.begin(); mapit!= independent.end(); mapit++)
1564  {
1565  sout << "[" << mapit->first << " , " << mapit->second << "] ";
1566  }
1567  if (logger->isDebugEnabled())
1568  {
1569  LOGPZ_DEBUG(logger, sout.str())
1570  }
1571  }
1572 #endif
1573  permute.Resize(0);
1574  permute.Resize(fConnectVec.NElements(),-1);
1575 
1576  int64_t count = 0;
1577  std::map<int64_t,int64_t>::iterator mapit;
1578  for(mapit=independent.begin(); mapit!=independent.end(); mapit++)
1579  {
1580  permute[mapit->first] = count++;
1581  }
1582 #ifdef LOG4CXX
1583  {
1584  std::stringstream sout;
1585  sout << "Index = " << Index() << " Permutation vector 1 " << permute;
1586  if (logger->isDebugEnabled())
1587  {
1588  LOGPZ_DEBUG(logger, sout.str())
1589  }
1590  }
1591 #endif
1592  std::map<int64_t,int64_t> seqmap;
1593  int64_t ind;
1594  for(ind=0; ind < fConnectVec.NElements(); ind++)
1595  {
1596  int64_t seqnum = fConnectVec[ind].SequenceNumber();
1597  if(seqnum == -1) continue;
1598  seqmap[seqnum]=ind;
1599  }
1600  for(mapit=seqmap.begin(); mapit!=seqmap.end(); mapit++)
1601  {
1602  if(permute[mapit->first] == -1) permute[mapit->first] = count++;
1603  }
1604 #ifdef LOG4CXX
1605  {
1606  std::stringstream sout;
1607  sout << "Index = " << Index() << " Permutation vector 2 " << permute;
1608  if (logger->isDebugEnabled())
1609  {
1610  LOGPZ_DEBUG(logger, sout.str())
1611  }
1612  }
1613 #endif
1614 }
1615 
1621 {
1622  this->ComputePermutationInternalFirst(permute);
1623  LOGPZ_DEBUG(logger, "Permuting")
1624  Permute(permute);
1625 }
1626 
1628  //compute number of internal nodes -> numinternal
1629  // TPZCompMesh::Print();
1630 
1631  ComputeNodElCon();
1632 
1633  int64_t i=0, numinternal=0, numconstraints = 0, numexternal=0;
1634  //int countinternal=0
1635  int64_t countconstraint=0;
1636  int64_t nconnects = fConnectVec.NElements();
1637  std::set<int64_t> internalseqnum;
1638  //std::cout << "fExternalLocIndex\n";
1639  //for(i=0; i<nconnects; i++) std::cout << fExternalLocIndex[i] << ' ';
1640  //std::cout << std::endl;
1641  for(i=0;i<nconnects; i++){
1642  if (fExternalLocIndex[i]==-1){
1643  // which are not free and which are not constrained
1644  TPZConnect &no = fConnectVec[i];
1645 
1646  if(no.NElConnected() == 0) continue;
1647  //se n�o tiver elemento conectado tambe'm
1648  if(no.HasDependency() || no.IsCondensed())
1649  {
1650  numconstraints++;
1651  }
1652  else
1653  {
1654  numinternal+= 1;
1655  internalseqnum.insert(no.SequenceNumber());
1656  }
1657  }
1658  else
1659  {
1660  numexternal++;
1661  }
1662  }
1663  countconstraint = numinternal+numexternal;
1664  // initialize a counter for internal nodes
1665  i=0;
1666  TPZManVector<int64_t> permute(nconnects);
1667  for (i=0;i<nconnects;i++) permute[i] = i;
1668  std::set<int64_t>::iterator it;
1669  int64_t seqnum = 0;
1670  for(it=internalseqnum.begin(); it!=internalseqnum.end(); it++)
1671  {
1672  permute[*it] = seqnum++;
1673  }
1674 #ifdef LOG4CXX
1675  if(logger->isDebugEnabled())
1676  {
1677  std::stringstream sout;
1678  sout << " numinternal " << numinternal << " numconstraints " << numconstraints << " numexternal " << numexternal << std::endl;
1679  // sout << " permute so far " << permute;
1680  LOGPZ_DEBUG(logger,sout.str())
1681  }
1682 #endif
1683  // loop over all nodes
1684  for (i=0;i<fConnectVec.NElements();i++){
1685  // take seqnum = the sequencenumber of the node
1686  //TPZConnect &no = fConnectIndex[fExternalLocIndex[i]];
1687  TPZConnect &no = fConnectVec[i];
1688  seqnum = no.SequenceNumber();
1689  // if the node is free or constrained
1690  // if (no.HasDependency() || no.NElConnected() == 0) {
1691  // //->set permute[sequnum] to itself
1692  // continue;
1693  // }
1694  // if the node is internal
1695  if (fExternalLocIndex[i] == -1){
1696  //-> set permute[sequnum] to counter
1697  // ->set permute[seqnum] = fExternalConnectIndex+numinternal
1698  if(no.HasDependency() || no.IsCondensed())
1699  {
1700  permute[seqnum] = countconstraint;
1701  countconstraint++;
1702  }
1703  else
1704  {
1705  }
1706  }
1707  // if the node is external
1708  else
1709  {
1710  permute [seqnum] = fExternalLocIndex[i]+numinternal;
1711  // end loop
1712  }
1713  }
1714 #ifdef LOG4CXX
1715  if(logger->isDebugEnabled())
1716  {
1717  std::stringstream sout;
1718  sout << "Index = " << " Permutations " << permute << std::endl;
1719  std::set<int64_t> permval;
1720  permval.insert(&permute[0], (&permute[permute.size()-1]+1));
1721  sout << " Number of distinct values in permute " << permval.size();
1722  LOGPZ_DEBUG(logger,sout.str())
1723  }
1724 #endif
1725  Permute(permute);
1726 }
1727 
1729 
1730  int64_t i=0;
1731  int64_t seqnumext;
1732  int64_t seqnumint;
1733  // int numinteq = NumInternalEquations();
1734  int size;
1735  TPZFMatrix<STATE> &sol = Mesh()->Solution();
1736 
1737  for (i=0;i<fConnectVec.NElements(); i++) {
1738  if (fExternalLocIndex[i] != -1) {
1740  TPZConnect &noint = fConnectVec[i];
1741  seqnumext = noext.SequenceNumber();
1742  size = (Mesh()->Block()).Size(seqnumext);
1743  seqnumint = noint.SequenceNumber();
1744  int64_t posext = Mesh()->Block().Position(seqnumext);
1745  int64_t posint = fBlock.Position(seqnumint);
1746  int l;
1747  for(l=0;l<size;l++) {
1748  fSolution(posint+l,0) = sol(posext+l,0);
1749  }
1750  }
1751  }
1752 // fSolution.Print(std::cout);
1753 
1756 }
1757 
1761 TPZVec<STATE> TPZSubCompMesh::IntegrateSolution(const std::string &varname, const std::set<int> &matids)
1762 {
1763  return TPZCompMesh::Integrate(varname,matids);
1764 }
1765 
1766 
1768 {
1769  int64_t nel = this->NElements();
1770 #ifdef LOG4CXX
1771  if (logger->isDebugEnabled()) {
1772  std::stringstream sout;
1773  fSolution.Print("SubMeshSol",sout);
1774  LOGPZ_DEBUG(logger, sout.str())
1775  }
1776 #endif
1777  for (int64_t iel = 0; iel < nel; iel++) {
1778  TPZCompEl *cel = this->Element(iel);
1779  if (!cel) {
1780  continue;
1781  }
1783  }
1784 }
1785 
1786 
1788  TPZCompMesh::Skyline(skyline);
1789  skyline.Resize(NumInternalEquations());
1790 }
1791 
1793  int64_t nmeshnodes = fConnectVec.NElements();
1794  int64_t numeq=0;
1795  //??
1796 
1797  int64_t i;
1798  for (i=0; i< nmeshnodes; i++){
1799  if(fExternalLocIndex[i] == -1) {
1800  TPZConnect &df = fConnectVec[i];
1801  if(df.HasDependency() || df.IsCondensed() || !df.NElConnected() || df.SequenceNumber() == -1) continue;
1802  int dfsize = df.NShape()*df.NState();
1803 #ifdef PZDEBUG
1804  int64_t seqnum = df.SequenceNumber();
1805  int blsize = Block().Size(seqnum);
1806  if (blsize != dfsize) {
1807  DebugStop();
1808  }
1809 #endif
1810  numeq += dfsize;
1811  }
1812  }
1813  return numeq;
1814 
1815 }
1816 
1817 
1818 REAL TPZSubCompMesh::CompareElement(int var, char *matname){
1819  return CompareMesh(var,matname);
1820 }
1821 
1822 
1824 {
1826 }
1827 
1828 /*
1829  void TPZSubCompMesh::GetExternalConnectIndex (TPZVec<int> &extconn){
1830  int i;
1831  int count = 0;
1832  for(i=0; i<fExternalLocIndex.NElements(); i++){
1833  if (fExternalLocIndex[i] != -1) count++;
1834  }
1835  extconn.Resize(count);
1836 
1837  count=0;
1838  for(i=0; i<fExternalLocIndex.NElements(); i++){
1839  if (fExternalLocIndex[i] != -1){
1840  extconn[count] = i;
1841  count++;
1842  }
1843  }
1844  return;
1845  }
1846  */
1847 
1848 
1853  return Hash("TPZSubCompMesh") ^ TPZCompMesh::ClassId() << 1 ^ TPZCompEl::ClassId() << 2;
1854 }
1855 
1856 #ifndef BORLAND
1857 template class TPZRestoreClass< TPZSubCompMesh>;
1858 #endif
1859 
1863 void TPZSubCompMesh::Write(TPZStream &buf, int withclassid) const
1864 {
1865  //std::map<int, TPZMaterial *> matmap = MaterialVec();
1866  //MaterialVec().clear();
1867  TPZCompEl::Write(buf,withclassid);
1868  TPZCompMesh::Write(buf,0);
1869  //MaterialVec() = matmap;//AQUIFRAN
1870  const std::map<int, TPZMaterial *> &matmap = fMaterialVec;
1871  TPZManVector<int> matindex(matmap.size(),-1);
1872  int count=0;
1873  for (std::map<int,TPZMaterial *>::const_iterator it = matmap.begin(); it != matmap.end(); it++) {
1874  matindex[count++] = it->first;
1875  }
1876  buf.Write( matindex);
1877  buf.Write(fConnectIndex);
1878  buf.Write(fExternalLocIndex);
1879  buf.Write(fFatherToLocal);
1880  buf.Write(&fSingularConnect,1);
1881 }
1882 
1886 void TPZSubCompMesh::Read(TPZStream &buf, void *context)
1887 {
1888  TPZCompEl::Read(buf,context);
1889  TPZCompMesh::Read(buf,Mesh()->Reference());
1890  TPZCompMesh *mesh = (TPZCompMesh *) context;
1891  TPZManVector<int> matindex;
1892  buf.Read( matindex);
1893  int sz = matindex.size();
1894  for (int im=0; im<sz; im++) {
1895  MaterialVec()[matindex[im]] = mesh->MaterialVec()[matindex[im]];
1896  }
1897  buf.Read(fConnectIndex);
1898  buf.Read(fExternalLocIndex);
1899  buf.Read( fFatherToLocal);
1900  buf.Read(&fSingularConnect,1);
1901 }
1902 
1904  TPZSolVec &sol, TPZGradSolVec &dsol,TPZFMatrix<REAL> &axes){
1905  PZError << __PRETTY_FUNCTION__ << " - ERROR! This method is not implemented\n";
1906 }
1907 
1909  const TPZFMatrix<REAL> &axes, TPZSolVec &sol, TPZGradSolVec &dsol){
1910  PZError << __PRETTY_FUNCTION__ << " - ERROR! This method is not implemented\n";
1911 }
1912 
1914  TPZVec<REAL> &normal,
1915  TPZSolVec &leftsol, TPZGradSolVec &dleftsol,TPZFMatrix<REAL> &leftaxes,
1916  TPZSolVec &rightsol, TPZGradSolVec &drightsol,TPZFMatrix<REAL> &rightaxes){
1917  PZError << __PRETTY_FUNCTION__ << " - ERROR! This method is not implemented\n";
1918 }
1919 
1920 
1928 {
1929  int64_t nel = fElementVec.NElements();
1930  int64_t iel;
1931  for(iel=0; iel<nel; iel++)
1932  {
1933  TPZCompEl *cel = fElementVec[iel];
1934  if(!cel) continue;
1935  cel->CreateGraphicalElement(graphmesh, dimension);
1936  }
1937 }
1938 
1942 bool TPZSubCompMesh::NeedsComputing(const std::set<int> &matids)
1943 {
1944  std::set<int> meshmatids;
1945  if(! matids.size())
1946  {
1947  return true;
1948  }
1949  int numtrue=0,numfalse=0;
1950  // loop over the elements
1951  int64_t iel, nelem = ElementVec().NElements();
1952  for(iel=0; iel<nelem; iel++)
1953  {
1954  TPZCompEl *cel = ElementVec()[iel];
1955  if(!cel) continue;
1956  TPZMaterial * mat = cel->Material();
1957  if(!mat)
1958  {
1959  TPZSubCompMesh *submesh = dynamic_cast<TPZSubCompMesh *> (cel);
1960  if(submesh)
1961  {
1962  bool result = submesh->NeedsComputing(matids);
1963  if(result) numtrue++;
1964  else numfalse++;
1965  }
1966  }
1967  else
1968  {
1969  int matid = mat->Id();
1970  meshmatids.insert(matid);
1971  if(matids.find(matid) != matids.end())
1972  {
1973  numtrue++;
1974  }
1975  else {
1976  numfalse++;
1977  }
1978 
1979  }
1980  }
1981  {
1982  std::stringstream sout;
1983  sout << "Material ids contained in the mesh ";
1984  std::set<int>::iterator it;
1985  for(it = meshmatids.begin(); it != meshmatids.end(); it++)
1986  {
1987  sout << *it << " ";
1988  }
1989  sout << std::endl << "Material ids which should be computed ";
1990  std::set<int>::const_iterator it2;
1991  for(it2= matids.begin(); it2 != matids.end(); it2++)
1992  {
1993  sout << *it2 << " ";
1994  }
1995  sout << std::endl;
1996  LOGPZ_DEBUG(logger, sout.str())
1997  }
1998  if(numtrue && numfalse)
1999  {
2000  std::stringstream sout;
2001  sout << "A substructure should have either all elements computable or not numtrue " << numtrue << " numfalse " << numfalse;
2002  LOGPZ_WARN(logger,sout.str())
2003  }
2004  if(numtrue)
2005  {
2006  return true;
2007  }
2008  else {
2009  return false;
2010  }
2011  return false;
2012 }
2013 
2015 {
2016  // all elements of fConnectIndex should be found in fFatherToSub map
2017  if (fConnectIndex.NElements() != fFatherToLocal.size()) {
2018  DebugStop();
2019  }
2020  int64_t numberexternal = fConnectIndex.NElements();
2021  int64_t i;
2022  for (i=0; i<numberexternal; i++) {
2023  if (fFatherToLocal.find(fConnectIndex[i]) == fFatherToLocal.end()) {
2024  DebugStop();
2025  }
2026  }
2027  // the number of external connects in the fExternalLocIndex should be size also
2028  int64_t nel = fExternalLocIndex.NElements();
2029  int64_t numext = 0;
2030  for (i=0; i<nel; i++) {
2031  if (fExternalLocIndex[i] != -1) {
2032  numext++;
2033  }
2034  }
2035  if (numext != numberexternal) {
2036  DebugStop();
2037  }
2038  std::map<int64_t,int64_t>::iterator it;
2039  for (it=fFatherToLocal.begin(); it!=fFatherToLocal.end(); it++) {
2040  if (fExternalLocIndex[it->second] == -1) {
2041  DebugStop();
2042  }
2043  }
2044  return true;
2045 }
2046 
2048 bool TPZSubCompMesh::HasMaterial(const std::set<int> &materialids) const
2049 {
2050  int nel = NElements();
2051  for (int el=0; el<nel ; el++) {
2052  TPZCompEl *cel = ElementVec()[el];
2053  if (!cel) {
2054  continue;
2055  }
2056  bool has_material_Q = cel->HasMaterial(materialids);
2057  if (has_material_Q) {
2058  return true;
2059  }
2060  }
2061  return false;
2062 }
2063 
2065 {
2066  if (fSingularConnect == -1) {
2067  return 0;
2068  }
2069  int64_t seqnum = fConnectVec[fSingularConnect].SequenceNumber();
2070  return fBlock.Size(seqnum);
2071 
2072 }
2073 
2074 
2076 void TPZSubCompMesh::SetNumberRigidBodyModes(int nrigid, unsigned char lagrange)
2077 {
2078  if (fSingularConnect == -1) {
2079  int nshape = nrigid;
2080  int nstate = 1;
2081  int order = 1;
2082  fSingularConnect = AllocateNewConnect(nshape,nstate,order);
2083  fConnectVec[fSingularConnect].IncrementElConnected();
2084  fConnectVec[fSingularConnect].SetLagrangeMultiplier(lagrange);
2085  int64_t extind = FatherMesh()->AllocateNewConnect(nshape,nstate,order);
2086  FatherMesh()->ConnectVec()[extind].IncrementElConnected();
2087  FatherMesh()->ConnectVec()[extind].SetLagrangeMultiplier(lagrange);
2088  int64_t next = fConnectIndex.NElements();
2089  fConnectIndex.Resize(next+1);
2090  fConnectIndex[next] = extind;
2092  fFatherToLocal[extind] = fSingularConnect;
2093  ExpandSolution();
2094  }
2095  else if(fSingularConnect != -1 && nrigid >0 ) {
2096  int64_t seqnum = fConnectVec[fSingularConnect].SequenceNumber();
2097  fConnectVec[fSingularConnect].SetLagrangeMultiplier(lagrange);
2098  fBlock.Set(seqnum,nrigid);
2099  ExpandSolution();
2100  int64_t extind = fExternalLocIndex[fSingularConnect];
2101  TPZCompMesh *fathermesh = FatherMesh();
2102  if (fathermesh && extind < 0) {
2103  DebugStop();
2104  }
2105  while (fathermesh && extind > 0) {
2106  seqnum = fathermesh->ConnectVec()[extind].SequenceNumber();
2107  fathermesh->ConnectVec()[extind].SetLagrangeMultiplier(lagrange);
2108  fathermesh->Block().Set(seqnum, nrigid);
2109  fathermesh->ExpandSolution();
2110  TPZSubCompMesh *subfather = dynamic_cast<TPZSubCompMesh *> (fathermesh);
2111  if (subfather) {
2112  extind = subfather->fExternalLocIndex[extind];
2113  }
2114  fathermesh = fathermesh->FatherMesh();
2115  }
2116  }
2117  else {
2118  // not implemented yet
2119  DebugStop();
2120  }
2121 }
2122 
2124 void TPZSubCompMesh::BuildCornerConnectList(std::set<int64_t> &connectindexes) const
2125 {
2126  int nel = NElements();
2127  for (int el=0; el<nel ; el++) {
2128  TPZCompEl *cel = ElementVec()[el];
2129  if (!cel) {
2130  continue;
2131  }
2132  std::set<int64_t> locconind;
2133  cel->BuildCornerConnectList(locconind);
2134  std::set<int64_t>::iterator it;
2135  for (it=locconind.begin(); it != locconind.end(); it++) {
2136  int64_t index = *it;
2137  int64_t extlocindex = fExternalLocIndex[index];
2138  if (extlocindex != -1) {
2139  int64_t cornerind = fConnectIndex[extlocindex];
2140  connectindexes.insert(cornerind);
2141  }
2142  }
2143  }
2144 }
2145 
2147 int64_t TPZSubCompMesh::InternalIndex(int64_t IndexinFather)
2148 {
2149  if (fFatherToLocal.find(IndexinFather) == fFatherToLocal.end()) {
2150  DebugStop();
2151  }
2152  return fFatherToLocal[IndexinFather];
2153 }
2154 
2155 void TPZSubCompMesh::EvaluateError(std::function<void(const TPZVec<REAL> &loc,TPZVec<STATE> &val,TPZFMatrix<STATE> &deriv)> fp,
2156  TPZVec<REAL> &errors, bool store_errors){
2157 
2158  fAnalysis->SetExact(fp);
2159  fAnalysis->PostProcessError(errors,store_errors);
2160  int NErrors = errors.size();
2161  if(store_errors)
2162  {
2163  int64_t index = Index();
2164  TPZFMatrix<STATE> &elvals = Mesh()->ElementSolution();
2165  if (elvals.Cols() < NErrors) {
2166  std::cout << "The element solution of the mesh should be resized before EvaluateError\n";
2167  DebugStop();
2168  }
2169  for (int ier=0; ier <NErrors; ier++) {
2170  elvals(index,ier) = errors[ier];
2171  }
2172  }
2173 
2174 }
int64_t NElements() const
Number of computational elements allocated.
Definition: pzcmesh.h:169
TPZFMatrix< STATE > & ElementSolution()
Access method for the element solution vectors.
Definition: pzcmesh.h:225
TPZGeoMesh * Reference() const
Returns a pointer to the geometrical mesh associated.
Definition: pzcmesh.h:209
Represents a graphical mesh used for post processing purposes. Post processing.
Definition: pzgraphmesh.h:34
virtual int64_t GetFromSuperMesh(int64_t superind, TPZCompMesh *super) override
Gets an external connection from the supermesh - Supermesh is one mesh who contains the analised subm...
Definition: pzsubcmesh.cpp:461
TPZAutoPointer< TPZAnalysis > fAnalysis
Pointer to submesh analysis object. Defines the resolution type.
Definition: pzsubcmesh.h:42
int AllocateNewElement()
Makes more room for new elements.
Definition: pzadmchunk.h:184
virtual TPZCompMesh * FatherMesh() const
Get the father meshes stack.
Definition: pzcmesh.h:337
virtual TPZCompMesh * RootMesh(int64_t local) override
Returns the rootmesh who have the specified connection.
Definition: pzsubcmesh.cpp:590
virtual void BuildCornerConnectList(std::set< int64_t > &connectindexes) const =0
adds the connect indexes associated with base shape functions to the set
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
TPZAutoPointer< TPZMatrix< STATE > > Matrix()
Definition: pzsmanal.h:39
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 SetName(const std::string &nm)
Set the mesh name.
Definition: pzcmesh.cpp:232
Definition: pzmatrix.h:52
Contains the TPZParSkylineStructMatrix class which defines parallel structural matrix for skyline mat...
TPZCompEl * LeftElement() const
Returns the left element from the element interface.
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 VisualMatrix(TPZFMatrix< TVar > &matrix, const std::string &outfilename)
This function calls the function that create a Data Explorer file or VTK file that allow to visualiz...
TPZFMatrix< STATE > fSolution
Solution vector.
Definition: pzcmesh.h:72
virtual void SetConnectIndex(int inode, int64_t index)=0
Set the index i to node inode.
virtual TPZVec< STATE > IntegrateSolution(const std::string &varname, const std::set< int > &matids) override
Compute the integral of a variable defined by the string if the material id is included in matids...
void Unwrap()
put the elements in the element group back in the mesh and delete the element group ...
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
Analysis procedure to computational sub mesh. Analysis.
Definition: pzsmanal.h:20
Contains the TPZFrontStructMatrix class which responsible for a interface among Finite Element Packag...
TPZAdmChunkVector< TPZCompEl * > fElementVec
List of pointers to elements.
Definition: pzcmesh.h:60
virtual void CalcResidual(TPZElementMatrix &ef) override
Computes the element right hand side.
TPZManVector< int64_t > fConnectIndex
Pointer to external location index of the connection.
Definition: pzsubcmesh.h:46
virtual int64_t AllocateNewConnect(int nshape, int nstate, int order) override
Virtual Method to allocate new connect.
Definition: pzsubcmesh.cpp:386
int64_t NEquations()
This computes the number of equations associated with non-restrained nodes.
Definition: pzcmesh.cpp:721
clarg::argInt dimension("-d", "Matrices dimension M x M", 1000)
TPZStack< int64_t > fConnect
Vector of pointers to TPZConnect objects.
Definition: pzelmat.h:39
virtual int64_t AllocateNewConnect(int nshape, int nstate, int order)
Returns an index to a new connect.
Definition: pzcmesh.h:693
virtual void CalcStiff(TPZElementMatrix &ek, TPZElementMatrix &ef) override
Calculates the submesh stiffness matrix.
void SaddlePermute()
Put the sequence number of the pressure connects after the seq number of the flux connects...
Definition: pzcmesh.cpp:2328
TPZFNMatrix< 1000, STATE > fMat
Pointer to a blocked matrix object.
Definition: pzelmat.h:41
Contains declaration of TPZCompEl class which defines the interface of a computational element...
Implements a matrix stored in a frontal decomposition scheme. Frontal.
Contains the declaration of the TPZElementGroup class, which implements an computational element whic...
void InitializeEF(TPZElementMatrix &ef)
Initialize the datastructure of ef.
virtual void CreateGraphicalElement(TPZGraphMesh &graphmesh, int dimension) override
Creates corresponding graphical element(s) if the dimension matches graphical elements are used to ge...
int ClassId() const override
Returns the unique identifier for reading/writing objects to streams.
Definition: pzcmesh.cpp:1969
Defines step solvers class. Solver.
Definition: pzmganalysis.h:17
void SetSolver(TPZMatrixSolver< STATE > &solver)
Set solver matrix.
virtual void Read(TPZStream &buf, void *context) override
Reads the element data from a stream.
virtual TPZEquationFilter & EquationFilter()
access method for the equation filter
void SetF(const TPZFMatrix< TVar > &F)
Copies the F vector in the internal data structure.
Definition: pzmatred.cpp:140
TPZVec< STATE > Integrate(const std::string &varname, const std::set< int > &matids)
Integrate the variable name over the mesh.
Definition: pzcmesh.cpp:2162
virtual int Dimension() const override
Virtual Method! See declaration in TPZCompEl class.
Definition: pzsubcmesh.cpp:350
virtual void Print(std::ostream &out=std::cout) const override
Prints the submesh information on the specified device/file out.
Definition: pzsubcmesh.cpp:498
virtual void Assemble()
Assemble the stiffness matrix and load vector.
Definition: pzanalysis.cpp:304
virtual int NStateVariables() const =0
Returns the number of state variables associated with the material.
virtual void ComputeSolution(TPZVec< REAL > &qsi, TPZSolVec &sol, TPZGradSolVec &dsol, TPZFMatrix< REAL > &axes) override
Computes solution and its derivatives in the local coordinate qsi. SHOULD NEVER BE CALLED...
virtual void CreateGraphicalElement(TPZGraphMesh &graphmesh, int dimension)
Creates corresponding graphical element(s) if the dimension matches graphical elements are used to ge...
Definition: pzcompel.h:785
int64_t fSingularConnect
Number of rigid body modes expected by the internal matrix inversion.
Definition: pzsubcmesh.h:55
Class which groups elements to characterize dense 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.
REAL val(STATE &number)
Returns value of the variable.
Definition: pzartdiff.h:23
virtual int64_t TransferElementFrom(TPZCompMesh *mesh, int64_t elindex)
Transfer one element from a specified mesh to the current submesh.
Definition: pzcmesh.h:356
#define LOGPZ_WARN(A, B)
Define log for warnings.
Definition: pzlog.h:91
void ComputeFillIn(int64_t resolution, TPZFMatrix< REAL > &fillin)
This method will fill the matrix passed as parameter with a representation of the fillin of the globa...
Definition: pzcmesh.cpp:1869
void SetGuiInterface(TPZAutoPointer< TPZGuiInterface > gui)
Defines gui interface object.
Definition: pzanalysis.h:108
TPZBlock< STATE > fBlock
Block structure associated with fMat.
Definition: pzelmat.h:43
virtual void TransferMultiphysicsElementSolution() override
Transfer multiphysics solution.
TPZMatrixSolver< STATE > & Solver()
Get the solver matrix.
Definition: pzanalysis.h:367
virtual void Resize(const int64_t newsize, const T &object)
Resizes the vector object.
Definition: pzmanvector.h:426
Contains the declaration of the VisualMatrix functions to VTK and DX packages.
void LoadReferences()
Map this grid in the geometric grid.
Definition: pzcmesh.cpp:482
virtual ~TPZSubCompMesh()
Destructor.
Definition: pzsubcmesh.cpp:194
void ReducedRightHandSide(TPZFMatrix< STATE > &rhs)
compute the reduced right hand side using the current stiffness. Abort if there is no stiffness compu...
Definition: pzsmanal.cpp:145
void SetMesh(TPZCompMesh *mesh)
Sets the grid of the element.
Definition: pzcompel.cpp:281
void SetExact(std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &result, TPZFMatrix< STATE > &deriv)> f)
Sets the pointer of the exact solution function.
Definition: pzanalysis.h:360
virtual void TransferMultiphysicsElementSolution()
Definition: pzcompel.h:415
void PermuteInternalFirst(TPZVec< int64_t > &permute)
Permutes the potentially internal connects to the first on the list Respect the previous order of th...
Contains the TPZFStructMatrix class which implements Full Structural Matrices.
virtual int ClassId() const override
Define the class id associated with the class.
Definition: pzcompel.cpp:1129
Implements SkyLine Structural Matrices. Structural Matrix.
TPZFMatrix< STATE > & Solution()
Returns the solution matrix.
Definition: pzanalysis.h:177
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...
virtual void LoadElementReference() override
This method will load the elements of the mesh in their corresponding geometric elements.
TPZAutoPointer< TPZStructMatrix > StructMatrix()
Returns a reference to the structural matrix.
Definition: pzanalysis.h:182
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
virtual void SetConnectIndex(int inode, int64_t index) override
Changes the current node index -inode- to the specified node- index.
Definition: pzsubcmesh.cpp:813
virtual TPZCompMesh * FatherMesh() const override
Returns the current submesh father mesh.
Definition: pzsubcmesh.cpp:260
void MakeExternal(int64_t local)
Changes an local internal connection to a external connection in the father mesh. ...
Definition: pzsubcmesh.cpp:424
void Print(std::ostream &out)
Definition: pzelmat.cpp:47
int Zero() override
Makes Zero all the elements.
Definition: pzfmatrix.h:651
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
int64_t TransferElementTo(TPZCompMesh *mesh, int64_t elindex) override
Transfers one element from a submesh to another mesh.
Definition: pzsubcmesh.cpp:945
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
TPZSubCompMesh()
Default constructor.
Definition: pzsubcmesh.cpp:189
TPZCompEl * Element(int64_t iel)
Definition: pzcmesh.h:185
void CondensedSolution(TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef)
CondensedSolution: returns the condensed stiffness matrix - ek - and the condensed solution vector - ...
Definition: pzsmanal.cpp:121
void LoadSolution(const TPZFMatrix< STATE > &sol)
Given the solution of the global system of equations, computes and stores the solution for the restri...
Definition: pzcmesh.cpp:441
void Resize(const int newsize)
Increase the size of the chunk vector.
Definition: pzadmchunk.h:280
void Unwrap()
unwrap the condensed element from the computational element and delete the condensed element ...
sin
Definition: fadfunc.h:63
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 Read(TPZStream &buf, void *context) override
Read the element data from a stream.
Definition: pzcompel.cpp:736
void Push(const T object)
Pushes a copy of the object on the stack.
Definition: pzstack.h:80
virtual void Skyline(TPZVec< int64_t > &skyline)
This method computes the skyline of the system of equations.
Definition: pzcmesh.cpp:787
virtual void SkylineInternal(TPZVec< int64_t > &skyline)
This method computes the skyline of the system of equations.
int64_t fDepConnectIndex
Definition: pzconnect.h:69
virtual void LoadSolution()
Load the solution into the computable grid.
Definition: pzanalysis.cpp:441
virtual void Write(const bool val)
Definition: TPZStream.cpp:8
Structure to reference dependency.
Definition: pzconnect.h:67
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
virtual void EvaluateError(std::function< void(const TPZVec< REAL > &loc, TPZVec< STATE > &val, TPZFMatrix< STATE > &deriv)> func, TPZVec< REAL > &errors, bool store_error) override
Performs an error estimate on the elemen.
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
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 Read(TPZStream &buf, void *context) override
Read the element data from a stream.
Definition: pzcmesh.cpp:2006
Contains declaration of TPZInterfaceElement class which computes the contribution over an interface b...
virtual void ExpandSolution()
Adapt the solution vector to new block dimensions.
Definition: pzcmesh.cpp:396
virtual void SetNumThreads(int n)
#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
TPZCompMesh * Mesh() const
Return a pointer to the grid of the element.
Definition: pzcompel.cpp:288
Contains the TPZSkylineStructMatrix class which implements SkyLine Structural Matrices.
Definition: pzmatrix.h:52
int Dimension() const
Returns the dimension of the simulation.
Definition: pzcmesh.h:148
virtual void ComputeNodElCon()
Compute the number of elements connected to each connect object.
Definition: pzcmesh.cpp:668
Contains the declaration of multiphysic interface class.
virtual void LoadSolution() override
Load the father mesh solution to all submesh connects - (internal and external).
std::map< int64_t, int64_t > fFatherToLocal
Maps indicating the correspondence between the connect index of the father mesh and de local connect ...
Definition: pzsubcmesh.h:52
int64_t InternalIndex(int64_t IndexinFather)
return the index in the subcompmesh of a connect with index within the father
TPZMaterial * FindMaterial(int id)
Find the material with identity id.
Definition: pzcmesh.cpp:297
virtual int64_t PutinSuperMesh(int64_t local, TPZCompMesh *super)
Put an local connection in the supermesh - Supermesh is one mesh who contains the analised submesh...
Definition: pzcmesh.cpp:1512
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
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
void MakeInternalFast(int64_t local)
Marks the connect to be local.
Definition: pzsubcmesh.cpp:566
virtual int ClassId() const override
Returns the unique identifier for reading/writing objects to streams.
virtual int64_t TransferElementTo(TPZCompMesh *mesh, int64_t elindex)
Transfer one element from a submesh to another mesh.
Definition: pzcmesh.h:363
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 TransferDependencies(int64_t local)
Transfer the dependency list of a connect. This will make the dependency disappear for the correspond...
Definition: pzsubcmesh.cpp:523
TPZCompMesh * fMesh
Definition: pzelmat.h:36
Contains TPZSkyline class which implements a skyline storage format.
TPZDepend * fNext
Definition: pzconnect.h:71
TPZAdmChunkVector< TPZConnect > & ConnectVec()
Return a reference to the connect pointers vector.
Definition: pzcmesh.h:198
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
virtual bool HasMaterial(const std::set< int > &materialids) const override
Verifies if the material associated with the element is contained in the set.
#define LOGPZ_ERROR(A, B)
Define log for errors (cout)
Definition: pzlog.h:93
Analysis for substructuring. Use a frontal matrix. Analysis.
int SetNBlocks(const int num_of_blocks)
Sets number of blocks on diagonal matrix.
Definition: pzblock.cpp:91
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
TPZBlock< STATE > fBlock
Block structure to right construction of the stiffness matrix and load vector.
Definition: pzcmesh.h:80
void PotentialInternal(std::list< int64_t > &connectindices) const
Puts the nodes which can be transferred in an ordered list.
Definition: pzsubcmesh.cpp:789
virtual TPZConnect & Connect(int i) const
Returns a pointer to the ith node.
Definition: pzcompel.cpp:298
void CompactDataStructure(CompactScheme type=CompactScheme::ALWAYS)
Sets the method to compact the data structure based on the.
Definition: pzadmchunk.h:222
virtual void IterativeProcess(std::ostream &out, REAL tol, int numiter, bool linesearch=false, bool checkconv=false)
It process a Newton&#39;s method to solve the non-linear problem.
void CondensedSolution(TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef)
CondensedSolution: returns the condensed stiffness matrix - ek - and the condensed solution vector - ...
virtual TPZFront< TVar > & GetFront()=0
void SetDimModel(int dim)
Set de dimension of the domain of the problem.
Definition: pzcmesh.h:139
unsigned char NState() const
Number of state variables associated with the connect.
Definition: pzconnect.h:146
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
bool NeedsComputing(const std::set< int > &matids) override
Verifies if any element needs to be computed corresponding to the material ids.
int32_t Hash(std::string str)
Definition: TPZHash.cpp:10
std::map< int,TPZMaterial *> & MaterialVec()
Returns a reference to the material pointers vector.
Definition: pzcmesh.h:203
void PermuteExternalConnects()
Optimizes the connections positions on block.
virtual int64_t ConnectIndex(int i) const =0
Returns the index of the ith connectivity of the element.
Computes the contribution over an interface between two discontinuous elements. Computational Element...
bool AmIKilled()
Returns if the process was canceled through gui interface.
Definition: pzanalysis.h:118
This class associates an element matrix with the coeficients of its contribution in the global stiffn...
Definition: pzelmat.h:30
static REAL angle
Angle in radians to test.
Definition: pzsubcmesh.cpp:53
virtual int NConnects() const =0
Returns the number of nodes of the element.
Implements Full Structural Matrices. Structural Matrix.
Definition: pzfstrmatrix.h:19
TPZFNMatrix< 50, REAL > fDepMatrix
Definition: pzconnect.h:70
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
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
virtual void BuildCornerConnectList(std::set< int64_t > &connectindexes) const override
adds the connect indexes associated with base shape functions to the set
virtual void CleanUpUnconnectedNodes()
Delete the nodes which have no elements connected to them.
Definition: pzcmesh.cpp:498
virtual void Assemble() override
Assemble the stiffness matrix in locally kept datastructure.
Implements a simple substructuring of a linear system of equations, composed of 4 submatrices...
Definition: pzmatred.h:34
int IsAllowedElement(TPZCompMesh *mesh, int64_t elindex)
Verifies the transfer possibility of the connection elindex from the mesh to the current submesh...
Definition: pzsubcmesh.cpp:994
virtual void AssembleResidual()
Assemble the load vector.
Definition: pzanalysis.cpp:288
TPZGeoEl * Reference() const
Return a pointer to the corresponding geometric element if such exists, return 0 otherwise.
Definition: pzcompel.cpp:1137
TPZCompEl * LeftElement() const
Returns the left element from the element interface.
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 REAL CompareMesh(int var, char *matname)
This method will initiate the comparison between the current computational mesh and the mesh which is...
Definition: pzcmesh.cpp:1522
void SetForcingFunction(TPZAutoPointer< TPZFunction< STATE > > fp)
Sets a procedure as source function for the material.
Definition: TPZMaterial.h:368
void SetMinMaxEq(int64_t mineq, int64_t maxeq)
Define as equacoes ativas de [mineq, maxeq)
Computes the contribution over an interface between two discontinuous elements. Computational Element...
virtual void ComputeNodElCon() override
Computes the number of elements connected to each connect object.
Definition: pzsubcmesh.cpp:298
void Reset(TPZCompMesh *mesh=NULL, MType type=Unknown)
Reset the data structure.
Definition: pzelmat.h:59
TPZDepend * FirstDepend()
Definition: pzconnect.h:296
Contains TPZSubMeshFrontalAnalysis class which implements the analysis for substructuring.
static int main()
Static function for validation tests.
Definition: pzsubcmesh.cpp:63
virtual void MakeAllInternal() override
Makes all mesh connections internal mesh connections.
Definition: pzsubcmesh.cpp:602
virtual REAL CompareElement(int var, char *matname) override
This method will initiate the comparison between the current computational mesh and the mesh which is...
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
Contains TPZSolver class which defines a abstract class of solvers which will be used by matrix class...
TPZAutoPointer< TPZMatrix< TVar > > Matrix() const
Returns a pointer to TPZMatrix<>
Definition: pzsolve.h:138
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
virtual void Print(std::ostream &out=std::cout) const
Prints element data.
Definition: pzcompel.cpp:322
void SetAnalysisFrontal(int numThreads, TPZAutoPointer< TPZGuiInterface > guiInterface)
Sets the analysis type.
TPZCompEl * RightElement() const
Returns the right element from the element interface.
int NElConnected() const
Returns fNElConnected.
Definition: pzconnect.h:264
int64_t NumInternalEquations()
Computes the number of internal equations.
int NumLoadCases()
returns the number of load cases for this material object
Definition: TPZMaterial.h:186
int Id() const
Definition: TPZMaterial.h:170
Contains TPZStepSolver class which defines step solvers class.
int64_t Cols() const
Returns number of cols.
Definition: pzmatrix.h:809
Derived class from TPZAnalysis implements non linear analysis (Newton&#39;s method). Analysis.
void SetFree(int index)
Indicate an element as free.
Definition: pzadmchunk.h:199
void Write(TPZStream &buf, int withclassid) const override
Save the element data to a stream.
Definition: pzcompel.cpp:726
pthread_cond_t cond
Definition: numatst.cpp:318
virtual void Print(std::ostream &out) const
Definition: pzmatrix.h:253
Contains TPZSubMeshAnalysis class which implements the analysis procedure to computational sub mesh...
std::map< int, TPZMaterial *> fMaterialVec
Map of pointers to materials.
Definition: pzcmesh.h:65
Defines the interface for saving and reading data. Persistency.
Definition: TPZStream.h:50
int64_t NElements() const
Returns the number of elements of the vector.
Definition: pzvec.h:190
TPZCompEl * RightElement() const
Returns the right element from the element interface.
static void Forcing(const TPZVec< REAL > &x, TPZVec< STATE > &disp)
Defining function force (external to material) .
Definition: pzsubcmesh.cpp:57
void SetNumberRigidBodyModes(int nrigid, unsigned char lagrange=0)
Set the number of rigid body modes associated with the internal degrees of freedom.
void SetAnalysisSkyline(int numThreads, int preconditioned, TPZAutoPointer< TPZGuiInterface > guiInterface)
Condense the internal equations using a skyline symetric matrix the preconditioned argument indicates...
TPZFMatrix< STATE > & Rhs()
Returns the load vector.
Definition: pzanalysis.h:174
void SetDirect(const DecomposeType decomp)
void InitializeBlock()
Resequence the block object, remove unconnected connect objects and reset the dimension of the soluti...
Definition: pzcmesh.cpp:391
bool VerifyDatastructureConsistency()
Method to verify that the datastructures are consistent.
virtual void MakeInternal(int64_t local) override
Makes a specified connection a internal mesh connection.
Definition: pzsubcmesh.cpp:550
void SetFront(TPZFront< STATE > &front)
Sets the front matrix.
virtual void Write(TPZStream &buf, int withclassid) const override
Saves the element data to a stream.
void ResetReference()
Resets all load references in elements and nodes.
Definition: pzgmesh.cpp:197
void Print(const TPZCompMesh &mesh, std::ostream &out=std::cout)
Print the information for the connect element.
Definition: pzconnect.cpp:67
Is a structural matrix with parallel techniques included. Structural Matrix Frontal.
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
TPZManVector< int64_t > fExternalLocIndex
Indexes of the external connections.
Definition: pzsubcmesh.h:50
TPZFlopCounter cos(const TPZFlopCounter &orig)
Returns the cosine in radians and increments the counter of the Cosine.
Definition: pzreal.h:514
virtual TPZCompMesh * CommonMesh(TPZCompMesh *mesh) override
Gives the commom father mesh of the specified mesh and the current submesh.
Definition: pzsubcmesh.cpp:265
virtual void PostProcessError(TPZVec< REAL > &, bool store_error=true, std::ostream &out=std::cout)
Compute the local error over all elements and global errors in several norms and print out without ca...
Definition: pzanalysis.cpp:573
Contains the declaration of the TPZCondensedCompEl class, which implements an computational element w...
void SetRenumber(TPZAutoPointer< TPZRenumbering > renumber)
Change the renumbering scheme.
Definition: pzanalysis.h:142
virtual void SetCompMesh(TPZCompMesh *mesh, bool mustOptimizeBandwidth)
Set the computational mesh of the analysis.
Definition: pzanalysis.cpp:115
TPZFMatrix< STATE > & Solution()
Access the solution vector.
Definition: pzcmesh.h:219
Contains the TPZParFrontStructMatrix class which is a structural matrix with parallel techniques incl...
virtual TPZMatrix< STATE > * Create() override
Definition: pzstrmatrix.cpp:65
static void GatherDependency(TPZCompMesh &cmesh, TPZConnect &start, std::set< int64_t > &dependency)
Definition: pzsubcmesh.cpp:576
void Write(TPZStream &buf, int withclassid) const override
Save the element data to a stream.
Definition: pzcmesh.cpp:1975
TPZAdmChunkVector< TPZConnect > fConnectVec
List of pointers to nodes.
Definition: pzcmesh.h:62
Implements an interface to register a class id and a restore function. Persistence.
Definition: TPZSavable.h:150
#define PZError
Defines the output device to error messages and the DebugStop() function.
Definition: pzerror.h:15
void Reset()
Reset method.
virtual int64_t PutinSuperMesh(int64_t local, TPZCompMesh *super) override
Puts an local connection in the supermesh - Supermesh is one mesh who contains the analised submesh...
Definition: pzsubcmesh.cpp:455
TPZMaterial * Material() const
Definition: pzbndcond.h:263
virtual void Print(std::ostream &out=std::cout) const
Prints mesh data.
Definition: pzcmesh.cpp:236
void SetGMRES(const int64_t numiterations, const int numvectors, const TPZMatrixSolver< TVar > &pre, const REAL tol, const int64_t FromCurrent)
virtual void Read(bool &val)
Definition: TPZStream.cpp:91
virtual int64_t GetFromSuperMesh(int64_t superind, TPZCompMesh *super)
Get an external connection from the supermesh - Supermesh is one mesh who contains the analised subme...
Definition: pzcmesh.cpp:1517
int64_t TransferElementFrom(TPZCompMesh *mesh, int64_t elindex) override
Transfers one element from a specified mesh to the current submesh.
Definition: pzsubcmesh.cpp:817
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
Contains TPZNonLinearAnalysis class which implements the non linear analysis.
This class implements a reference counter mechanism to administer a dynamically allocated object...