ParMooN
 All Classes Functions Variables Friends Pages
import3Dmesh3.h
1 
2 void ExampleFile()
3 {
4 #define __ROBINBC__
5  OutPut("Example: import3Dmesh.h" << endl);
6 }
7 
8 void GridBoundCondition(double x, double y, double z, BoundCond &cond)
9 {
10  cond = DIRICHLET;
11 }
12 
13 void SolidBoundCondition(double x, double y, double z, BoundCond &cond)
14 {
15  cond = DIRICHLET;
16 }
17 
18 void FluidBoundCondition(double x, double y, double z, BoundCond &cond)
19 {
20  cond = DIRICHLET;
21 }
22 // exact solution
23 void Exact(double x, double y, double z, double *values)
24 {
25  values[0] = 0;
26  values[1] = 0;
27  values[2] = 0;
28  values[3] = 0;
29  values[4] = 0;
30 }
31 
32 // initial conditon
33 void InitialCondition(double x, double y, double z, double *values)
34 {
35  values[0] = TDatabase::ParamDB->P0;
36  values[0] = x;
37 }
38 
39 // kind of boundary condition
40 void BoundCondition(double x, double y, double z, BoundCond &cond)
41 {
42 
43 
44  if(TDatabase::ParamDB->P13==0)
45  { //Lside
46  if(x<225 && y>76 )
47  { cond = DIRICHLET; }
48  else
49  { cond = NEUMANN; }
50  }
51  else
52  {
53  if(x< 240 && y<-72. && z>191)
54  { cond = DIRICHLET; }
55  else
56  { cond = NEUMANN; }
57  }
58 }
59 
60 // value of boundary condition
61 void BoundValue(double x, double y, double z, double &value)
62 {
63  value = 0;
64 }
65 
66 void BilinearCoeffs(int n_points, double *X, double *Y, double *Z, double **parameters, double **coeffs)
67 {
68  double eps=0;
69  int i, Region_ID;
70  double *coeff;
71  double x, y, z, hK;
72 
73 
74  double rhsfact, alpha, char_L, beam_r, DimlessBeam_r;
75  double xp=0., yp=0., zp=0., SourceCoord, Sourceradius;
76 
77  if(TDatabase::ParamDB->P13==0)
78  { // center point is (291.977 -1.8733333333333333 159.89 ), old, wrong
79  // 2.922130000000000e+02 -1.8733333333333333e+00 1.583800000000000e+02 1
80  xp=2.922130000000000e+02;
81  zp=1.583800000000000e+02;
82  }
83  else
84  { // side 300.565 69.66666666666666 157.038
85  xp= 300.565;
86  zp=157.038;
87  }
88 
89 
90  beam_r = TDatabase::ParamDB->P11;
91  char_L = TDatabase::ParamDB->P12;
92 
93  hK = coeffs[0][2];
94  Region_ID = coeffs[0][1];
95  DimlessBeam_r = beam_r/char_L;
96 
97 
98  // cout << "Region_ID = " << Region_ID << endl;
99 
100  if(Region_ID==0) //0-Fat/adipost
101  {
102  eps = TDatabase::ParamDB->P1/(1100.*4483*char_L*char_L); ;
103  alpha = TDatabase::ParamDB->P7;
104  rhsfact = (alpha*TDatabase::ParamDB->P4)/(1100.0*4483*(22./7.)*beam_r*beam_r);
105  }
106  else
107  {
108  eps = TDatabase::ParamDB->P2/(1040.0*3500.*char_L*char_L);
109  alpha = TDatabase::ParamDB->P8;
110  rhsfact = alpha*TDatabase::ParamDB->P4/(1040.0*3500.*(22./7.)*beam_r*beam_r);
111  }
112 
113 // //2-CSF
114 // eps = TDatabase::ParamDB->P2/(997.0*3710*char_L*char_L);
115 // alpha = TDatabase::ParamDB->P8;
116 // rhsfact = alpha*TDatabase::ParamDB->P4/(997.0*3710*(22./7.)*beam_r*beam_r);
117 //
118  for(i=0;i<n_points;i++)
119  {
120  coeff = coeffs[i];
121 
122  x = X[i];
123  y = Y[i];
124  z = Z[i];
125 
126  // diffusion
127  coeff[0] = eps;
128  // convection in x direction
129  coeff[1] = 0;
130  // convection in y direction
131  coeff[2] = 0;
132  // convection in z direction
133  coeff[3] = 0;
134  // reaction
135  coeff[4] = 0;
136  // rhs
137 
138  if(TDatabase::ParamDB->P13==0)
139  { //Lside
140  Sourceradius = sqrt( (x-xp)*(x-xp) + (z-zp)*(z-zp));
141 
142  if(y<0.) y=0.;
143  SourceCoord = -y;
144  }
145  else
146  { // side
147  Sourceradius = sqrt( (x-xp)*(x-xp) + (z-zp)*(z-zp));
148 
149  if(y>0.) y=0.;
150  SourceCoord = y;
151  }
152 
153 
154  if(Sourceradius<=DimlessBeam_r)
155  {
156  coeff[5] = rhsfact*exp(alpha*SourceCoord*char_L);// f
157 
158 // cout << "Sourceradius : " << Sourceradius << endl;
159  if(TDatabase::ParamDB->P14< hK) TDatabase::ParamDB->P14=hK;
160  }
161  else
162  {coeff[5] = 0; }// f
163 
164  coeff[6] = 0;
165 
166  }
167 }
168 
169 void ReadAdditionalNModes(char *NODE, tetgenio &InAddPts)
170 {
171  int i, N_Vertices, index, N_Atribude, BDMarker;
172  char line[100];
173  std::ifstream dat(NODE);
174 
175  if (!dat)
176  {
177  cerr << "cannot open '" << NODE << "' for input" << endl;
178  exit(-1);
179  }
180 
181  dat.getline (line, 99); // first line is command line
182  dat >> N_Vertices >> N_Atribude >> BDMarker;
183  dat.getline (line, 99);
184 
185  InAddPts.numberofpoints = N_Vertices;
186  InAddPts.pointlist = new double[3*N_Vertices];
187 
188  //read from file
189  for(i=0;i<N_Vertices; i++)
190  {
191  dat.getline (line, 99);
192  dat >> index >> InAddPts.pointlist[3*i] >> InAddPts.pointlist[3*i+1] >> InAddPts.pointlist[3*i+2];
193  cout<< i << " vert X: " <<InAddPts.pointlist[3*i] << " vert Y: " <<InAddPts.pointlist[3*i+1] <<" vert Z: " <<InAddPts.pointlist[3*i+2] <<endl;
194 
195  }
196 
197  dat.close();
198 } // ReadAdditionalNModes
199 
200 void ReadMeditMeshWithCells(char *SMESH, tetgenio &In)
201 {
202  int i, j, k, dimension, N_FVert, N_Faces, N_Vertices, N_CVert, N_Cells;
203  int BDComp_Min=10000000;
204  int *plist;
205 
206  double x;
207 
208  char line[100];
209  tetgenio::facet *F;
210  tetgenio::polygon *P;
211 
212  std::ifstream dat(SMESH);
213 
214  if (!dat)
215  {
216  cerr << "cannot open '" << SMESH << "' for input" << endl;
217  exit(-1);
218  }
219 
220  // check the dimension
221  while (!dat.eof())
222  {
223  dat >> line;
224 
225  if ( (!strcmp(line, "Dimension")) || (!strcmp(line, "dimension")) || (!strcmp(line, "DIMENSION")))
226  {
227  dat.getline (line, 99);
228  dat >> dimension;
229  break;
230  }
231  // read until end of line
232  dat.getline (line, 99);
233  }
234 
235  if(dimension!=3)
236  {
237  cerr << "dimension: " << dimension << endl;
238  exit(-1);
239  }
240 
241  // find the N_Vertices
242  while (!dat.eof())
243  {
244  dat >> line;
245 
246  if ( (!strcmp(line, "Vertices")) || (!strcmp(line, "vertices")) || (!strcmp(line, "VERTICES")) )
247  {
248  dat.getline (line, 99);
249  dat >> N_Vertices;
250  break;
251  }
252  // read until end of line
253  dat.getline (line, 99);
254  }
255 
256  In.numberofpoints = N_Vertices;
257  In.pointlist = new double[3*N_Vertices];
258 
259  //read from file
260  for(i=0;i<N_Vertices; i++)
261  {
262  dat.getline (line, 99);
263  dat >> x >> In.pointlist[3*i+1] >> In.pointlist[3*i+2];
264 
265  In.pointlist[3*i] = x - 3.606596999999999830777142; //center point
266 
267 // if(i==1063 || i==1061 || i==1062 || i==162 || i==1175 || i==517 || i==1064 || i==518 ) // i==1063 as center
268  if(i==1062 || i==916 || i==341 || i==914 || i==162 || i==1063 || i==1061) // i==1062 as center
269  cout<< i << " vert X: " <<In.pointlist[3*i] << " vert Y: " <<In.pointlist[3*i+1] <<" vert Z: " <<In.pointlist[3*i+2] <<endl;
270  }
271 
272  // find the N_Triangles
273  // find the N_Vertices
274  while (!dat.eof())
275  {
276  dat >> line;
277 
278  if ( (!strcmp(line, "Triangles")) || (!strcmp(line, "triangles")) || (!strcmp(line, "TRIANGLES")) )
279  {
280  N_FVert = 3;
281  dat.getline (line, 99);
282  dat >> N_Faces;
283  break;
284  }
285  else if ( (!strcmp(line, "Quadrilaterals")) || (!strcmp(line, "quadrilaterals")) || (!strcmp(line, "QUADRILATERALS")) )
286  {
287  N_FVert = 4;
288  dat.getline (line, 99);
289  dat >> N_Faces;
290  break;
291  }
292 
293  // read until end of line
294  dat.getline (line, 99);
295  }
296 
297  In.numberoffacets = N_Faces;
298  In.facetlist = new tetgenio::facet[In.numberoffacets];
299  In.facetmarkerlist = new int[In.numberoffacets];
300 
301  for(i=0;i<N_Faces; i++)
302  {
303  dat.getline (line, 99);
304 
305  F = &In.facetlist[i];
306  tetgenio::init(F);
307  F->numberofpolygons = 1;
308  F->polygonlist = new tetgenio::polygon[F->numberofpolygons];
309  F->numberofholes = 0;
310  F->holelist = NULL;
311  P = &F->polygonlist[0];
312  tetgenio::init(P);
313  P->numberofvertices = N_FVert;
314  P->vertexlist = new int[P->numberofvertices];
315 
316  for(j=0;j<N_FVert;j++)
317  {
318  dat >> k;
319  P->vertexlist[j] = k-1; // c numbering
320  }
321 
322  dat >> In.facetmarkerlist[i];
323 
324  if(BDComp_Min > In.facetmarkerlist[i])
325  BDComp_Min=In.facetmarkerlist[i];
326 
327  // cout << i<< " P->vertexlist[j]: " << In.facetmarkerlist[i] << endl;
328  } // for(i=0;i<N_Faces; i++)
329 
330  BDComp_Min--;
331 
332  for(i=0;i<N_Faces; i++)
333  In.facetmarkerlist[i] -= BDComp_Min;
334 
335 
337  while (!dat.eof())
338  {
339  dat >> line;
340 
341  if ( (!strcmp(line, "Tetrahedron")) || (!strcmp(line, "tetrahedron")) || (!strcmp(line, "TETRAHEDRON")) )
342  {
343  N_CVert = 4;
344  dat.getline (line, 99);
345  dat >> N_Cells;
346  break;
347  }
348  // read until end of line
349  dat.getline (line, 99);
350  }
351 
352  In.numberoftetrahedra = N_Cells;
353  In.numberofcorners = N_CVert;
354  In.numberoftetrahedronattributes = 1;
355  In.tetrahedronlist = new int[N_Cells * 4];
356  In.tetrahedronattributelist = new REAL[N_Cells];
357  In.tetrahedronvolumelist = new REAL[N_Cells];
358 
359  for(i=0; i<N_Cells; i++)
360  {
361  dat.getline (line, 99);
362  plist = &(In.tetrahedronlist[i * 4]);
363 
364  for(j=0;j<N_CVert;j++)
365  {
366  dat >> k;
367  plist[j] = k -1;
368  }// for(j=0;j<N_CV
369 
370 
371  dat >> In.tetrahedronattributelist[i];
372 // In.tetrahedronvolumelist [i] = 100.;
373  } // for(i=0; i<N_
374 
375 
376 
377 // cout << i<< " N_Cells : " << N_Cells << endl;
378 // exit(0);
379 
380  dat.close();
381 // exit(0);
382 } // ReadMeditMeshWithCells
383 
384 void ReadMeditMesh(char *SMESH, tetgenio &In)
385 {
386  int i, j, k, dimension, N_FVert, N_Faces, N_Vertices;
387  int BDComp_Min=10000000;
388 
389  char line[100];
390  tetgenio::facet *F;
391  tetgenio::polygon *P;
392 
393  std::ifstream dat(SMESH);
394 
395  if (!dat)
396  {
397  cerr << "cannot open '" << SMESH << "' for input" << endl;
398  exit(-1);
399  }
400 
401  // check the dimension
402  while (!dat.eof())
403  {
404  dat >> line;
405 
406  if ( (!strcmp(line, "Dimension")) || (!strcmp(line, "dimension")) || (!strcmp(line, "DIMENSION")))
407  {
408  dat.getline (line, 99);
409  dat >> dimension;
410  break;
411  }
412  // read until end of line
413  dat.getline (line, 99);
414  }
415 
416  if(dimension!=3)
417  {
418  cerr << "dimension: " << dimension << endl;
419  exit(-1);
420  }
421 
422  // find the N_Vertices
423  while (!dat.eof())
424  {
425  dat >> line;
426 
427  if ( (!strcmp(line, "Vertices")) || (!strcmp(line, "vertices")) || (!strcmp(line, "VERTICES")) )
428  {
429  dat.getline (line, 99);
430  dat >> N_Vertices;
431  break;
432  }
433  // read until end of line
434  dat.getline (line, 99);
435  }
436 
437  In.numberofpoints = N_Vertices;
438  In.pointlist = new double[3*N_Vertices];
439 
440  //read from file
441  for(i=0;i<N_Vertices; i++)
442  {
443  dat.getline (line, 99);
444  dat >> In.pointlist[3*i] >> In.pointlist[3*i+1] >> In.pointlist[3*i+2];
445  //cout<< i << " vert X: " <<In.pointlist[3*i] << " vert Y: " <<In.pointlist[3*i+1] <<endl;
446  }
447 
448  // find the N_Triangles
449  // find the N_Vertices
450  while (!dat.eof())
451  {
452  dat >> line;
453 
454  if ( (!strcmp(line, "Triangles")) || (!strcmp(line, "triangles")) || (!strcmp(line, "TRIANGLES")) )
455  {
456  N_FVert = 3;
457  dat.getline (line, 99);
458  dat >> N_Faces;
459  break;
460  }
461  else if ( (!strcmp(line, "Quadrilaterals")) || (!strcmp(line, "quadrilaterals")) || (!strcmp(line, "QUADRILATERALS")) )
462  {
463  N_FVert = 4;
464  dat.getline (line, 99);
465  dat >> N_Faces;
466  break;
467  }
468 
469  // read until end of line
470  dat.getline (line, 99);
471  }
472 
473  In.numberoffacets = N_Faces;
474  In.facetlist = new tetgenio::facet[In.numberoffacets];
475  In.facetmarkerlist = new int[In.numberoffacets];
476 
477  for(i=0;i<N_Faces; i++)
478  {
479  dat.getline (line, 99);
480 
481  F = &In.facetlist[i];
482  tetgenio::init(F);
483  F->numberofpolygons = 1;
484  F->polygonlist = new tetgenio::polygon[F->numberofpolygons];
485  F->numberofholes = 0;
486  F->holelist = NULL;
487  P = &F->polygonlist[0];
488  tetgenio::init(P);
489  P->numberofvertices = N_FVert;
490  P->vertexlist = new int[P->numberofvertices];
491 
492  for(j=0;j<N_FVert;j++)
493  {
494  dat >> k;
495  P->vertexlist[j] = k-1; // c numbering
496  }
497 
498  dat >> In.facetmarkerlist[i];
499 
500  if(BDComp_Min > In.facetmarkerlist[i])
501  BDComp_Min=In.facetmarkerlist[i];
502 
503  // cout << i<< " P->vertexlist[j]: " << In.facetmarkerlist[i] << endl;
504  } // for(i=0;i<N_Faces; i++)
505 
506  BDComp_Min--;
507 
508  for(i=0;i<N_Faces; i++)
509  In.facetmarkerlist[i] -= BDComp_Min;
510 
511 // cout << i<< " P->vertexlist[j]: " << BDComp_Min << endl;
512 
513 
514  dat.close();
515 // exit(0);
516 } // ReadMeditMesh
517 
518 void DeleteDomain(TDomain *&Domain)
519 {
520  int i, j, k, N_Cells, N_RootCells;
521  int CurrVertex, N_Faces, N_Vertices, ID;
522  TBaseCell **CellTree, **SurfCellTree, *cell;
523  TGridCell **DelCell;
524  TVertex **VertexDel;
525  TCollection *coll;
526 
527  Domain->GetTreeInfo(CellTree,N_RootCells);
528  coll = Domain->GetCollection(It_Finest, 0);
529  N_Cells = coll->GetN_Cells();
530 
531  VertexDel = new TVertex*[8*N_RootCells];
532  CurrVertex = 0;
533 
534  for(i=0;i<N_Cells;i++)
535  {
536  cell = coll->GetCell(i);
537  N_Faces = cell->GetN_Faces();
538  N_Vertices = cell->GetN_Vertices();
539  for(j=0;j<N_Faces;j++)
540  {
541  if(CurrVertex==0)
542  {
543  VertexDel[CurrVertex++] = cell->GetVertex(j);
544  }
545  else
546  {
547  ID = 0;
548  for(k=0;k<CurrVertex;k++)
549  if(VertexDel[k]==cell->GetVertex(j))
550  {
551  ID = 1; break;
552  }
553  if(ID!=1)
554  {
555  VertexDel[CurrVertex++] = cell->GetVertex(j);
556  }
557  } // else if(CurrVertex==0)
558 
559  ID = 0;
560  for(k=0;k<CurrVertex;k++)
561  if(VertexDel[k]==cell->GetVertex((j+1)%N_Vertices))
562  {
563  ID = 1; break;
564  }
565  if(ID!=1)
566  {
567  VertexDel[CurrVertex++] = cell->GetVertex((j+1)%N_Vertices);
568  }
569 
570  ID = 0;
571  for(k=0;k<CurrVertex;k++)
572  if(VertexDel[k]==cell->GetVertex((j+2)%N_Vertices))
573  {
574  ID = 1; break;
575  }
576  if(ID!=1)
577  {
578  VertexDel[CurrVertex++] = cell->GetVertex((j+2)%N_Vertices);
579  }
580  if(N_Faces==6) // If old cell is hexahedrol
581  {
582  ID = 0;
583  for(k=0;k<CurrVertex;k++)
584  if(VertexDel[k]==cell->GetVertex((j+4)%N_Vertices))
585  {
586  ID = 1; break;
587  }
588  if(ID!=1)
589  {
590  VertexDel[CurrVertex++] = cell->GetVertex((j+4)%N_Vertices);
591  }
592  }
593  } // for j
594  } // for i
595 
596  for(i=0;i<CurrVertex;i++)
597  delete VertexDel[i];
598  delete [] VertexDel;
599  OutPut(CurrVertex<<" vertices were deleted"<<endl);
600 
601  // remove all existing cells and joints
602  for(i=0;i<N_RootCells;i++)
603  delete (TGridCell*)CellTree[i];
604  delete [] CellTree;
605  OutPut(N_RootCells<<" cells were deleted"<<endl);
606 }
607 
608 void getboundarytag(int a, int b, int c, int NBdfaces, int *Bdfaces, int &CompID)
609 {
610  int i,a1,b1,c1,marker;
611  for(i=0;i<NBdfaces;i++)
612  {
613  a1 = Bdfaces[4*i];
614  b1 = Bdfaces[4*i+1];
615  c1 = Bdfaces[4*i+2];
616  marker = Bdfaces[4*i+3];
617  if((a==a1||a==b1||a==c1)&&(b==a1||b==b1||b==c1)&&(c==a1||c==b1||c==c1))
618  {
619  CompID = marker;
620  break;
621  }
622  }
623 }
624 
625 void TetrameshCreate(TDomain *&Domain)
626 {
627  int i, j, k, l, dimension, N_Vertices, NBdfaces, *Bdfaces, BoundaryMarker;
628  int N_Faces, N_RootCells;
629  int v1, v2, v3, v4, CellMarker, RefLevel=0, *Tetrahedrals, *PointNeighb, maxEpV, *Tetrahedrals_loc;
630  int a, b, c, Neib[2], Neighb_tmp, CurrNeib, len1, len2, len3, CurrComp=0, N_Points ;
631  int N_Cells, MaxLen, jj, CompID=0;
632  const int *EdgeVertex, *TmpFV, *TmpLen;
633  double X, Y, Z, DispX, DispY, DispZ;
634  double Xmin = 1e10, Xmax = -1e10, Ymin = 1e10, Ymax = -1e10;
635  double Zmin = 1e10, Zmax = -1e10, T[4]={0,0,0,0}, S[4]={0,0,0,0};
636  double StartX, StartY, StartZ, BoundX, BoundY, BoundZ;
637  char *MESH, line[100];
638  MESH = TDatabase::ParamDB->SMESHFILE;
639  TVertex **NewVertices;
640  TBaseCell **CellTree, *cell;
641  TJoint *Joint;
642  TBoundComp3D *bdface;
643  TShapeDesc *ShapeDesc;
644  TCollection *coll;
645  TBoundPart *someboundpart;
646  std::ifstream dat(MESH);
647 
648  if (!dat)
649  {
650  cerr << "cannot open '" << MESH << "' for input" << endl;
651  exit(-1);
652  }
653 
654  while (!dat.eof())
655  {
656  dat >> line;
657  if ( (!strcmp(line, "Dimension")) || (!strcmp(line, "dimension")) || (!strcmp(line, "DIMENSION")))
658  {
659  dat.getline (line, 99);
660  dat >> dimension;
661  break;
662  }
663  dat.getline (line, 99);
664  }
665 
666  if(dimension!=3)
667  {
668  cerr << "dimension: " << dimension << endl;
669  exit(-1);
670  }
671 
672  while (!dat.eof())
673  {
674  dat >> line;
675  if ( (!strcmp(line, "Vertices")) || (!strcmp(line, "vertices")) || (!strcmp(line, "VERTICES")) )
676  {
677  dat.getline (line, 99);
678  dat >> N_Vertices;
679  break;
680  }
681  dat.getline (line, 99);
682  }
683  cout <<"N_Vertices "<<N_Vertices<<endl;
684  NewVertices = new TVertex*[N_Vertices];
685 
686  for(i=0;i<N_Vertices; i++)
687  {
688  dat.getline (line, 99);
689  dat >> X >> Y >> Z;
690  NewVertices[i] = new TVertex(X, Y, Z);
691  if (X > Xmax) Xmax = X;
692  if (X < Xmin) Xmin = X;
693  if (Y > Ymax) Ymax = Y;
694  if (Y < Ymin) Ymin = Y;
695  if (Z > Zmax) Zmax = Z;
696  if (Z < Zmin) Zmin = Z;
697  }
698  // set bounding box
699  StartX = Xmin;
700  StartY = Ymin;
701  StartZ = Zmin;
702  BoundX = Xmax - Xmin;
703  BoundY = Ymax - Ymin;
704  BoundZ = Zmax - Zmin;
705 
706  Domain->SetBoundBox(StartX, StartY, StartZ, BoundX, BoundY, BoundZ);
707  //========================
708  while (!dat.eof())
709  {
710  dat >> line;
711  if ( (!strcmp(line, "Triangles")) || (!strcmp(line, "triangles")) || (!strcmp(line, "TRIANGLES")) )
712  {
713  dat.getline (line, 99);
714  dat >> NBdfaces;
715  break;
716  }
717  // read until end of line
718  dat.getline (line, 99);
719  }
720  cout<<"number of boundary faces "<<NBdfaces<<endl;
721  Bdfaces = new int[4*NBdfaces];
722 
723  for(i=0;i<NBdfaces;i++)
724  {
725  dat.getline (line, 99);
726  dat >> v1 >> v2 >> v3 >> BoundaryMarker;
727  Bdfaces[4*i ] = v1;
728  Bdfaces[4*i + 1] = v2;
729  Bdfaces[4*i + 2] = v3;
730  Bdfaces[4*i + 3] = BoundaryMarker;
731  }
732  //================================
733 
734 
735  while (!dat.eof())
736  {
737  dat >> line;
738 
739  if ( (!strcmp(line, "Tetrahedra")) || (!strcmp(line, "tetrahedra")) || (!strcmp(line, "TETRAHEDRA")) )
740  {
741  dat.getline (line, 99);
742  dat >> N_RootCells;
743  break;
744  }
745  // read until end of line
746  dat.getline (line, 99);
747  }
748  // generate new cells
749  cout << "number of root cells "<<N_RootCells <<endl;
750  CellTree = new TBaseCell*[N_RootCells];
751  Tetrahedrals = new int[4*N_RootCells];
752 
753  for (i=0;i<N_RootCells;i++)
754  {
755  dat.getline (line, 99);
756  dat >> v1 >> v2 >> v3 >> v4 >> CellMarker;
757  Tetrahedrals[4*i ] = v1-1;
758  Tetrahedrals[4*i + 1] = v2-1;
759  Tetrahedrals[4*i + 2] = v3-1;
760  Tetrahedrals[4*i + 3] = v4-1;
761 
762  CellTree[i] = new TMacroCell(TDatabase::RefDescDB[Tetrahedron], RefLevel);
763  CellTree[i]->SetRegionID(CellMarker);
764  CellTree[i]->SetAsLayerCell(1);
765  CellTree[i]->SetVertex(0, NewVertices[Tetrahedrals[4*i ]]);
766  CellTree[i]->SetVertex(1, NewVertices[Tetrahedrals[4*i + 1]]);
767  CellTree[i]->SetVertex(2, NewVertices[Tetrahedrals[4*i + 2]]);
768  CellTree[i]->SetVertex(3, NewVertices[Tetrahedrals[4*i + 3]]);
769  CellTree[i]->SetClipBoard(i);
770  ((TMacroCell *) CellTree[i])->SetSubGridID(0);
771  }
772  dat.close();
773 
774  Domain->SetTreeInfo(CellTree, N_RootCells);
775  // initialize iterators
776  TDatabase::IteratorDB[It_EQ]->SetParam(Domain);
777  TDatabase::IteratorDB[It_LE]->SetParam(Domain);
778  TDatabase::IteratorDB[It_Finest]->SetParam(Domain);
779  TDatabase::IteratorDB[It_Between]->SetParam(Domain);
780  TDatabase::IteratorDB[It_OCAF]->SetParam(Domain);
781 
782  // search neighbours
783  PointNeighb = new int[N_Vertices];
784  memset(PointNeighb, 0, N_Vertices*SizeOfInt);
785 
786  for (i=0;i<4*N_RootCells;i++)
787  PointNeighb[Tetrahedrals[i]]++;
788 
789  maxEpV = 0;
790  for (i=0;i<N_Vertices;i++)
791  if (PointNeighb[i] > maxEpV) maxEpV = PointNeighb[i];
792  delete [] PointNeighb;
793  cout<<"maximum edges per vertex "<< maxEpV<<endl;
794  PointNeighb = new int[++maxEpV * N_Vertices];
795  memset(PointNeighb, 0, maxEpV*N_Vertices*SizeOfInt);
796 
797  // every vertex contains "maxEpV" columns
798  // for every vertex at first colomn contains the number of cells containing this vertex
799  // at further columns we set the index of corresponding cells containing this vertex
800  // cout<<"maxEpV*N_Vertices "<<maxEpV*N_Vertices<<endl;
801  for(i=0;i<4*N_RootCells;i++)
802  {
803  j = Tetrahedrals[i]*maxEpV;
804  PointNeighb[j]++;
805  PointNeighb[j + PointNeighb[j]] = i / 4;
806  }
807 
808  coll=Domain->GetCollection(It_Finest, 0);
809  N_Cells = coll->GetN_Cells();
810 
811  for(i=0;i<N_Cells;i++)
812  {
813  cell = coll->GetCell(i);
814  ShapeDesc= cell->GetShapeDesc();
815  ShapeDesc->GetFaceVertex(TmpFV, TmpLen, MaxLen);
816  // which vertices are on each face(local indexing)
817  // number of vertices on each face
818  // maximum number of vertices per face
819  N_Faces = cell->GetN_Faces();
820  Tetrahedrals_loc = Tetrahedrals+4*i;
821  for(jj=0;jj<N_Faces;jj++)
822  if(cell->GetJoint(jj) == NULL)
823  {
824  N_Points = TmpLen[jj];
825  if(N_Points!=3)
826  {
827  cerr << "Only Tria faces are allowed!!! N_FVert: " << N_Points << endl;
828  exit(-1);
829  }
830  a = Tetrahedrals_loc[TmpFV[jj*MaxLen]];
831  b = Tetrahedrals_loc[TmpFV[jj*MaxLen + 1]];
832  c = Tetrahedrals_loc[TmpFV[jj*MaxLen + 2]];
833  Neib[0] = -1;
834  Neib[1] = -1;
835  CurrNeib = 0;
836  len1 = PointNeighb[a*maxEpV];
837  len2 = PointNeighb[b*maxEpV];
838  len3 = PointNeighb[c*maxEpV];
839  // find the index of the cells containing current face with node(vertex) indices a,b,c
840  for (j=1;j<=len1;j++)
841  {
842  Neighb_tmp = PointNeighb[a*maxEpV + j];
843  for (k=1;k<=len2;k++)
844  {
845  if(Neighb_tmp == PointNeighb[b*maxEpV + k])
846  {
847  for(l=1;l<=len3;l++)
848  if(Neighb_tmp == PointNeighb[c*maxEpV + l])
849  {
850  Neib[CurrNeib++] = Neighb_tmp;
851  break;
852  }
853  }
854  }
855  if (CurrNeib == 2) break; //break from j<=len1 loop
856  } // for (j=1;j<=len1;j++)
857  if(CurrNeib==1)
858  {
859  getboundarytag(a,b,c, NBdfaces, Bdfaces, CompID);
860 // cout << "compid"<<CompID<<endl;
861  bdface = new TBdPlane(CompID);
862  if(bdface->GetTSofXYZ(NewVertices[a]->GetX(), NewVertices[a]->GetY(), NewVertices[a]->GetZ(), T[1], S[1])
863  ||
864  bdface->GetTSofXYZ(NewVertices[b]->GetX(), NewVertices[b]->GetY(), NewVertices[b]->GetZ(), T[2], S[2])
865  ||
866  bdface->GetTSofXYZ(NewVertices[c]->GetX(), NewVertices[c]->GetY(), NewVertices[c]->GetZ(), T[3], S[3]) )
867  {
868  cout <<"chek 11.5 "<< i<<endl;
869  cerr<<"Error: could not set parameter values"<<endl;
870  OutPut(NewVertices[a]<<endl);
871  OutPut(NewVertices[b]<<endl);
872  OutPut(NewVertices[c]<<endl);
873  exit(0);
874  }
875  Joint = new TBoundFace(bdface, T, S);
876  }
877  else
878  {
879  if (CurrNeib != 2)
880  cerr << "Error !!!!!!!! not enough neighbours!" << endl;
881  Joint = new TJointEqN(CellTree[Neib[0]], CellTree[Neib[1]]);
882  }
883  // First element containing the current face
884  // find the local index for the point 'a' on the cell
885  for (j=0;j<4;j++)
886  if (Tetrahedrals[4*Neib[0]+j] == a) break;
887  // find the local index for the point 'b' on the cell
888  for (k=0;k<4;k++)
889  if (Tetrahedrals[4*Neib[0]+k] == b) break;
890  // find the local index for the point 'c' on the cell
891  for (l=0;l<4;l++)
892  if (Tetrahedrals[4*Neib[0]+l] == c) break;
893  l = l*100 + k*10 + j;
894 
895  switch (l) // j will contain the local index for the current face
896  {
897  case 210: case 21: case 102:
898  case 120: case 12: case 201:
899  j = 0;
900  break;
901  case 310: case 31: case 103:
902  case 130: case 13: case 301:
903  j = 1;
904  break;
905  case 321: case 132: case 213:
906  case 231: case 123: case 312:
907  j = 2;
908  break;
909  case 230: case 23: case 302:
910  case 320: case 32: case 203:
911  j = 3;
912  break;
913 
914  default:
915  Error("Unable to set the face !!!!!!!!!!!!" << endl);
916  exit(0);
917  } // switch (l)
918  CellTree[Neib[0]]->SetJoint(j, Joint);
919 
920  // second cell containing the current face
921  if (Neib[1] != -1) // second element containing the current face
922  {
923  // find the local index for the point 'a' on the cell
924  for (j=0;j<4;j++)
925  if (Tetrahedrals[4*Neib[1]+j] == a) break;
926 
927  // find the local index for the point 'b' on the cell
928  for (k=0;k<4;k++)
929  if (Tetrahedrals[4*Neib[1]+k] == b) break;
930 
931  // find the local index for the point 'c' on the cell
932  for (l=0;l<4;l++)
933  if (Tetrahedrals[4*Neib[1]+l] == c) break;
934 
935  l = l*100 + k*10 + j;
936  switch (l) // j will contain the local index for the current face
937  {
938  case 210: case 21: case 102:
939  case 120: case 12: case 201:
940  j = 0;
941  break;
942  case 310: case 31: case 103:
943  case 130: case 13: case 301:
944  j = 1;
945  break;
946  case 321: case 132: case 213:
947  case 231: case 123: case 312:
948  j = 2;
949  break;
950  case 230: case 23: case 302:
951  case 320: case 32: case 203:
952  j = 3;
953  break;
954 
955  default:
956  Error("Unable to set the face !!!!!!!!!!!!" << endl);
957  exit(0);
958  }
959  CellTree[Neib[1]]->SetJoint(j, Joint);
960  } // if (Neib[1] != -1)
961 
962  if (Joint->GetType() == InterfaceJoint3D || Joint->GetType() == IsoInterfaceJoint3D)
963  {
964  ((TInterfaceJoint3D*)Joint)->SetMapType();
965  ((TInterfaceJoint3D*)(Joint))->CheckOrientation();
966  }
967  else
968  if (Joint->GetType() == JointEqN)
969  Joint->SetMapType();
970  } // if(cell->GetJoint(j)
971  } // for(i=0;i<N_Cells;i++)
972 delete [] PointNeighb;
973  cout<< "tetrameshcreatesuccessful " <<"\n";
974 }
975 
976 void TetrameshGen(TDomain *&Domain)
977 {
978  //======================================================================
979  // Tetgen for grid generation begin
980  //======================================================================
981  int i, j, k, l, N_Coord, *N_FVerts, N_Faces, *Facets;
982  int N, N_RootCells, N_Cells, CurrVertex, N_Vertices, ID, N_G, RefLevel=0;
983  int CurrNeib, len1, len2, len3, maxEpV = 0, a, b, c, Neib[2], Neighb_tmp, CurrComp;
984  int *Tetrahedrals, *PointNeighb, dimension;
985  int *Facelist, *Facemarkerlist;
986 
987  double *Coordinates, N_x, N_y, N_z;
988  double *Vertices;
989  double Xmin = 1e10, Xmax = -1e10, Ymin = 1e10, Ymax = -1e10;
990  double Zmin = 1e10, Zmax = -1e10, T[4]={0,0,0,0}, S[4]={0,0,0,0};
991  double StartX, StartY, StartZ, BoundX, BoundY, BoundZ;
992 
993  tetgenio In, Out;
994 
995  TBaseCell **CellTree, **SurfCellTree;
996  TGridCell **DelCell;
997  TVertex **VertexDel, **NewVertices, **NewSurfVertices;
998  TBoundPart *BoundPart;
999  TBdPlane **UpdateFaceParams;
1000  TJoint *Joint;
1001  TBoundComp3D *bdcomp;
1002  TCollection *coll, *SurfColl;
1003  TBaseCell *cell;
1004  TBdSphere *UpdateParam;
1005 
1006  char *SMESH, line[100];
1007  SMESH = TDatabase::ParamDB->SMESHFILE;
1008 
1009 #ifdef _MPI
1010  int rank, size;
1011  MPI_Comm Comm = TDatabase::ParamDB->Comm;
1012 
1013  MPI_Comm_rank(Comm, &rank);
1014  MPI_Comm_size(Comm, &size);
1015 
1016 #endif
1017 
1019 #ifdef _MPI
1020  if(rank==0)
1021 #endif
1022  {
1023 
1024 // In.initialize();
1025 // Out.initialize();
1026 
1027  std::ostringstream opts;
1028  opts << " ";
1029 
1030  opts.seekp(std::ios::beg);
1031  opts<<'p'; // Tetrahedralize the PLC. Switches are chosen to read a PLC (p)
1032 // opts<<'r'; // -r Reconstructs a previously generated mesh
1033 // opts<<"q"<<1.49; // quality mesh generation(q) with a specified quality bound
1034 // opts<<"a"<<0.1; // maximum volume constraint
1035 // opts<<'i'; // Inserts a list of additional points into mesh.
1036  opts<<'z'; // numbers all output items starting from zero
1037 // opts<<'d'; // Detect intersections of PLC facets.
1038  opts<<'f'; // Outputs all faces (including non-boundary)
1039  opts<<'e'; // Outputs a list of edges of the triangulation
1040 // opts<<'I'; // Suppresses mesh iteration numbers.
1041  opts<<'C'; // Checks the consistency of the final mesh.
1042 // opts<<'Q'; // Quiet: No terminal output except errors.
1043 // opts<<'g'; // Outputs mesh to .mesh file for viewing by Medit
1044  opts<<'Y'; // Suppresses boundary facets/segments splitting
1045 // opts<<'V'; //verbose mode
1046  opts<<ends;
1047 
1048 // cout << " SMESHFILE is " << SMESH << endl;
1050  ReadMeditMesh(SMESH, In);
1051 
1052 // for(i=0;i<In.numberofpoints;i++)
1053 // OutPut(i<<" (x, y, z) = "<<
1054 // In.pointlist[3*i]<<' '<<In.pointlist[3*i+1]<<' '<<In.pointlist[3*i+2]<<endl);
1055 
1056  // Calling tetrahedralize function of 3dtetgen mesh generator
1057  tetrahedralize((char*)opts.str().c_str(), &In, &Out);
1058 
1059  } // if(rank==0)
1060 
1061 
1062  // output: coordinates of all vertices
1063 // for(i=0;i<Out.numberofpoints;i++)
1064 // OutPut(" (x, y, z) = "<< Out.pointlist[3*i]<<' '<<Out.pointlist[3*i+1]<<' '<<Out.pointlist[3*i+2]<<endl);
1065 
1066  Domain->GetTreeInfo(CellTree,N_RootCells);
1067  coll = Domain->GetCollection(It_Finest, 0);
1068  N_Cells = coll->GetN_Cells();
1069 
1070 // cout<<"N_RootCells: "<<N_RootCells<<endl;
1071  // remove all existing vertices and joints
1072  VertexDel = new TVertex*[8*N_RootCells];
1073  CurrVertex = 0;
1074 
1075  for(i=0;i<N_Cells;i++)
1076  {
1077  cell = coll->GetCell(i);
1078  N_Faces = cell->GetN_Faces();
1079  N_Vertices = cell->GetN_Vertices();
1080  for(j=0;j<N_Faces;j++)
1081  {
1082  if(CurrVertex==0)
1083  {
1084  VertexDel[CurrVertex++] = cell->GetVertex(j);
1085  }
1086  else
1087  {
1088  ID = 0;
1089  for(k=0;k<CurrVertex;k++)
1090  if(VertexDel[k]==cell->GetVertex(j))
1091  {
1092  ID = 1; break;
1093  }
1094  if(ID!=1)
1095  {
1096  VertexDel[CurrVertex++] = cell->GetVertex(j);
1097  }
1098  } // else if(CurrVertex==0)
1099 
1100  ID = 0;
1101  for(k=0;k<CurrVertex;k++)
1102  if(VertexDel[k]==cell->GetVertex((j+1)%N_Vertices))
1103  {
1104  ID = 1; break;
1105  }
1106  if(ID!=1)
1107  {
1108  VertexDel[CurrVertex++] = cell->GetVertex((j+1)%N_Vertices);
1109  }
1110 
1111  ID = 0;
1112  for(k=0;k<CurrVertex;k++)
1113  if(VertexDel[k]==cell->GetVertex((j+2)%N_Vertices))
1114  {
1115  ID = 1; break;
1116  }
1117  if(ID!=1)
1118  {
1119  VertexDel[CurrVertex++] = cell->GetVertex((j+2)%N_Vertices);
1120  }
1121  if(N_Faces==6) // If old cell is hexahedrol
1122  {
1123  ID = 0;
1124  for(k=0;k<CurrVertex;k++)
1125  if(VertexDel[k]==cell->GetVertex((j+4)%N_Vertices))
1126  {
1127  ID = 1; break;
1128  }
1129  if(ID!=1)
1130  {
1131  VertexDel[CurrVertex++] = cell->GetVertex((j+4)%N_Vertices);
1132  }
1133  }
1134  } // for j
1135  } // for i
1136 
1137  for(i=0;i<CurrVertex;i++)
1138  delete VertexDel[i];
1139  delete [] VertexDel;
1140  OutPut(CurrVertex<<" vertices were deleted"<<endl);
1141 
1142  // remove all existing cells and joints
1143 
1144  for(i=0;i<N_RootCells;i++)
1145  delete (TGridCell*)CellTree[i];
1146  delete [] CellTree;
1147  OutPut(N_RootCells<<" cells were deleted"<<endl);
1148 
1149 
1150 #ifdef _MPI
1151  if(rank==0)
1152 #endif
1153  {
1154  N_RootCells = Out.numberoftetrahedra;
1155  N_G = Out.numberofpoints;
1156 
1157  // allocate auxillary fields
1158  Coordinates = Out.pointlist;
1159  Tetrahedrals = Out.tetrahedronlist;
1160  }
1161 
1162 
1163 
1164 #ifdef _MPI
1165  MPI_Bcast(&N_RootCells, 1, MPI_INT, 0, Comm);
1166  MPI_Bcast(&N_G, 1, MPI_INT, 0, Comm);
1167 
1168 
1169  if(rank!=0)
1170  {
1171  Coordinates = new double [3*N_G];
1172  Tetrahedrals = new int[4*N_RootCells];
1173  }
1174 
1175  MPI_Bcast(Coordinates, 3*N_G, MPI_DOUBLE, 0, Comm);
1176  MPI_Bcast(Tetrahedrals, 4*N_RootCells, MPI_INT, 0, Comm);
1177 #endif
1178 
1179  // generate new vertices
1180  NewVertices = new TVertex*[N_G];
1181 
1182  for (i=0;i<N_G;i++)
1183  {
1184  NewVertices[i] = new TVertex(Coordinates[3*i], Coordinates[3*i+1], Coordinates[3*i+2]);
1185 
1186  // set bounding box
1187  if (Coordinates[3*i] > Xmax) Xmax = Coordinates[3*i];
1188  if (Coordinates[3*i] < Xmin) Xmin = Coordinates[3*i];
1189  if (Coordinates[3*i+1] > Ymax) Ymax = Coordinates[3*i+1];
1190  if (Coordinates[3*i+1] < Ymin) Ymin = Coordinates[3*i+1];
1191  if (Coordinates[3*i+2] > Zmax) Zmax = Coordinates[3*i+2];
1192  if (Coordinates[3*i+2] < Zmin) Zmin = Coordinates[3*i+2];
1193  }
1194  // set bounding box
1195  StartX = Xmin;
1196  StartY = Ymin;
1197  StartZ = Zmin;
1198  BoundX = Xmax - Xmin;
1199  BoundY = Ymax - Ymin;
1200  BoundZ = Zmax - Zmin;
1201 
1202  Domain->SetBoundBox(StartX, StartY, StartZ, BoundX, BoundY, BoundZ);
1203 // cout<<Xmin <<" "<<Ymin <<" "<<Zmin<<endl;
1204 // cout<<Xmax <<" "<<Ymax <<" "<<Zmax<<endl;
1205 
1206  CellTree = new TBaseCell*[N_RootCells];
1207  // output of each tetraheron vertex indices (four vertices for each)
1208 // for (i=0;i<N_RootCells;i++)
1209 // cout<< Tetrahedrals[4*i]<<" "<<Tetrahedrals[4*i + 1]<<" "
1210 // <<Tetrahedrals[4*i + 2]<<" "<<Tetrahedrals[4*i + 3]<<endl;
1211 
1212 
1213  for (i=0;i<N_RootCells;i++)
1214  {
1215  CellTree[i] = new TMacroCell(TDatabase::RefDescDB[Tetrahedron],RefLevel);
1216 
1217  CellTree[i]->SetVertex(0, NewVertices[Tetrahedrals[4*i ]]);
1218  CellTree[i]->SetVertex(1, NewVertices[Tetrahedrals[4*i + 1]]);
1219  CellTree[i]->SetVertex(2, NewVertices[Tetrahedrals[4*i + 2]]);
1220  CellTree[i]->SetVertex(3, NewVertices[Tetrahedrals[4*i + 3]]);
1221 
1222  CellTree[i]->SetClipBoard(i);
1223  ((TMacroCell *) CellTree[i])->SetSubGridID(0);
1224  }
1225 
1226  Domain->SetTreeInfo(CellTree, N_RootCells);
1227 
1228  // initialize iterators
1229  TDatabase::IteratorDB[It_EQ]->SetParam(Domain);
1230  TDatabase::IteratorDB[It_LE]->SetParam(Domain);
1231  TDatabase::IteratorDB[It_Finest]->SetParam(Domain);
1232  TDatabase::IteratorDB[It_Between]->SetParam(Domain);
1233  TDatabase::IteratorDB[It_OCAF]->SetParam(Domain);
1234 
1235 
1236  // search neighbours
1237  PointNeighb = new int[N_G];
1238 #ifdef _MPI
1239  if(rank==0)
1240 #endif
1241  cout<<"numberofpoints "<<N_G<<endl;
1242  memset(PointNeighb, 0, N_G *SizeOfInt);
1243 
1244  for (i=0;i<4*N_RootCells;i++)
1245  PointNeighb[Tetrahedrals[i]]++;
1246 
1247  for (i=0;i<N_G;i++)
1248  if (PointNeighb[i] > maxEpV) maxEpV = PointNeighb[i];
1249  delete [] PointNeighb;
1250 
1251 #ifdef _MPI
1252  if(rank==0)
1253 #endif
1254  cout<<"maxEpV "<< maxEpV<<endl;
1255 
1256  PointNeighb = new int[++maxEpV * N_G];
1257 
1258  memset(PointNeighb, 0, maxEpV*N_G*SizeOfInt);
1259 
1260  // every vertex contains "maxEpV" columns
1261  // for every vertex at first colomn contains the number of cells containing this vertex
1262  // at further columns we set the index of corresponding cells containing this vertex
1263  // cout<<"maxEpV*N_G "<<maxEpV*N_G<<endl;
1264 
1265  for(i=0;i<4*N_RootCells;i++)
1266  {
1267  j = Tetrahedrals[i]*maxEpV;
1268  PointNeighb[j]++;
1269  //cout<<"j + PointNeighb[j] " << j <<endl;
1270  PointNeighb[j + PointNeighb[j]] = i / 4;
1271  }
1272  // output of PointNeighb columns for each point
1273 // for (i=0;i<N_G;i++)
1274 // {
1275 // for (j=0;j<maxEpV;j++)
1276 // cout<<" "<< PointNeighb[i*maxEpV+j];
1277 // cout<<endl;
1278 // }
1279 
1280  // generate new faces
1281 
1282 #ifdef _MPI
1283  if(rank==0)
1284 #endif
1285  {
1286  N_G = Out.numberoftrifaces;
1287  Facelist = Out.trifacelist;
1288  Facemarkerlist = Out.trifacemarkerlist;
1289  }
1290 
1291 #ifdef _MPI
1292  MPI_Bcast(&N_G, 1, MPI_INT, 0, Comm);
1293 
1294  if(rank!=0)
1295  {
1296  Facelist = new int [3*N_G];
1297  Facemarkerlist = new int[N_G];
1298  }
1299 
1300 
1301  MPI_Bcast(Facelist, 3*N_G, MPI_INT, 0, Comm);
1302  MPI_Bcast(Facemarkerlist, N_G, MPI_INT, 0, Comm);
1303 
1304  if(rank==0)
1305 #endif
1306  cout<<"numberoftrifaces "<<N_G<<endl;
1307 
1308  for (i=0;i<N_G;i++)
1309  {
1310  a = Facelist[3*i];
1311  b = Facelist[3*i+1];
1312  c = Facelist[3*i+2];
1313 
1314 // cout<<" "<< a<<" "<< b<<" "<< c<<endl;
1315 
1316  Neib[0] = -1;
1317  Neib[1] = -1;
1318  CurrNeib = 0;
1319 
1320  len1 = PointNeighb[a*maxEpV];
1321  len2 = PointNeighb[b*maxEpV];
1322  len3 = PointNeighb[c*maxEpV];
1323 
1324  // find the index of the cells containing current face with point indices a,b,c
1325  for (j=1;j<=len1;j++)
1326  {
1327  Neighb_tmp = PointNeighb[a*maxEpV + j];
1328  for (k=1;k<=len2;k++)
1329  {
1330  if (Neighb_tmp == PointNeighb[b*maxEpV + k])
1331  {
1332  for (l=1;l<=len3;l++)
1333  if (Neighb_tmp == PointNeighb[c*maxEpV + l])
1334  {
1335  Neib[CurrNeib++] = Neighb_tmp;
1336  break;
1337  }
1338  }
1339  }
1340  if (CurrNeib == 2) break;
1341  }
1342  // cout<<"CurrNeib " << CurrNeib <<endl;
1343 // cout<<"Facemarkerlist[i] : "<<Facemarkerlist[i]<<endl;
1344  if (Facemarkerlist[i]) // 0 for inner edges and Boundcomp+1 for Boundedge respect
1345  {
1346 
1347  CurrComp = Facemarkerlist[i] - 1;
1348 // cout<<"Boundary face CurrComp: "<<CurrComp<<endl;
1349 // exit(0);
1350  CurrComp = 0; // not yet implemented fully
1351 
1352 
1353 
1354  bdcomp = Domain->GetBdPart(0)->GetBdComp(CurrComp);
1355 
1356 
1357  if(bdcomp->GetTSofXYZ(NewVertices[a]->GetX(), NewVertices[a]->GetY(),
1358  NewVertices[a]->GetY(), T[1], S[1]) ||
1359  bdcomp->GetTSofXYZ(NewVertices[b]->GetX(), NewVertices[b]->GetY(),
1360  NewVertices[b]->GetY(), T[2], S[2]) ||
1361  bdcomp->GetTSofXYZ(NewVertices[c]->GetX(), NewVertices[c]->GetY(),
1362  NewVertices[c]->GetY(), T[3], S[3]) )
1363  {
1364  cerr<<"Error: could not set parameter values"<<endl;
1365  OutPut(NewVertices[a]<<endl);
1366  OutPut(NewVertices[b]<<endl);
1367  OutPut(NewVertices[c]<<endl);
1368  exit(0);
1369  }
1370 
1371  if (CurrNeib == 2)
1372  {
1373  if(bdcomp->IsFreeBoundary())
1374  Joint = new TIsoInterfaceJoint3D(bdcomp, T, S,
1375  CellTree[Neib[0]], CellTree[Neib[1]] );
1376  else
1377  Joint = new TInterfaceJoint3D(bdcomp, T, S,
1378  CellTree[Neib[0]], CellTree[Neib[1]] );
1379  }
1380  else
1381  {
1382  if(bdcomp->IsFreeBoundary())
1383  Joint = new TIsoBoundFace(bdcomp, T, S);
1384  else
1385  Joint = new TBoundFace(bdcomp, T, S);
1386  }
1387  }
1388  else
1389  {
1390 // // cout<<"Inner face"<<endl;
1391  if (CurrNeib != 2)
1392  cerr << "Error !!!!!!!! not enough neighbours!" << endl;
1393 
1394  Joint = new TJointEqN(CellTree[Neib[0]], CellTree[Neib[1]]);
1395 
1396  }
1397 
1398  // First element containing the current face
1399  // find the local index for the point 'a' on the cell
1400  for (j=0;j<4;j++)
1401  if (Tetrahedrals[4*Neib[0]+j] == a) break;
1402 
1403  // find the local index for the point 'b' on the cell
1404  for (k=0;k<4;k++)
1405  if (Tetrahedrals[4*Neib[0]+k] == b) break;
1406 
1407  // find the local index for the point 'c' on the cell
1408  for (l=0;l<4;l++)
1409  if (Tetrahedrals[4*Neib[0]+l] == c) break;
1410 
1411  l = l*100 + k*10 + j;
1412 
1413  switch (l) // j will contain the local index for the current face
1414  {
1415  case 210: case 21: case 102:
1416  case 120: case 12: case 201:
1417  j = 0;
1418  break;
1419  case 310: case 31: case 103:
1420  case 130: case 13: case 301:
1421  j = 1;
1422  break;
1423  case 321: case 132: case 213:
1424  case 231: case 123: case 312:
1425  j = 2;
1426  break;
1427  case 230: case 23: case 302:
1428  case 320: case 32: case 203:
1429  j = 3;
1430  break;
1431 
1432  default:
1433  Error("Unable to set the face !!!!!!!!!!!!" << endl);
1434  exit(0);
1435  }
1436  CellTree[Neib[0]]->SetJoint(j, Joint);
1437 
1438  if (Neib[1] != -1) // second element containing the current face
1439  {
1440  // find the local index for the point 'a' on the cell
1441  for (j=0;j<4;j++)
1442  if (Tetrahedrals[4*Neib[1]+j] == a) break;
1443 
1444  // find the local index for the point 'b' on the cell
1445  for (k=0;k<4;k++)
1446  if (Tetrahedrals[4*Neib[1]+k] == b) break;
1447 
1448  // find the local index for the point 'c' on the cell
1449  for (l=0;l<4;l++)
1450  if (Tetrahedrals[4*Neib[1]+l] == c) break;
1451 
1452  l = l*100 + k*10 + j;
1453 
1454 // cout<<""<< l <<endl;
1455 
1456  switch (l) // j will contain the local index for the current face
1457  {
1458  case 210: case 21: case 102:
1459  case 120: case 12: case 201:
1460  j = 0;
1461  break;
1462  case 310: case 31: case 103:
1463  case 130: case 13: case 301:
1464  j = 1;
1465  break;
1466  case 321: case 132: case 213:
1467  case 231: case 123: case 312:
1468  j = 2;
1469  break;
1470  case 230: case 23: case 302:
1471  case 320: case 32: case 203:
1472  j = 3;
1473  break;
1474 
1475  default:
1476  Error("Unable to set the face !!!!!!!!!!!!" << endl);
1477  exit(0);
1478  }
1479  CellTree[Neib[1]]->SetJoint(j, Joint);
1480  }
1481 
1482  if (Joint->GetType() == InterfaceJoint3D ||
1483  Joint->GetType() == IsoInterfaceJoint3D)
1484  {
1485  ((TInterfaceJoint3D*)Joint)->SetMapType();
1486  ((TInterfaceJoint3D*)(Joint))->CheckOrientation();
1487  }
1488  else
1489  if (Joint->GetType() == JointEqN)
1490  Joint->SetMapType();
1491 
1492  }
1493 
1494 #ifdef _MPI
1495  if(rank==0)
1496 #endif
1497  {
1498  }
1499 
1500 #ifdef _MPI
1501  if(rank!=0)
1502  {
1503  delete [] Tetrahedrals ;
1504  delete [] Coordinates;
1505  delete [] Facelist;
1506  delete [] Facemarkerlist;
1507  }
1508 #endif
1509  delete [] PointNeighb;
1510 }
1511 
1513 void TetraGen(TDomain *&Domain)
1514 {
1515  //======================================================================
1516  // Tetgen for grid generation begin
1517  //======================================================================
1518  int i, j, k, l, N_Coord, *N_FVerts, N_Faces, *Facets;
1519  int N, N_RootCells, N_Cells, CurrVertex, N_Vertices, ID, N_G, RefLevel=0;
1520  int CurrNeib, len1, len2, len3, maxEpV = 0, a, b, c, Neib[2], Neighb_tmp, CurrComp;
1521  int *Tetrahedrals, *PointNeighb, dimension;
1522  int jj, N_Points, MaxLen, *Tetrahedrals_loc;
1523  const int *EdgeVertex, *TmpFV, *TmpLen;
1524 
1525  double *Coordinates, N_x, N_y, N_z;
1526  double *Vertices;
1527  double Xmin = 1e10, Xmax = -1e10, Ymin = 1e10, Ymax = -1e10;
1528  double Zmin = 1e10, Zmax = -1e10, T[4]={0,0,0,0}, S[4]={0,0,0,0};
1529  double StartX, StartY, StartZ, BoundX, BoundY, BoundZ;
1530 
1531  tetgenio In, Out;
1532 
1533  TBaseCell **CellTree, **SurfCellTree;
1534  TGridCell **DelCell;
1535  TVertex **VertexDel, **NewVertices, **NewSurfVertices;
1536  TBoundPart *BoundPart;
1537  TBdPlane **UpdateFaceParams;
1538  TJoint *Joint;
1539  TBoundComp3D *bdcomp;
1540  TCollection *coll, *SurfColl;
1541  TBaseCell *cell;
1542  TBdSphere *UpdateParam;
1543  TShapeDesc *ShapeDesc;
1544 
1545  char *SMESH, line[100];
1546  SMESH = TDatabase::ParamDB->SMESHFILE;
1547 
1548 #ifdef _MPI
1549  int rank, size;
1550  MPI_Comm Comm = TDatabase::ParamDB->Comm;
1551  MPI_Comm_rank(Comm, &rank);
1552  MPI_Comm_size(Comm, &size);
1553 #endif
1554 
1556 #ifdef _MPI
1557  if(rank==0)
1558 #endif
1559  {
1560 // In.initialize();
1561 // Out.initialize();
1562  std::ostringstream opts;
1563  opts << " ";
1564 
1565  opts.seekp(std::ios::beg);
1566  opts<<'p'; // Tetrahedralize the PLC. Switches are chosen to read a PLC (p)
1567 // opts<<'r'; // -r Reconstructs a previously generated mesh
1568  opts<<"q"<<1.20; // quality mesh generation(q) with a specified quality bound
1569 // opts<<"a"<<0.1; // maximum volume constraint
1570 // opts<<'i'; // Inserts a list of additional points into mesh.
1571  opts<<'z'; // numbers all output items starting from zero
1572 // opts<<'d'; // Detect intersections of PLC facets.
1573  opts<<'f'; // Outputs all faces (including non-boundary)
1574  opts<<'e'; // Outputs a list of edges of the triangulation
1575 // opts<<'I'; // Suppresses mesh iteration numbers.
1576  opts<<'C'; // Checks the consistency of the final mesh.
1577 // opts<<'Q'; // Quiet: No terminal output except errors.
1578 // opts<<'g'; // Outputs mesh to .mesh file for viewing by Medit
1579  opts<<'Y'; // Suppresses boundary facets/segments splitting
1580 // opts<<'V'; //verbose mode
1581  opts<<ends;
1582 
1583 // cout << " SMESHFILE is " << SMESH << endl;
1585  ReadMeditMesh(SMESH, In);
1586 
1587 // for(i=0;i<In.numberofpoints;i++)
1588 // OutPut(i<<" (x, y, z) = "<<
1589 // In.pointlist[3*i]<<' '<<In.pointlist[3*i+1]<<' '<<In.pointlist[3*i+2]<<endl);
1590 
1591  // Calling tetrahedralize function of 3dtetgen mesh generator
1592  tetrahedralize((char*)opts.str().c_str(), &In, &Out);
1593 
1594  // output: coordinates of all vertices
1595 // for(i=0;i<Out.numberofpoints;i++)
1596 // OutPut(" (x, y, z) = "<< Out.pointlist[3*i]<<' '<<Out.pointlist[3*i+1]<<' '<<Out.pointlist[3*i+2]<<endl);
1597 
1598  } // if(rank==0)
1599 
1601  DeleteDomain(Domain);
1602 
1603 #ifdef _MPI
1604  if(rank==0)
1605 #endif
1606  {
1607  N_RootCells = Out.numberoftetrahedra;
1608  N_G = Out.numberofpoints;
1609  // allocate auxillary fields
1610  Coordinates = Out.pointlist;
1611  Tetrahedrals = Out.tetrahedronlist;
1612  }
1613 
1614 #ifdef _MPI
1615  MPI_Bcast(&N_RootCells, 1, MPI_INT, 0, Comm);
1616  MPI_Bcast(&N_G, 1, MPI_INT, 0, Comm);
1617  if(rank!=0)
1618  {
1619  Coordinates = new double [3*N_G];
1620  Tetrahedrals = new int[4*N_RootCells];
1621  }
1622  MPI_Bcast(Coordinates, 3*N_G, MPI_DOUBLE, 0, Comm);
1623  MPI_Bcast(Tetrahedrals, 4*N_RootCells, MPI_INT, 0, Comm);
1624 #endif
1625 
1626  // generate new vertices
1627  NewVertices = new TVertex*[N_G];
1628 
1629  for (i=0;i<N_G;i++)
1630  {
1631  NewVertices[i] = new TVertex(Coordinates[3*i], Coordinates[3*i+1], Coordinates[3*i+2]);
1632  // set bounding box
1633  if (Coordinates[3*i] > Xmax) Xmax = Coordinates[3*i];
1634  if (Coordinates[3*i] < Xmin) Xmin = Coordinates[3*i];
1635  if (Coordinates[3*i+1] > Ymax) Ymax = Coordinates[3*i+1];
1636  if (Coordinates[3*i+1] < Ymin) Ymin = Coordinates[3*i+1];
1637  if (Coordinates[3*i+2] > Zmax) Zmax = Coordinates[3*i+2];
1638  if (Coordinates[3*i+2] < Zmin) Zmin = Coordinates[3*i+2];
1639  }
1640 
1641  // set bounding box
1642  StartX = Xmin;
1643  StartY = Ymin;
1644  StartZ = Zmin;
1645  BoundX = Xmax - Xmin;
1646  BoundY = Ymax - Ymin;
1647  BoundZ = Zmax - Zmin;
1648 
1649  Domain->SetBoundBox(StartX, StartY, StartZ, BoundX, BoundY, BoundZ);
1650 // cout<<Xmin <<" "<<Ymin <<" "<<Zmin<<endl;
1651 // cout<<Xmax <<" "<<Ymax <<" "<<Zmax<<endl;
1652  CellTree = new TBaseCell*[N_RootCells];
1653  // output of each tetraheron vertex indices (four vertices for each)
1654 // for (i=0;i<N_RootCells;i++)
1655 // cout<< Tetrahedrals[4*i]<<" "<<Tetrahedrals[4*i + 1]<<" "
1656 // <<Tetrahedrals[4*i + 2]<<" "<<Tetrahedrals[4*i + 3]<<endl;
1657 
1658 
1659  for (i=0;i<N_RootCells;i++)
1660  {
1661  CellTree[i] = new TMacroCell(TDatabase::RefDescDB[Tetrahedron], RefLevel);
1662 
1663  CellTree[i]->SetVertex(0, NewVertices[Tetrahedrals[4*i ]]);
1664  CellTree[i]->SetVertex(1, NewVertices[Tetrahedrals[4*i + 1]]);
1665  CellTree[i]->SetVertex(2, NewVertices[Tetrahedrals[4*i + 2]]);
1666  CellTree[i]->SetVertex(3, NewVertices[Tetrahedrals[4*i + 3]]);
1667  CellTree[i]->SetClipBoard(i);
1668  ((TMacroCell *) CellTree[i])->SetSubGridID(0);
1669  }
1670 
1671  Domain->SetTreeInfo(CellTree, N_RootCells);
1672  // initialize iterators
1673  TDatabase::IteratorDB[It_EQ]->SetParam(Domain);
1674  TDatabase::IteratorDB[It_LE]->SetParam(Domain);
1675  TDatabase::IteratorDB[It_Finest]->SetParam(Domain);
1676  TDatabase::IteratorDB[It_Between]->SetParam(Domain);
1677  TDatabase::IteratorDB[It_OCAF]->SetParam(Domain);
1678  // search neighbours
1679  PointNeighb = new int[N_G];
1680 #ifdef _MPI
1681  if(rank==0)
1682 #endif
1683  cout<<"numberofpoints "<<N_G<<endl;
1684  memset(PointNeighb, 0, N_G *SizeOfInt);
1685 
1686  for (i=0;i<4*N_RootCells;i++)
1687  PointNeighb[Tetrahedrals[i]]++;
1688 
1689  for (i=0;i<N_G;i++)
1690  if (PointNeighb[i] > maxEpV) maxEpV = PointNeighb[i];
1691  delete [] PointNeighb;
1692 
1693 #ifdef _MPI
1694  if(rank==0)
1695 #endif
1696  cout<<"maxEpV "<< maxEpV<<endl;
1697  PointNeighb = new int[++maxEpV * N_G];
1698  memset(PointNeighb, 0, maxEpV*N_G*SizeOfInt);
1699  // every vertex contains "maxEpV" columns
1700  // for every vertex at first colomn contains the number of cells containing this vertex
1701  // at further columns we set the index of corresponding cells containing this vertex
1702  // cout<<"maxEpV*N_G "<<maxEpV*N_G<<endl;
1703 
1704  for(i=0;i<4*N_RootCells;i++)
1705  {
1706  j = Tetrahedrals[i]*maxEpV;
1707  PointNeighb[j]++;
1708  //cout<<"j + PointNeighb[j] " << j <<endl;
1709  PointNeighb[j + PointNeighb[j]] = i / 4;
1710  }
1711  // generate new faces
1712  coll=Domain->GetCollection(It_Finest, 0);
1713  N_Cells = coll->GetN_Cells();
1714  N_G = 0;
1715  for (i=0;i<N_Cells;i++)
1716  {
1717  cell = coll->GetCell(i);
1718  ShapeDesc= cell->GetShapeDesc();
1719  ShapeDesc->GetFaceVertex(TmpFV, TmpLen, MaxLen);
1720  N_Faces = cell->GetN_Faces();
1721  Tetrahedrals_loc = Tetrahedrals+4*i;
1722 
1723  for(jj=0;jj<N_Faces;jj++)
1724  if(cell->GetJoint(jj) == NULL)
1725  {
1726  N_G++;
1727  N_Points = TmpLen[jj];
1728 
1729  if(N_Points!=3)
1730  {
1731  cerr << "Only Tria faces are allowed!!! N_FVert: " << N_Points << endl;
1732  exit(-1);
1733  }
1734 
1735  //printf(" TetrameshGen Null joint \n" );
1736  a = Tetrahedrals_loc[TmpFV[jj*MaxLen]];
1737  b = Tetrahedrals_loc[TmpFV[jj*MaxLen + 1]];
1738  c = Tetrahedrals_loc[TmpFV[jj*MaxLen + 2]];
1739 // cout<<" "<< a<<" "<< b<<" "<< c<<endl;
1740  Neib[0] = -1;
1741  Neib[1] = -1;
1742  CurrNeib = 0;
1743 
1744  len1 = PointNeighb[a*maxEpV];
1745  len2 = PointNeighb[b*maxEpV];
1746  len3 = PointNeighb[c*maxEpV];
1747 
1748  // find the index of the cells containing current face with point indices a,b,c
1749  for (j=1;j<=len1;j++)
1750  {
1751  Neighb_tmp = PointNeighb[a*maxEpV + j];
1752  for (k=1;k<=len2;k++)
1753  {
1754  if (Neighb_tmp == PointNeighb[b*maxEpV + k])
1755  {
1756  for (l=1;l<=len3;l++)
1757  if (Neighb_tmp == PointNeighb[c*maxEpV + l])
1758  {
1759  Neib[CurrNeib++] = Neighb_tmp;
1760  break;
1761  }
1762  }
1763  }
1764  if (CurrNeib == 2) break;
1765  }// for (j=1;j<=len1;j++)
1766  // cout<<"CurrNeib " << CurrNeib <<endl;
1767  // cout<<"Facemarkerlist[i] : "<<Facemarkerlist[i]<<endl;
1768  if (CurrNeib == 1) // 0 for inner edges and Boundcomp+1 for Boundedge respect
1769  {
1770  CurrComp = 0; // not yet implemented fully
1771 
1772  bdcomp = Domain->GetBdPart(0)->GetBdComp(CurrComp);
1773 
1774  if(bdcomp->GetTSofXYZ(NewVertices[a]->GetX(), NewVertices[a]->GetY(),
1775  NewVertices[a]->GetY(), T[1], S[1]) ||
1776  bdcomp->GetTSofXYZ(NewVertices[b]->GetX(), NewVertices[b]->GetY(),
1777  NewVertices[b]->GetY(), T[2], S[2]) ||
1778  bdcomp->GetTSofXYZ(NewVertices[c]->GetX(), NewVertices[c]->GetY(),
1779  NewVertices[c]->GetY(), T[3], S[3]) )
1780  {
1781  cerr<<"Error: could not set parameter values"<<endl;
1782  OutPut(NewVertices[a]<<endl);
1783  OutPut(NewVertices[b]<<endl);
1784  OutPut(NewVertices[c]<<endl);
1785  exit(0);
1786  }
1787 
1788  if (CurrNeib == 2)
1789  {
1790  if(bdcomp->IsFreeBoundary())
1791  Joint = new TIsoInterfaceJoint3D(bdcomp, T, S,
1792  CellTree[Neib[0]], CellTree[Neib[1]] );
1793  else
1794  Joint = new TInterfaceJoint3D(bdcomp, T, S,
1795  CellTree[Neib[0]], CellTree[Neib[1]] );
1796  }
1797  else
1798  {
1799  if(bdcomp->IsFreeBoundary())
1800  Joint = new TIsoBoundFace(bdcomp, T, S);
1801  else
1802  Joint = new TBoundFace(bdcomp, T, S);
1803  }
1804  }
1805  else
1806  {
1807  // cout<<"Inner face"<<endl;
1808  if (CurrNeib != 2)
1809  cerr << "Error !!!!!!!! not enough neighbours!" << endl;
1810 
1811  Joint = new TJointEqN(CellTree[Neib[0]], CellTree[Neib[1]]);
1812  }
1813 
1814  // First element containing the current face
1815  // find the local index for the point 'a' on the cell
1816  for (j=0;j<4;j++)
1817  if (Tetrahedrals[4*Neib[0]+j] == a) break;
1818 
1819  // find the local index for the point 'b' on the cell
1820  for (k=0;k<4;k++)
1821  if (Tetrahedrals[4*Neib[0]+k] == b) break;
1822 
1823  // find the local index for the point 'c' on the cell
1824  for (l=0;l<4;l++)
1825  if (Tetrahedrals[4*Neib[0]+l] == c) break;
1826 
1827  l = l*100 + k*10 + j;
1828 
1829 // cout<<""<< l <<endl;
1830 
1831  switch (l) // j will contain the local index for the current face
1832  {
1833  case 210: case 21: case 102:
1834  case 120: case 12: case 201:
1835  j = 0;
1836  break;
1837  case 310: case 31: case 103:
1838  case 130: case 13: case 301:
1839  j = 1;
1840  break;
1841  case 321: case 132: case 213:
1842  case 231: case 123: case 312:
1843  j = 2;
1844  break;
1845  case 230: case 23: case 302:
1846  case 320: case 32: case 203:
1847  j = 3;
1848  break;
1849 
1850  default:
1851  Error("Unable to set the face !!!!!!!!!!!!" << endl);
1852  exit(0);
1853  }
1854  CellTree[Neib[0]]->SetJoint(j, Joint);
1855 
1856  if (Neib[1] != -1) // second element containing the current face
1857  {
1858  // find the local index for the point 'a' on the cell
1859  for (j=0;j<4;j++)
1860  if (Tetrahedrals[4*Neib[1]+j] == a) break;
1861 
1862  // find the local index for the point 'b' on the cell
1863  for (k=0;k<4;k++)
1864  if (Tetrahedrals[4*Neib[1]+k] == b) break;
1865 
1866  // find the local index for the point 'c' on the cell
1867  for (l=0;l<4;l++)
1868  if (Tetrahedrals[4*Neib[1]+l] == c) break;
1869 
1870  l = l*100 + k*10 + j;
1871 // cout<<""<< l <<endl;
1872  switch (l) // j will contain the local index for the current face
1873  {
1874  case 210: case 21: case 102:
1875  case 120: case 12: case 201:
1876  j = 0;
1877  break;
1878  case 310: case 31: case 103:
1879  case 130: case 13: case 301:
1880  j = 1;
1881  break;
1882  case 321: case 132: case 213:
1883  case 231: case 123: case 312:
1884  j = 2;
1885  break;
1886  case 230: case 23: case 302:
1887  case 320: case 32: case 203:
1888  j = 3;
1889  break;
1890 
1891  default:
1892  Error("Unable to set the face !!!!!!!!!!!!" << endl);
1893  exit(0);
1894  }
1895  CellTree[Neib[1]]->SetJoint(j, Joint);
1896  } // if (Neib[1] != -1)
1897 
1898  if (Joint->GetType() == InterfaceJoint3D ||
1899  Joint->GetType() == IsoInterfaceJoint3D)
1900  {
1901  ((TInterfaceJoint3D*)Joint)->SetMapType();
1902  ((TInterfaceJoint3D*)(Joint))->CheckOrientation();
1903  }
1904  else
1905  if (Joint->GetType() == JointEqN)
1906  Joint->SetMapType();
1907  } // if(cell->GetJoint(jj) == NULL)
1908  } // for (i=0;i<N_Cells;i++)
1909 
1910 #ifdef _MPI
1911  if(rank==0)
1912 #endif
1913  cout<<"numberoftrifaces after "<<N_G<<endl;
1914 
1915 #ifdef _MPI
1916  if(rank==0)
1917 #endif
1918  {
1919 // In.deinitialize();
1920 // Out.deinitialize();
1921  }
1922 #ifdef _MPI
1923  if(rank!=0)
1924  {
1925  delete [] Tetrahedrals ;
1926  delete [] Coordinates;
1927 // delete [] Facelist;
1928 // delete [] Facemarkerlist;
1929  }
1930 #endif
1931  delete [] PointNeighb;
1932 // #ifdef _MPI
1933 // if(rank==3)
1934 // printf(" TetrameshGen %d \n", N_G );
1935 // MPI_Finalize();
1936 // #endif
1937 // exit(0);
1938 //
1939 }
1940 
1942 void TetraGenWithInputCells(TDomain *&Domain)
1943 {
1944  //======================================================================
1945  // Tetgen for grid generation begin
1946  //======================================================================
1947  int i, j, k, l, N_Coord, *N_FVerts, N_Faces, *Facets;
1948  int N, N_RootCells, N_Cells, CurrVertex, N_Vertices, ID, N_G, RefLevel=0;
1949  int CurrNeib, len1, len2, len3, maxEpV = 0, a, b, c, Neib[2], Neighb_tmp, CurrComp;
1950  int *Tetrahedrals, *PointNeighb, dimension;
1951  int jj, N_Points, MaxLen, *Tetrahedrals_loc;
1952  const int *EdgeVertex, *TmpFV, *TmpLen;
1953  double *TetRegionlist;
1954 
1955  double *Coordinates, N_x, N_y, N_z;
1956  double *Vertices;
1957  double Xmin = 1e10, Xmax = -1e10, Ymin = 1e10, Ymax = -1e10;
1958  double Zmin = 1e10, Zmax = -1e10, T[4]={0,0,0,0}, S[4]={0,0,0,0};
1959  double StartX, StartY, StartZ, BoundX, BoundY, BoundZ;
1960 
1961  tetgenio In, Out, InAddPts;
1962 
1963  TBaseCell **CellTree, **SurfCellTree;
1964  TGridCell **DelCell;
1965  TVertex **VertexDel, **NewVertices, **NewSurfVertices;
1966  TBoundPart *BoundPart;
1967  TBdPlane **UpdateFaceParams;
1968  TJoint *Joint;
1969  TBoundComp3D *bdcomp;
1970  TCollection *coll, *SurfColl;
1971  TBaseCell *cell;
1972  TBdSphere *UpdateParam;
1973  TShapeDesc *ShapeDesc;
1974 
1975  char *SMESH, *NODE, line[100];
1976  SMESH = TDatabase::ParamDB->SMESHFILE;
1977  NODE = TDatabase::ParamDB->GRAPEBASENAME;
1978 
1979 #ifdef _MPI
1980  int rank, size;
1981  MPI_Comm Comm = TDatabase::ParamDB->Comm;
1982  MPI_Comm_rank(Comm, &rank);
1983  MPI_Comm_size(Comm, &size);
1984 #endif
1985 
1987 #ifdef _MPI
1988  if(rank==0)
1989 #endif
1990  {
1991 // In.initialize();
1992 // Out.initialize();
1993  std::ostringstream opts;
1994  opts << " ";
1995 
1996  opts.seekp(std::ios::beg);
1997 // opts<<'p'; // Tetrahedralize the PLC. Switches are chosen to read a PLC (p)
1998  opts<<'r'; // -r Reconstructs a previously generated mesh
1999 // opts<<'q'<<1.2; // quality mesh generation(q) with a specified quality bound
2000 // opts<<'a'<<100; // maximum volume constraint
2001 // opts<<'a'; // maximum volume constraint
2002 // opts<<'i'; // Inserts a list of additional points into mesh.
2003  opts<<'z'; // numbers all output items starting from zero
2004 // opts<<'d'; // Detect intersections of PLC facets.
2005  opts<<'f'; // Outputs all faces (including non-boundary)
2006  opts<<'e'; // Outputs a list of edges of the triangulation
2007 // opts<<'I'; // Suppresses mesh iteration numbers.
2008  opts<<'C'; // Checks the consistency of the final mesh.
2009 // opts<<'Q'; // Quiet: No terminal output except errors.
2010 // opts<<'g'; // Outputs mesh to .mesh file for viewing by Medit
2011  opts<<'Y'; // Suppresses boundary facets/segments splitting
2012 // opts<<'A'; // output region attribute list
2013 // opts<<'V'; //verbose mode
2014 // opts<<'V';
2015 // opts<<'V';
2016 // opts<<'V';
2017 // opts<<'V';
2018  opts<<ends;
2019 
2020  In.initialize();
2021  Out.initialize();
2022 // InAddPts.initialize();
2023 // cout << " SMESHFILE is " << SMESH << endl;
2025  ReadMeditMeshWithCells(SMESH, In);
2026 // ReadAdditionalNModes(NODE, InAddPts);
2027 // In.load_medit(SMESH, 1);
2028 
2029 // InAddPts.load_node(NODE);
2030 // cout<< " InAddPts.numberofpoints " << InAddPts.numberofpoints <<endl;
2031 // exit(0);
2032 // for(i=0;i<In.numberofpoints;i++)
2033 // OutPut(i<<" (x, y, z) = "<<
2034 // In.pointlist[3*i]<<' '<<In.pointlist[3*i+1]<<' '<<In.pointlist[3*i+2]<<endl);
2035 // exit(0);
2036  // Calling tetrahedralize function of 3dtetgen mesh generator
2037  tetrahedralize((char*)opts.str().c_str(), &In, &Out);
2038 // tetrahedralize((char*)opts.str().c_str(), &In, &Out, &InAddPts);
2039 // exit(0);
2040  // output: coordinates of all vertices
2041 // for(i=0;i<Out.numberofpoints;i++)
2042 // OutPut(" (x, y, z) = "<< Out.pointlist[3*i]<<' '<<Out.pointlist[3*i+1]<<' '<<Out.pointlist[3*i+2]<<endl);
2043 
2044  } // if(rank==0)
2045 
2047  DeleteDomain(Domain);
2048 
2049 #ifdef _MPI
2050  if(rank==0)
2051 #endif
2052  {
2053  N_RootCells = Out.numberoftetrahedra;
2054  N_G = Out.numberofpoints;
2055  // allocate auxillary fields
2056  Coordinates = Out.pointlist;
2057  Tetrahedrals = Out.tetrahedronlist;
2058  TetRegionlist = Out.tetrahedronattributelist;
2059  }
2060 
2061 #ifdef _MPI
2062  MPI_Bcast(&N_RootCells, 1, MPI_INT, 0, Comm);
2063  MPI_Bcast(&N_G, 1, MPI_INT, 0, Comm);
2064 
2065  if(rank!=0)
2066  {
2067  Coordinates = new double [3*N_G];
2068  Tetrahedrals = new int[4*N_RootCells];
2069  TetRegionlist = new double [N_RootCells];
2070  }
2071 
2072  MPI_Bcast(Coordinates, 3*N_G, MPI_DOUBLE, 0, Comm);
2073  MPI_Bcast(Tetrahedrals, 4*N_RootCells, MPI_INT, 0, Comm);
2074  MPI_Bcast(TetRegionlist, N_RootCells, MPI_DOUBLE, 0, Comm);
2075 #endif
2076 
2077  // generate new vertices
2078  NewVertices = new TVertex*[N_G];
2079 
2080  for (i=0;i<N_G;i++)
2081  {
2082  NewVertices[i] = new TVertex(Coordinates[3*i], Coordinates[3*i+1], Coordinates[3*i+2]);
2083 
2084  // set bounding box
2085  if (Coordinates[3*i] > Xmax) Xmax = Coordinates[3*i];
2086  if (Coordinates[3*i] < Xmin) Xmin = Coordinates[3*i];
2087  if (Coordinates[3*i+1] > Ymax) Ymax = Coordinates[3*i+1];
2088  if (Coordinates[3*i+1] < Ymin) Ymin = Coordinates[3*i+1];
2089  if (Coordinates[3*i+2] > Zmax) Zmax = Coordinates[3*i+2];
2090  if (Coordinates[3*i+2] < Zmin) Zmin = Coordinates[3*i+2];
2091  }
2092 
2093  // set bounding box
2094  StartX = Xmin;
2095  StartY = Ymin;
2096  StartZ = Zmin;
2097  BoundX = Xmax - Xmin;
2098  BoundY = Ymax - Ymin;
2099  BoundZ = Zmax - Zmin;
2100  Domain->SetBoundBox(StartX, StartY, StartZ, BoundX, BoundY, BoundZ);
2101 // cout<<Xmin <<" "<<Ymin <<" "<<Zmin<<endl;
2102 // cout<<Xmax <<" "<<Ymax <<" "<<Zmax<<endl;
2103  CellTree = new TBaseCell*[N_RootCells];
2104  // output of each tetraheron vertex indices (four vertices for each)
2105 // for (i=0;i<N_RootCells;i++)
2106 // cout<< Tetrahedrals[4*i]<<" "<<Tetrahedrals[4*i + 1]<<" "
2107 // <<Tetrahedrals[4*i + 2]<<" "<<Tetrahedrals[4*i + 3]<<endl;
2108 
2109 
2110  for (i=0;i<N_RootCells;i++)
2111  {
2112  CellTree[i] = new TMacroCell(TDatabase::RefDescDB[Tetrahedron], RefLevel);
2113  CellTree[i]->SetVertex(0, NewVertices[Tetrahedrals[4*i ]]);
2114  CellTree[i]->SetVertex(1, NewVertices[Tetrahedrals[4*i + 1]]);
2115  CellTree[i]->SetVertex(2, NewVertices[Tetrahedrals[4*i + 2]]);
2116  CellTree[i]->SetVertex(3, NewVertices[Tetrahedrals[4*i + 3]]);
2117  CellTree[i]->SetClipBoard(i);
2118  ((TMacroCell *) CellTree[i])->SetSubGridID(0);
2119 // cout << i<< " TetRegionlist : " << TetRegionlist[i] << endl;
2120  CellTree[i]->SetRegionID((int)TetRegionlist[i]);
2121  }
2122 // exit(0);
2123  Domain->SetTreeInfo(CellTree, N_RootCells);
2124  // initialize iterators
2125  TDatabase::IteratorDB[It_EQ]->SetParam(Domain);
2126  TDatabase::IteratorDB[It_LE]->SetParam(Domain);
2127  TDatabase::IteratorDB[It_Finest]->SetParam(Domain);
2128  TDatabase::IteratorDB[It_Between]->SetParam(Domain);
2129  TDatabase::IteratorDB[It_OCAF]->SetParam(Domain);
2130 
2131  // search neighbours
2132  PointNeighb = new int[N_G];
2133 #ifdef _MPI
2134  if(rank==0)
2135 #endif
2136  cout<<"numberofpoints "<<N_G<<endl;
2137  memset(PointNeighb, 0, N_G *SizeOfInt);
2138 
2139  for (i=0;i<4*N_RootCells;i++)
2140  PointNeighb[Tetrahedrals[i]]++;
2141 
2142  for (i=0;i<N_G;i++)
2143  if (PointNeighb[i] > maxEpV) maxEpV = PointNeighb[i];
2144  delete [] PointNeighb;
2145 
2146 #ifdef _MPI
2147  if(rank==0)
2148 #endif
2149  cout<<"maxEpV "<< maxEpV<<endl;
2150 
2151  PointNeighb = new int[++maxEpV * N_G];
2152 
2153  memset(PointNeighb, 0, maxEpV*N_G*SizeOfInt);
2154  // every vertex contains "maxEpV" columns
2155  // for every vertex at first colomn contains the number of cells containing this vertex
2156  // at further columns we set the index of corresponding cells containing this vertex
2157  // cout<<"maxEpV*N_G "<<maxEpV*N_G<<endl;
2158 
2159  for(i=0;i<4*N_RootCells;i++)
2160  {
2161  j = Tetrahedrals[i]*maxEpV;
2162  PointNeighb[j]++;
2163  //cout<<"j + PointNeighb[j] " << j <<endl;
2164  PointNeighb[j + PointNeighb[j]] = i / 4;
2165  }
2166 
2167  // generate new faces
2168  coll=Domain->GetCollection(It_Finest, 0);
2169  N_Cells = coll->GetN_Cells();
2170  N_G = 0;
2171  for (i=0;i<N_Cells;i++)
2172  {
2173  cell = coll->GetCell(i);
2174  ShapeDesc= cell->GetShapeDesc();
2175  ShapeDesc->GetFaceVertex(TmpFV, TmpLen, MaxLen);
2176  N_Faces = cell->GetN_Faces();
2177  Tetrahedrals_loc = Tetrahedrals+4*i;
2178 
2179  for(jj=0;jj<N_Faces;jj++)
2180  if(cell->GetJoint(jj) == NULL)
2181  {
2182  N_G++;
2183  N_Points = TmpLen[jj];
2184 
2185  if(N_Points!=3)
2186  {
2187  cerr << "Only Tria faces are allowed!!! N_FVert: " << N_Points << endl;
2188  exit(-1);
2189  }
2190 
2191  //printf(" TetrameshGen Null joint \n" );
2192  a = Tetrahedrals_loc[TmpFV[jj*MaxLen]];
2193  b = Tetrahedrals_loc[TmpFV[jj*MaxLen + 1]];
2194  c = Tetrahedrals_loc[TmpFV[jj*MaxLen + 2]];
2195 
2196 // cout<<" "<< a<<" "<< b<<" "<< c<<endl;
2197 
2198  Neib[0] = -1;
2199  Neib[1] = -1;
2200  CurrNeib = 0;
2201 
2202  len1 = PointNeighb[a*maxEpV];
2203  len2 = PointNeighb[b*maxEpV];
2204  len3 = PointNeighb[c*maxEpV];
2205 
2206  // find the index of the cells containing current face with point indices a,b,c
2207  for (j=1;j<=len1;j++)
2208  {
2209  Neighb_tmp = PointNeighb[a*maxEpV + j];
2210  for (k=1;k<=len2;k++)
2211  {
2212  if (Neighb_tmp == PointNeighb[b*maxEpV + k])
2213  {
2214  for (l=1;l<=len3;l++)
2215  if (Neighb_tmp == PointNeighb[c*maxEpV + l])
2216  {
2217  Neib[CurrNeib++] = Neighb_tmp;
2218  break;
2219  }
2220  }
2221  }
2222  if (CurrNeib == 2) break;
2223  }// for (j=1;j<=len1;j++)
2224  // cout<<"CurrNeib " << CurrNeib <<endl;
2225  // cout<<"Facemarkerlist[i] : "<<Facemarkerlist[i]<<endl;
2226  if (CurrNeib == 1) // 0 for inner edges and Boundcomp+1 for Boundedge respect
2227  {
2228  CurrComp = 0; // not yet implemented fully
2229 
2230  bdcomp = Domain->GetBdPart(0)->GetBdComp(CurrComp);
2231 
2232  if(bdcomp->GetTSofXYZ(NewVertices[a]->GetX(), NewVertices[a]->GetY(),
2233  NewVertices[a]->GetY(), T[1], S[1]) ||
2234  bdcomp->GetTSofXYZ(NewVertices[b]->GetX(), NewVertices[b]->GetY(),
2235  NewVertices[b]->GetY(), T[2], S[2]) ||
2236  bdcomp->GetTSofXYZ(NewVertices[c]->GetX(), NewVertices[c]->GetY(),
2237  NewVertices[c]->GetY(), T[3], S[3]) )
2238  {
2239  cerr<<"Error: could not set parameter values"<<endl;
2240  OutPut(NewVertices[a]<<endl);
2241  OutPut(NewVertices[b]<<endl);
2242  OutPut(NewVertices[c]<<endl);
2243  exit(0);
2244  }
2245 
2246  if (CurrNeib == 2)
2247  {
2248  if(bdcomp->IsFreeBoundary())
2249  Joint = new TIsoInterfaceJoint3D(bdcomp, T, S,
2250  CellTree[Neib[0]], CellTree[Neib[1]] );
2251  else
2252  Joint = new TInterfaceJoint3D(bdcomp, T, S,
2253  CellTree[Neib[0]], CellTree[Neib[1]] );
2254  }
2255  else
2256  {
2257  if(bdcomp->IsFreeBoundary())
2258  Joint = new TIsoBoundFace(bdcomp, T, S);
2259  else
2260  Joint = new TBoundFace(bdcomp, T, S);
2261  }
2262  }
2263  else
2264  {
2265  // cout<<"Inner face"<<endl;
2266  if (CurrNeib != 2)
2267  cerr << "Error !!!!!!!! not enough neighbours!" << endl;
2268 
2269  Joint = new TJointEqN(CellTree[Neib[0]], CellTree[Neib[1]]);
2270  }
2271 
2272  // First element containing the current face
2273  // find the local index for the point 'a' on the cell
2274  for (j=0;j<4;j++)
2275  if (Tetrahedrals[4*Neib[0]+j] == a) break;
2276 
2277  // find the local index for the point 'b' on the cell
2278  for (k=0;k<4;k++)
2279  if (Tetrahedrals[4*Neib[0]+k] == b) break;
2280 
2281  // find the local index for the point 'c' on the cell
2282  for (l=0;l<4;l++)
2283  if (Tetrahedrals[4*Neib[0]+l] == c) break;
2284 
2285  l = l*100 + k*10 + j;
2286 // cout<<""<< l <<endl;
2287 
2288  switch (l) // j will contain the local index for the current face
2289  {
2290  case 210: case 21: case 102:
2291  case 120: case 12: case 201:
2292  j = 0;
2293  break;
2294  case 310: case 31: case 103:
2295  case 130: case 13: case 301:
2296  j = 1;
2297  break;
2298  case 321: case 132: case 213:
2299  case 231: case 123: case 312:
2300  j = 2;
2301  break;
2302  case 230: case 23: case 302:
2303  case 320: case 32: case 203:
2304  j = 3;
2305  break;
2306 
2307  default:
2308  Error("Unable to set the face !!!!!!!!!!!!" << endl);
2309  exit(0);
2310  }
2311  CellTree[Neib[0]]->SetJoint(j, Joint);
2312 
2313  if (Neib[1] != -1) // second element containing the current face
2314  {
2315  // find the local index for the point 'a' on the cell
2316  for (j=0;j<4;j++)
2317  if (Tetrahedrals[4*Neib[1]+j] == a) break;
2318 
2319  // find the local index for the point 'b' on the cell
2320  for (k=0;k<4;k++)
2321  if (Tetrahedrals[4*Neib[1]+k] == b) break;
2322 
2323  // find the local index for the point 'c' on the cell
2324  for (l=0;l<4;l++)
2325  if (Tetrahedrals[4*Neib[1]+l] == c) break;
2326 
2327  l = l*100 + k*10 + j;
2328 
2329 // cout<<""<< l <<endl;
2330 
2331  switch (l) // j will contain the local index for the current face
2332  {
2333  case 210: case 21: case 102:
2334  case 120: case 12: case 201:
2335  j = 0;
2336  break;
2337  case 310: case 31: case 103:
2338  case 130: case 13: case 301:
2339  j = 1;
2340  break;
2341  case 321: case 132: case 213:
2342  case 231: case 123: case 312:
2343  j = 2;
2344  break;
2345  case 230: case 23: case 302:
2346  case 320: case 32: case 203:
2347  j = 3;
2348  break;
2349 
2350  default:
2351  Error("Unable to set the face !!!!!!!!!!!!" << endl);
2352  exit(0);
2353  }
2354  CellTree[Neib[1]]->SetJoint(j, Joint);
2355  } // if (Neib[1] != -1)
2356 
2357  if (Joint->GetType() == InterfaceJoint3D ||
2358  Joint->GetType() == IsoInterfaceJoint3D)
2359  {
2360  ((TInterfaceJoint3D*)Joint)->SetMapType();
2361  ((TInterfaceJoint3D*)(Joint))->CheckOrientation();
2362  }
2363  else
2364  if (Joint->GetType() == JointEqN)
2365  Joint->SetMapType();
2366  } // if(cell->GetJoint(jj) == NULL)
2367  } // for (i=0;i<N_Cells;i++)
2368 
2369 #ifdef _MPI
2370  if(rank==0)
2371 #endif
2372  cout<<"numberoftrifaces after "<<N_G<<endl;
2373 
2374 #ifdef _MPI
2375  if(rank==0)
2376 #endif
2377  {
2378  }
2379 
2380 #ifdef _MPI
2381  if(rank!=0)
2382  {
2383  delete [] Tetrahedrals ;
2384  delete [] Coordinates;
2385  delete [] TetRegionlist;
2386  }
2387 #endif
2388 
2389  delete [] PointNeighb;
2390 }
void SetMapType()
Definition: Joint.C:107
Definition: BdPlane.h:18
void GetTreeInfo(TBaseCell **&celltree, int &N_rootcells)
get tree of cells
Definition: Domain.h:176
virtual int GetTSofXYZ(double X, double Y, double Z, double &T, double &S)=0
JointType GetType()
Definition: Joint.h:75
TBaseCell * GetCell(int i) const
return Cell with index i in Cells-array
Definition: Collection.h:50
int GetN_Faces()
return the number of faces of the cell
Definition: BaseCell.h:190
TBoundComp2D * GetBdComp(int i)
Definition: BoundPart.h:49
int SetJoint(int J_i, TJoint *J)
set the pointer to face J_i to J
Definition: BaseCell.h:168
double P0
Definition: Database.h:533
int GetFaceVertex(const int *&TmpFV, const int *&TmpLen, int &MaxLen)
Definition: ShapeDesc.h:124
contains the boundary description, the virtual cell tree and macro grid
Definition: Domain.h:36
void SetClipBoard(int value)
set value in ClipBoard
Definition: BaseCell.h:263
int SetParam(TDomain *domain)
Definition: Iterator.C:17
static TIterator ** IteratorDB
Definition: Database.h:1131
Definition: IsoBoundFace.h:18
TShapeDesc * GetShapeDesc() const
return shape descriptor of refinement descriptor
Definition: BaseCell.h:134
store cells in an array, used by cell iterators
Definition: Collection.h:18
virtual TVertex * GetVertex(int Vert_i)=0
return the pointer to vertex with number i
represent a unit of the macro grid
Definition: MacroCell.h:15
Definition: JointEqN.h:20
int GetN_Cells() const
return number of cells
Definition: Collection.h:46
Definition: ShapeDesc.h:29
TJoint * GetJoint(int J_i)
return the pointer to face with number i
Definition: BaseCell.h:175
MPI_Comm Comm
Definition: Database.h:978
Definition: Joint.h:48
TBoundPart * GetBdPart(int i)
get i-th boundary part
Definition: Domain.h:172
Definition: BoundComp3D.h:17
double GetY() const
Definition: Vertex.h:99
int GetN_Vertices()
return the number of vertices of the cell
Definition: BaseCell.h:179
Definition: BoundPart.h:21
void SetTreeInfo(TBaseCell **celltree, int N_rootcells)
set tree of cells
Definition: Domain.h:183
bool IsFreeBoundary() const
Definition: BoundComp.h:57
TCollection * GetCollection(Iterators it, int level)
produce a collection with all cells returned by iterator it
Definition: Domain.C:1982
Definition: IsoInterfaceJoint3D.h:19
information for finite element data structure
Definition: BaseCell.h:25
Definition: Vertex.h:19
static TRefDesc ** RefDescDB
Definition: Database.h:1125
void SetRegionID(int val)
set region number to this cell
Definition: BaseCell.h:347
Definition: BoundFace.h:18
double GetZ() const
Definition: Vertex.h:103
double GetX() const
Definition: Vertex.h:96
virtual int SetVertex(int Vert_i, TVertex *Vert)=0
set the pointer of vertex Vert_i to Vert
represent geometric information of the cell
Definition: GridCell.h:15
void SetAsLayerCell(int val)
set as LayerCell cell
Definition: BaseCell.h:355
Definition: BdSphere.h:18
static TParamDB * ParamDB
Definition: Database.h:1134
Definition: InterfaceJoint3D.h:18