NeoPZ
TPZAcademicGeoMesh.cpp
Go to the documentation of this file.
1 //
2 // TPZAcademicGeoMesh.cpp
3 // PZ
4 //
5 // Created by Philippe Devloo on 5/29/16.
6 //
7 //
8 
9 #include "TPZAcademicGeoMesh.h"
10 
11 #include "pzshapelinear.h"
12 
13 #include "pzgengrid.h"
14 #include "TPZExtendGridDimension.h"
15 #include "TPZVTKGeoMesh.h"
16 
17 #include "pzcheckmesh.h"
18 
19 #include "TPZMaterial.h"
20 #include "pzbndcond.h"
21 #include "pzpoisson3d.h"
22 
23 #include "pzanalysis.h"
24 #include "pzstepsolver.h"
25 
26 #include "TPZRefPatternTools.h"
27 
28 
29 #include <stdio.h>
30 
31 #include "pzlog.h"
32 
33 #include "pzgeoelbc.h"
34 
36 #ifdef LOG4CXX
37 static LoggerPtr logger(Logger::getLogger("pz.academicmesh"));
38 #endif
39 
40 
41 
42 TPZAcademicGeoMesh::TPZAcademicGeoMesh() : fMeshType(EHexa), fBCNumbers(6,-1), fShouldDeform(false), fNumberElements(1)
43 {
44  REAL coord[8][3] = {
45  // {0,0,0},
46  // {1,0,0},
47  // {1,1,0},
48  // {0,1,0},
49  // {0,0,1},
50  // {1,0,1},
51  // {1,1,1},
52  // {0,1,1}
53  {-0.5,-0.5,-0.5},
54  {2,-1,-1},
55  {1.1,1.1,-0.1},
56  {-1,2,-1},
57  {-1,-1,2},
58  {1.2,-0.2,1.2},
59  {2,2,2},
60  {-0.5,1.5,1.5}
61  };
63  TPZManVector<int64_t,8> indices(8);
64  for (int i=0; i<8; i++) {
65  indices[i] = i;
66  for (int c=0; c<3; c++) {
67  fDeformed.NodeVec()[i].SetCoord(c, coord[i][c]);
68  }
69  }
70  int64_t index;
71  fDeformed.CreateGeoElement(ECube, indices, 1, index);
72 }
73 
76 {
77  int64_t nnodes = gmesh.NodeVec().NElements();
78  TPZManVector<REAL,3> xbefore(3),xafter(3);
79  for (int64_t nod=0; nod<nnodes; nod++) {
80  gmesh.NodeVec()[nod].GetCoordinates(xbefore);
81  for (int i=0; i<3; i++) {
82  xbefore[i] = 2.*xbefore[i]-1.;
83  }
84  fDeformed.ElementVec()[0]->X(xbefore, xafter);
85  gmesh.NodeVec()[nod].SetCoord(xafter);
86  }
87 }
88 
89 static int pyramid[2][5]=
90 {
91  {0,1,2,3,4},
92  {4,5,6,7,2}
93 };
94 static int tetraedra[2][4]=
95 {
96  {1,2,5,4},
97  {4,7,3,2}
98 };
99 static int tetraedra_2[6][4]=
100 {
101  {1,2,5,4},
102  {4,7,3,2},
103  {0,1,2,4},
104  {0,2,3,4},
105  {4,5,6,2},
106  {4,6,7,2}
107 };
108 
109 static int prism[2][6] =
110 {
111  {0,1,2,4,5,6},
112  {0,2,3,4,6,7}
113 };
114 
115 
117 {
118  int64_t nelem = fNumberElements;
119  gmesh->NodeVec().Resize((nelem+1)*(nelem+1)*(nelem+1));
120  for (int64_t i=0; i<=nelem; i++) {
121  for (int64_t j=0; j<=nelem; j++) {
122  for (int64_t k=0; k<=nelem; k++) {
124  x[0] = k*1./nelem;
125  x[1] = j*1./nelem;
126  x[2] = i*1./nelem;
127  gmesh->NodeVec()[i*(nelem+1)*(nelem+1)+j*(nelem+1)+k].Initialize(x, *gmesh);
128  }
129  }
130  }
131 }
132 
134 {
135  int64_t nelem = fNumberElements;
136  int MaterialId = fMaterialId;
137  TPZGeoMesh *gmesh = new TPZGeoMesh;
138  const int dim = 3;
139  gmesh->SetDimension(dim);
140  GenerateNodes(gmesh);
141 
142  for (int64_t i=0; i<nelem; i++) {
143  for (int64_t j=0; j<nelem; j++) {
144  for (int64_t k=0; k<nelem; k++) {
145  TPZManVector<int64_t,8> nodes(8,0);
146  nodes[0] = k*(nelem+1)*(nelem+1)+j*(nelem+1)+i;
147  nodes[1] = k*(nelem+1)*(nelem+1)+j*(nelem+1)+i+1;
148  nodes[2] = k*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i+1;
149  nodes[3] = k*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i;
150  nodes[4] = (k+1)*(nelem+1)*(nelem+1)+j*(nelem+1)+i;
151  nodes[5] = (k+1)*(nelem+1)*(nelem+1)+j*(nelem+1)+i+1;
152  nodes[6] = (k+1)*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i+1;
153  nodes[7] = (k+1)*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i;
154 #ifdef LOG4CXX
155  if(logger->isDebugEnabled())
156  {
157  std::stringstream sout;
158  sout << "Pyramid and tetrahedral nodes " << nodes;
159  LOGPZ_DEBUG(logger, sout.str())
160  }
161 #endif
162  for (int el=0; el<2; el++)
163  {
164  TPZManVector<int64_t,5> elnodes(5);
165  for (int il=0; il<5; il++) {
166  elnodes[il] = nodes[pyramid[el][il]];
167  }
168  int64_t index;
169  gmesh->CreateGeoElement(EPiramide, elnodes, MaterialId, index);
170  elnodes.resize(4);
171  for (int il=0; il<4; il++) {
172  elnodes[il] = nodes[tetraedra[el][il]];
173  }
174  gmesh->CreateGeoElement(ETetraedro, elnodes, MaterialId, index);
175  }
176  }
177  }
178  }
179  gmesh->BuildConnectivity();
180 // AddBoundaryElements(gmesh);
182  if (fShouldDeform) {
183  DeformGMesh(*gmesh);
184  }
185  return gmesh;
186 }
187 
189 
190  int64_t n = fNumberElements;
191  int MaterialId = fMaterialId;
192  TPZGeoMesh *gmesh = new TPZGeoMesh;
193  GenerateNodes(gmesh);
194  gmesh->SetDimension(3);
195 
196  int64_t s = n + 1;
197  REAL dir = 1.0;
198  REAL dx = 0.5/n;
199  int64_t n_nodes = gmesh->NNodes();
200  TPZStack<int> perm;
201  perm.push_back(0);
202  perm.push_back(1);
203  perm.push_back(3);
204  perm.push_back(2);
205  perm.push_back(4);
206  perm.push_back(5);
207  perm.push_back(7);
208  perm.push_back(6);
209 
210  TPZManVector<int,8> indexes(8);
211  for (int k = 0; k < n; k++) {
212  for (int j = 0; j < n; j++) {
213  for (int i = 0; i < n; i++) {
214 
215  int stride = j*s + k*s*s;
216  indexes[0] = i + stride;
217  indexes[1] = i + 1 + stride;
218  indexes[2] = i + s + stride;
219  indexes[3] = i + 1 + s + stride;
220  for (int l = 0; l < 4; l++) { // Vertical stride
221  indexes[l+4] = indexes[l]+s*s;
222  }
223 
224  bool is_even_Q = (i+j+k)%2==0;
225  if (is_even_Q) {
226 
227  TPZManVector<int64_t,8> nodes(8);
228  for (int inode=0; inode<8; inode++) {
229  nodes[inode] = gmesh->NodeVec()[indexes[perm[inode]]].Id();
230  }
231 
232  // Adding center point
234  gmesh->NodeVec()[nodes[0]].GetCoordinates(x);
235  for (int i = 0; i < 3; i++) {
236  x[i]+=dx*dir;
237  }
238 
239  gmesh->NodeVec().Resize(n_nodes+1);
240  gmesh->NodeVec()[n_nodes].Initialize(x, *gmesh);
241  int center_id = gmesh->NodeVec()[n_nodes].Id();
242  n_nodes++;
243 
244  // inserting Pyramids
245  int64_t index;
246  TPZManVector<int64_t,5> el_nodes(5);
247 
248  el_nodes[0] = nodes[0];
249  el_nodes[1] = nodes[1];
250  el_nodes[2] = nodes[2];
251  el_nodes[3] = nodes[3];
252  el_nodes[4] = center_id;
253  gmesh->CreateGeoElement(EPiramide, el_nodes, MaterialId, index);
254 
255  el_nodes[0] = nodes[0];
256  el_nodes[1] = nodes[3];
257  el_nodes[2] = nodes[7];
258  el_nodes[3] = nodes[4];
259  el_nodes[4] = center_id;
260  gmesh->CreateGeoElement(EPiramide, el_nodes, MaterialId, index);
261 
262  el_nodes[0] = nodes[0];
263  el_nodes[1] = nodes[1];
264  el_nodes[2] = nodes[5];
265  el_nodes[3] = nodes[4];
266  el_nodes[4] = center_id;
267  gmesh->CreateGeoElement(EPiramide, el_nodes, MaterialId, index);
268 
269  el_nodes[0] = nodes[1];
270  el_nodes[1] = nodes[2];
271  el_nodes[2] = nodes[6];
272  el_nodes[3] = nodes[5];
273  el_nodes[4] = center_id;
274  gmesh->CreateGeoElement(EPiramide, el_nodes, MaterialId, index);
275 
276  el_nodes[0] = nodes[3];
277  el_nodes[1] = nodes[2];
278  el_nodes[2] = nodes[6];
279  el_nodes[3] = nodes[7];
280  el_nodes[4] = center_id;
281  gmesh->CreateGeoElement(EPiramide, el_nodes, MaterialId, index);
282 
283  el_nodes[0] = nodes[4];
284  el_nodes[1] = nodes[5];
285  el_nodes[2] = nodes[6];
286  el_nodes[3] = nodes[7];
287  el_nodes[4] = center_id;
288  gmesh->CreateGeoElement(EPiramide, el_nodes, MaterialId, index);
289 
290 
291  }else{
292  TPZManVector<int64_t,8> el_nodes(8);
293  int64_t index;
294  for (int inode=0; inode<8; inode++) {
295  el_nodes[inode] = gmesh->NodeVec()[indexes[perm[inode]]].Id();
296  }
297  gmesh->CreateGeoElement(ECube, el_nodes, MaterialId, index);
298  }
299 
300  }
301  }
302  }
303 
304  gmesh->BuildConnectivity();
305  AddBoundaryElements(gmesh);
306  if (fShouldDeform) {
307  DeformGMesh(*gmesh);
308  }
309 
310  return gmesh;
311 }
312 
314 {
315  int64_t nelem = fNumberElements;
316  int MaterialId = fMaterialId;
317  TPZGeoMesh *gmesh = new TPZGeoMesh;
318  GenerateNodes(gmesh);
319  gmesh->SetDimension(3);
320 
321  for (int64_t i=0; i<nelem; i++) {
322  for (int64_t j=0; j<nelem; j++) {
323  for (int64_t k=0; k<nelem; k++) {
324  TPZManVector<int64_t,8> nodes(8,0);
325  nodes[0] = k*(nelem+1)*(nelem+1)+j*(nelem+1)+i;
326  nodes[1] = k*(nelem+1)*(nelem+1)+j*(nelem+1)+i+1;
327  nodes[2] = k*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i+1;
328  nodes[3] = k*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i;
329  nodes[4] = (k+1)*(nelem+1)*(nelem+1)+j*(nelem+1)+i;
330  nodes[5] = (k+1)*(nelem+1)*(nelem+1)+j*(nelem+1)+i+1;
331  nodes[6] = (k+1)*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i+1;
332  nodes[7] = (k+1)*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i;
333 #ifdef LOG4CXX
334  if(logger->isDebugEnabled())
335  {
336  std::stringstream sout;
337  sout << "Tetrahedral nodes " << nodes;
338  LOGPZ_DEBUG(logger, sout.str())
339  }
340 #endif
341  for (int el=0; el<6; el++)
342  {
343  TPZManVector<int64_t,4> elnodes(4);
344  int64_t index;
345  for (int il=0; il<4; il++) {
346  elnodes[il] = nodes[tetraedra_2[el][il]];
347  }
348  gmesh->CreateGeoElement(ETetraedro, elnodes, MaterialId, index);
349  }
350  }
351  }
352  }
353  gmesh->BuildConnectivity();
354  AddBoundaryElements(gmesh);
355  if (fShouldDeform) {
356  DeformGMesh(*gmesh);
357  }
358  return gmesh;
359 }
360 
361 
363 {
364  int64_t nelem = fNumberElements;
365  int MaterialId = fMaterialId;
366  TPZGeoMesh *gmesh = new TPZGeoMesh;
367  gmesh->SetDimension(3);
368  GenerateNodes(gmesh);
369 
370  for (int64_t k=0; k<nelem; k++) {
371  for (int64_t j=0; j<nelem; j++) {
372  for (int64_t i=0; i<nelem; i++) {
373  TPZManVector<int64_t,8> nodes(8,0);
374  nodes[0] = k*(nelem+1)*(nelem+1)+j*(nelem+1)+i;
375  nodes[1] = k*(nelem+1)*(nelem+1)+j*(nelem+1)+i+1;
376  nodes[2] = k*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i+1;
377  nodes[3] = k*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i;
378  nodes[4] = (k+1)*(nelem+1)*(nelem+1)+j*(nelem+1)+i;
379  nodes[5] = (k+1)*(nelem+1)*(nelem+1)+j*(nelem+1)+i+1;
380  nodes[6] = (k+1)*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i+1;
381  nodes[7] = (k+1)*(nelem+1)*(nelem+1)+(j+1)*(nelem+1)+i;
382 #ifdef LOG4CXX
383  if (logger->isDebugEnabled())
384  {
385  std::stringstream sout;
386  sout << "Cube nodes " << nodes;
387  LOGPZ_DEBUG(logger, sout.str())
388  }
389 #endif
390  int64_t index;
391  gmesh->CreateGeoElement(ECube, nodes, MaterialId, index);
392  }
393  }
394  }
395  gmesh->BuildConnectivity();
396  AddBoundaryElements(gmesh);
397  if (fShouldDeform) {
398  DeformGMesh(*gmesh);
399  }
400  return gmesh;
401 }
402 
405 {
406  int64_t nel = mesh->NElements();
407  for(int64_t el=0; el<nel; el++) {
408  TPZGeoEl *gel = mesh->ElementVec()[el];
409  int nsides = gel->NSides();
410  for(int is=0; is<nsides; is++) {
411  TPZGeoElSide gelside(gel,is);
412  if(gelside.Dimension() != 2) {
413  continue;
414  }
415  if(gelside.Neighbour() != gelside) {
416  continue;
417  }
418  TPZManVector<REAL,2> xi(2,0.);
419  gelside.CenterPoint(xi);
420  TPZFNMatrix<6,REAL> axes(2,3);
421  TPZFNMatrix<4,REAL> jac(2,2),jacinv(2,2);
422  REAL detjac;
423  gelside.Jacobian(xi, jac, axes, detjac, jacinv);
424  TPZManVector<REAL,3> x(3,0.);
425  gelside.X(xi, x);
426  TPZManVector<REAL,3> normal(3);
427  normal[0] = fabs(axes(0,1)*axes(1,2)-axes(0,2)*axes(1,1));
428  normal[1] = fabs(-axes(0,0)*axes(1,2)+axes(0,2)*axes(1,0));
429  normal[2] = fabs(axes(0,0)*axes(1,1)-axes(0,1)*axes(1,0));
430  REAL tol = 1.e-6;
431  REAL xmin = 1., xmax = 0.;
432  int numtol = 0;
433  for(int i=0; i<3; i++) {
434  if(xmin > x[i]) xmin = x[i];
435  if(xmax < x[i]) {
436  xmax = x[i];
437  }
438  if(normal[i] > tol) {
439  numtol++;
440  }
441  }
442  if(numtol != 1) {
443  DebugStop();
444  }
445  if(xmin > tol && xmax < 1.-tol) {
446  DebugStop();
447  }
448  }
449  }
450 }
451 
452 
454 {
455  int64_t nel = gmesh->NElements();
456  for(int64_t el=0; el<nel; el++) {
457  TPZGeoEl *gel = gmesh->ElementVec()[el];
458  int nsides = gel->NSides();
459  for(int is=0; is<nsides; is++) {
460  TPZGeoElSide gelside(gel,is);
461  if(gelside.Dimension() != 2) {
462  continue;
463  }
464  if(gelside.Neighbour() != gelside) {
465  continue;
466  }
467  TPZManVector<REAL,2> xi(2,0.);
468  //gelside.CenterPoint(xi);
469  TPZFNMatrix<6,REAL> axes(2,3);
470  TPZFNMatrix<4,REAL> jac(2,2),jacinv(2,2);
471  REAL detjac;
472  gelside.Jacobian(xi, jac, axes, detjac, jacinv);
473  TPZManVector<REAL,3> x(3,0.);
474  //gelside.X(xi, x);
475  TPZManVector<REAL,3> normal(3);
476  normal[0] = axes(0,1)*axes(1,2)-axes(0,2)*axes(1,1);
477  normal[1] = -axes(0,0)*axes(1,2)+axes(0,2)*axes(1,0);
478  normal[2] = axes(0,0)*axes(1,1)-axes(0,1)*axes(1,0);
479  REAL tol = 1.e-6;
480  REAL xmin = 1., xmax = 0.;
481  int numfound = 0;
482  int dir = -5;
483  for (int i=0; i<3; i++) {
484  if (fabs(normal[i] - 1.) < tol) {
485  dir = 1+i;
486  numfound++;
487  }
488  if (fabs(normal[i] + 1.) < tol) {
489  dir = -1-i;
490  numfound++;
491  }
492  }
493  if (dir == -5) {
494  DebugStop();
495  }
496  if(numfound != 1) {
497  DebugStop();
498  }
499  switch (dir) {
500  case -3:
501  TPZGeoElBC(gelside,fBCNumbers[0]);
502  break;
503  case 1:
504  TPZGeoElBC(gelside,fBCNumbers[1]);
505  break;
506  case 2:
507  TPZGeoElBC(gelside,fBCNumbers[2]);
508  break;
509  case -1:
510  TPZGeoElBC(gelside,fBCNumbers[3]);
511  break;
512  case -2:
513  TPZGeoElBC(gelside,fBCNumbers[4]);
514  break;
515  case 3:
516  TPZGeoElBC(gelside,fBCNumbers[5]);
517  break;
518  default:
519  DebugStop();
520  break;
521  }
522  }
523  }
524  return 0;
525 }
526 
528 {
529  int64_t nel = gmesh->NElements();
530  for(int64_t el=0; el<nel; el++) {
531  TPZGeoEl *gel = gmesh->ElementVec()[el];
532  int nsides = gel->NSides();
533  for(int is=0; is<nsides; is++) {
534  TPZGeoElSide gelside(gel,is);
535  if(gelside.Dimension() != 2) {
536  continue;
537  }
538  if(gelside.Neighbour() != gelside) {
539  continue;
540  }
541  TPZManVector<REAL,2> xi(2,0.);
542  TPZFNMatrix<6,REAL> axes(2,3);
543  TPZFNMatrix<4,REAL> jac(2,2),jacinv(2,2);
544  REAL detjac;
545  gelside.Jacobian(xi, jac, axes, detjac, jacinv);
546  TPZManVector<REAL,3> x(3,0.);
547 
548  REAL tol = 1.e-6;
549  REAL coordmin = 0., coordmax = 0.;
550  const int nnodes = gelside.NSideNodes();
551  int nxbot = 0, nybot = 0, nzbot = 0, nxtop = 0, nytop = 0, nztop = 0;
552  for (int i = 0; i < nnodes; i++) {
553  const int nodeindex = gelside.SideNodeIndex(i);
554  TPZManVector<REAL,3> coord(3,0);
555  gmesh->NodeVec()[nodeindex].GetCoordinates(coord);
556  const REAL diffxbot = fabs(coord[0]);
557  const REAL diffybot = fabs(coord[1]);
558  const REAL diffzbot = fabs(coord[2]);
559  const REAL diffxtop = fabs(coord[0]-1.);
560  const REAL diffytop = fabs(coord[1]-1.);
561  const REAL diffztop = fabs(coord[2]-1.);
562  if(diffxbot == 0){
563  nxbot++;
564  }
565  if(diffybot == 0){
566  nybot++;
567  }
568  if(diffzbot == 0){
569  nzbot++;
570  }
571  if(diffxtop == 0){
572  nxtop++;
573  }
574  if(diffytop == 0){
575  nytop++;
576  }
577  if(diffztop == 0){
578  nztop++;
579  }
580  }
581 
582  int numbcsfound = 0;
583  if (nxbot >= 3) {
584  TPZGeoElBC(gelside,fBCNumbers[0]);
585  numbcsfound++;
586  }
587  if (nybot >= 3) {
588  TPZGeoElBC(gelside,fBCNumbers[1]);
589  numbcsfound++;
590  }
591  if (nzbot >= 3) {
592  TPZGeoElBC(gelside,fBCNumbers[2]);
593  numbcsfound++;
594  }
595  if (nxtop >= 3) {
596  TPZGeoElBC(gelside,fBCNumbers[3]);
597  numbcsfound++;
598  }
599  if (nytop >= 3) {
600  TPZGeoElBC(gelside,fBCNumbers[4]);
601  numbcsfound++;
602  }
603  if (nztop >= 3) {
604  TPZGeoElBC(gelside,fBCNumbers[5]);
605  numbcsfound++;
606  }
607 
608  if (numbcsfound != 1) {
609  DebugStop();
610  }
611  }
612  }
613  return 0;
614 }
615 
616 
617 
618 
int NSideNodes() const
TPZAcademicGeoMesh()
sets up parameters to create de hexahedral regular mesh
expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ expr_ fabs
Definition: tfadfunc.h:140
Contains TPZAnalysis class which implements the sequence of actions to perform a finite element analy...
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.
Implements a vector class which allows to use external storage provided by the user. Utility.
Definition: pzquad.h:16
static int pyramid[2][5]
TPZGeoMesh * RedBlackPyramidalAndHexagonalMesh()
virtual TPZGeoEl * CreateGeoElement(MElementType type, TPZVec< int64_t > &cornerindexes, int matid, int64_t &index, int reftype=1)
Generic method for creating a geometric element. Putting this method centrally facilitates the modifi...
Definition: pzgmesh.cpp:1296
virtual void resize(const int64_t newsize)
Definition: pzvec.h:213
int64_t NElements() const
Number of elements of the mesh.
Definition: pzgmesh.h:129
Utility class which represents an element with its side. The Geometric approximation classes Geometry...
Definition: pzgeoelside.h:83
Contains declaration of TPZGeoElBC class, it is a structure to help the construction of geometric ele...
int64_t NElements() const
Access method to query the number of elements of the vector.
TPZGeoElSide Neighbour() const
Definition: pzgeoel.h:754
virtual int NSides() const =0
Returns the number of connectivities of the element.
void CenterPoint(TPZVec< REAL > &center) const
return the coordinates of the center in master element space (associated with the side) ...
void SetDimension(int dim)
Set Dimension.
Definition: pzgmesh.h:285
Contains the TPZMatPoisson3d class.
TPZManVector< int, 6 > fBCNumbers
indices of the boundary conditions
Contains TPZShapeLinear class which implements the shape functions of a linear one-dimensional elemen...
bool fShouldDeform
whether the mesh should be deformed or not
TPZGeoMesh * HexahedralMesh()
static int tetraedra[2][4]
Contains the TPZBndCond class which implements a boundary condition for TPZMaterial objects...
void Resize(const int newsize)
Increase the size of the chunk vector.
Definition: pzadmchunk.h:280
static const double tol
Definition: pzgeoprism.cpp:23
expr_ dx(i) *cos(expr_.val())
Defines the behaviour of all geometric elements. GeometryTPZGeoEl is the common denominator for all g...
Definition: pzgeoel.h:43
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
int64_t NNodes() const
Number of nodes of the mesh.
Definition: pzgmesh.h:126
#define LOGPZ_DEBUG(A, B)
Define log for debug info.
Definition: pzlog.h:87
void X(TPZVec< REAL > &loc, TPZVec< REAL > &result) const
X coordinate of a point loc of the side.
TPZAdmChunkVector< TPZGeoNode > & NodeVec()
Definition: pzgmesh.h:140
int fMaterialId
material index for volume elements
Contains the TPZGenGrid class which implements the generation of a multilayered geometric grid (two-d...
Contains declaration of TPZCheckMesh class which verifies the consistency of the datastructure of a T...
Contains the TPZExtendGridDimension class which generates a three dimensional mesh as an extension of...
int64_t fNumberElements
number of elements in any direction
TPZGeoMesh * PyramidalAndTetrahedralMesh()
static int prism[2][6]
Structure to help the construction of geometric elements along side of a given geometric element...
Definition: pzgeoelbc.h:21
void BuildConnectivity()
Build the connectivity of the grid.
Definition: pzgmesh.cpp:967
int AddBoundaryElementsByCoord(TPZGeoMesh *gmesh)
This class implements a geometric mesh for the pz environment. Geometry.
Definition: pzgmesh.h:48
int Dimension() const
the dimension associated with the element/side
void Jacobian(TPZVec< REAL > &param, TPZFMatrix< REAL > &jacobian, TPZFMatrix< REAL > &axes, REAL &detjac, TPZFMatrix< REAL > &jacinv) const
Jacobian associated with the side of the element.
Definition: pzeltype.h:61
static int tetraedra_2[6][4]
Contains TPZStepSolver class which defines step solvers class.
void GenerateNodes(TPZGeoMesh *gmesh)
put the geometric nodes in the geometric mesh
void DeformGMesh(TPZGeoMesh &gmesh)
Deform the geometric mesh according to the coordinates of fDeformed.
TPZGeoMesh * TetrahedralMesh()
int AddBoundaryElements(TPZGeoMesh *gmesh)
void CheckConsistency(TPZGeoMesh *gmesh)
verify if the faces without neighbour should be orthogonal to the main planes
Contains the TPZVTKGeoMesh class which implements the graphical mesh to VTK environment to geometric ...
void push_back(const T object)
Definition: pzstack.h:48
int64_t SideNodeIndex(int nodenum) const
Returns the index of the nodenum node of side.
Contains the TPZRefPatternTools class which defines tools of pattern.
Non abstract class which implements full matrices with preallocated storage with (N+1) entries...
Definition: pzfmatrix.h:716
TPZAdmChunkVector< TPZGeoEl * > & ElementVec()
Methods for handling pzlists.
Definition: pzgmesh.h:138