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