Salome HOME
0e68f0214a9752fb4d5ab133d88b4213ab6afd3b
[modules/smesh.git] / src / StdMeshers / StdMeshers_Hexa_3D.cxx
1 //  Copyright (C) 2007-2008  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 //  Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 //  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 //  This library is free software; you can redistribute it and/or
7 //  modify it under the terms of the GNU Lesser General Public
8 //  License as published by the Free Software Foundation; either
9 //  version 2.1 of the License.
10 //
11 //  This library is distributed in the hope that it will be useful,
12 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
13 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 //  Lesser General Public License for more details.
15 //
16 //  You should have received a copy of the GNU Lesser General Public
17 //  License along with this library; if not, write to the Free Software
18 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 //  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 //  SMESH SMESH : implementaion of SMESH idl descriptions
23 //  File   : StdMeshers_Hexa_3D.cxx
24 //           Moved here from SMESH_Hexa_3D.cxx
25 //  Author : Paul RASCLE, EDF
26 //  Module : SMESH
27 //
28 #include "StdMeshers_Hexa_3D.hxx"
29
30 #include "StdMeshers_CompositeHexa_3D.hxx"
31 #include "StdMeshers_FaceSide.hxx"
32 #include "StdMeshers_HexaFromSkin_3D.hxx"
33 #include "StdMeshers_Penta_3D.hxx"
34 #include "StdMeshers_Prism_3D.hxx"
35 #include "StdMeshers_Quadrangle_2D.hxx"
36
37 #include "SMESH_Gen.hxx"
38 #include "SMESH_Mesh.hxx"
39 #include "SMESH_subMesh.hxx"
40 #include "SMESH_Comment.hxx"
41
42 #include "SMDS_MeshElement.hxx"
43 #include "SMDS_MeshNode.hxx"
44 #include "SMDS_FacePosition.hxx"
45 #include "SMDS_VolumeTool.hxx"
46 #include "SMDS_VolumeOfNodes.hxx"
47
48 #include <TopExp.hxx>
49 #include <TopExp_Explorer.hxx>
50 #include <TopTools_IndexedDataMapOfShapeListOfShape.hxx>
51 #include <TopTools_ListIteratorOfListOfShape.hxx>
52 #include <TopTools_ListOfShape.hxx>
53 #include <TopTools_SequenceOfShape.hxx>
54 #include <TopTools_MapOfShape.hxx>
55 #include <TopoDS.hxx>
56 #include <gp_Pnt2d.hxx>
57
58 #include "utilities.h"
59 #include "Utils_ExceptHandlers.hxx"
60
61 typedef SMESH_Comment TComm;
62
63 using namespace std;
64
65 static SMESH_ComputeErrorPtr ComputePentahedralMesh(SMESH_Mesh &,
66                                                     const TopoDS_Shape &);
67
68 static bool EvaluatePentahedralMesh(SMESH_Mesh &, const TopoDS_Shape &,
69                                     MapShapeNbElems &);
70
71 //=============================================================================
72 /*!
73  *  
74  */
75 //=============================================================================
76
77 StdMeshers_Hexa_3D::StdMeshers_Hexa_3D(int hypId, int studyId, SMESH_Gen * gen)
78   :SMESH_3D_Algo(hypId, studyId, gen)
79 {
80   MESSAGE("StdMeshers_Hexa_3D::StdMeshers_Hexa_3D");
81   _name = "Hexa_3D";
82   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);       // 1 bit /shape type
83   _requireShape = false;
84 }
85
86 //=============================================================================
87 /*!
88  *  
89  */
90 //=============================================================================
91
92 StdMeshers_Hexa_3D::~StdMeshers_Hexa_3D()
93 {
94   MESSAGE("StdMeshers_Hexa_3D::~StdMeshers_Hexa_3D");
95 }
96
97 //================================================================================
98 /*!
99  * \brief Clear fields and return the argument
100   * \param res - the value to return
101   * \retval bool - the argument value
102  */
103 //================================================================================
104
105 bool StdMeshers_Hexa_3D::ClearAndReturn(FaceQuadStruct* theQuads[6], const bool res)
106 {
107   for (int i = 0; i < 6; i++) {
108     delete theQuads[i];
109     theQuads[i] = NULL;
110   }
111   return res;
112 }
113
114
115 //=============================================================================
116 /*!
117  *  
118  */
119 //=============================================================================
120
121 bool StdMeshers_Hexa_3D::CheckHypothesis
122                          (SMESH_Mesh&                          aMesh,
123                           const TopoDS_Shape&                  aShape,
124                           SMESH_Hypothesis::Hypothesis_Status& aStatus)
125 {
126   // check nb of faces in the shape
127 /*  PAL16229
128   aStatus = SMESH_Hypothesis::HYP_BAD_GEOMETRY;
129   int nbFaces = 0;
130   for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next())
131     if ( ++nbFaces > 6 )
132       break;
133   if ( nbFaces != 6 )
134     return false;
135 */
136   aStatus = SMESH_Hypothesis::HYP_OK;
137   return true;
138 }
139
140 //=======================================================================
141 //function : isCloser
142 //purpose  : 
143 //=======================================================================
144
145 inline bool isCloser(const int i, const int j, const int nbhoriz,
146                      const FaceQuadStruct* quad, const gp_Pnt2d uv,
147                      double & minDist)
148 {
149   int ij = j * nbhoriz + i;
150   gp_Pnt2d uv2( quad->uv_grid[ij].u, quad->uv_grid[ij].v );
151   double dist = uv.SquareDistance( uv2 );
152   if ( dist < minDist ) {
153     minDist = dist;
154     return true;
155   }
156   return false;
157 }
158
159 //=======================================================================
160 //function : findIJ
161 //purpose  : return i,j of the node
162 //=======================================================================
163
164 static bool findIJ (const SMDS_MeshNode* node, const FaceQuadStruct * quad, int& I, int& J)
165 {
166   const SMDS_FacePosition* fpos =
167     static_cast<const SMDS_FacePosition*>(node->GetPosition().get());
168   if ( ! fpos ) return false;
169   gp_Pnt2d uv( fpos->GetUParameter(), fpos->GetVParameter() );
170
171   double minDist = DBL_MAX;
172   const int nbhoriz  = quad->side[0]->NbPoints();
173   const int nbvertic = quad->side[1]->NbPoints();
174   I = nbhoriz/2; J = nbvertic/2;
175   int oldI, oldJ;
176   do {
177     oldI = I; oldJ = J;
178     while ( I + 2 < nbhoriz &&  isCloser( I + 1, J, nbhoriz, quad, uv, minDist ))
179       I += 1;
180     if ( I == oldI )
181       while ( I - 1 > 0     &&  isCloser( I - 1, J, nbhoriz, quad, uv, minDist ))
182         I -= 1;
183     if ( minDist < DBL_MIN )
184       break;
185
186     while ( J + 2 < nbvertic && isCloser( I, J + 1, nbhoriz, quad, uv, minDist ))
187       J += 1;
188     if ( J == oldJ )
189       while ( J - 1 > 0      && isCloser( I, J - 1, nbhoriz, quad, uv, minDist ))
190         J -= 1;
191     if ( minDist < DBL_MIN )
192       break;
193
194   } while ( I != oldI || J != oldJ );
195
196   if ( minDist > DBL_MIN ) {
197     for (int i = 1; i < nbhoriz - 1; i++)
198       for (int j = 1; j < nbvertic - 1; j++)
199         if ( isCloser( i, j, nbhoriz, quad, uv, minDist ))
200           I = i, J = j;
201   }
202   return true;
203 }
204
205
206 //=============================================================================
207 /*!
208  * Hexahedron mesh on hexaedron like form
209  * -0.  - shape and face mesh verification
210  * -1.  - identify faces and vertices of the "cube"
211  * -2.  - Algorithm from:
212  * "Application de l'interpolation transfinie à la création de maillages
213  *  C0 ou G1 continus sur des triangles, quadrangles, tetraedres, pentaedres
214  *  et hexaedres déformés."
215  * Alain PERONNET - 8 janvier 1999
216  */
217 //=============================================================================
218
219 bool StdMeshers_Hexa_3D::Compute(SMESH_Mesh &         aMesh,
220                                  const TopoDS_Shape & aShape)// throw(SALOME_Exception)
221 {
222   // PAL14921. Enable catching std::bad_alloc and Standard_OutOfMemory outside
223   //Unexpect aCatch(SalomeException);
224   MESSAGE("StdMeshers_Hexa_3D::Compute");
225   SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
226
227   // 0.  - shape and face mesh verification
228   // 0.1 - shape must be a solid (or a shell) with 6 faces
229
230   vector < SMESH_subMesh * >meshFaces;
231   for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
232     SMESH_subMesh *aSubMesh = aMesh.GetSubMeshContaining(exp.Current());
233     ASSERT(aSubMesh);
234     meshFaces.push_back(aSubMesh);
235   }
236   if (meshFaces.size() != 6) {
237     //return error(COMPERR_BAD_SHAPE, TComm(meshFaces.size())<<" instead of 6 faces in a block");
238     static StdMeshers_CompositeHexa_3D compositeHexa(-10, 0, aMesh.GetGen());
239     if ( !compositeHexa.Compute( aMesh, aShape ))
240       return error( compositeHexa.GetComputeError() );
241     return true;
242   }
243
244   // 0.2 - is each face meshed with Quadrangle_2D? (so, with a wire of 4 edges)
245
246   // tool for working with quadratic elements
247   SMESH_MesherHelper aTool (aMesh);
248   _quadraticMesh = aTool.IsQuadraticSubMesh(aShape);
249
250   // cube structure
251   typedef struct cubeStruct
252   {
253     TopoDS_Vertex V000;
254     TopoDS_Vertex V001;
255     TopoDS_Vertex V010;
256     TopoDS_Vertex V011;
257     TopoDS_Vertex V100;
258     TopoDS_Vertex V101;
259     TopoDS_Vertex V110;
260     TopoDS_Vertex V111;
261     faceQuadStruct* quad_X0;
262     faceQuadStruct* quad_X1;
263     faceQuadStruct* quad_Y0;
264     faceQuadStruct* quad_Y1;
265     faceQuadStruct* quad_Z0;
266     faceQuadStruct* quad_Z1;
267     Point3DStruct* np; // normalised 3D coordinates
268   } CubeStruct;
269
270   CubeStruct aCube;
271
272   // bounding faces
273   FaceQuadStruct* aQuads[6];
274   for (int i = 0; i < 6; i++)
275     aQuads[i] = 0;
276
277   for (int i = 0; i < 6; i++)
278   {
279     TopoDS_Shape aFace = meshFaces[i]->GetSubShape();
280     SMESH_Algo *algo = _gen->GetAlgo(aMesh, aFace);
281     string algoName = algo->GetName();
282     bool isAllQuad = false;
283     if (algoName == "Quadrangle_2D") {
284       SMESHDS_SubMesh * sm = meshDS->MeshElements( aFace );
285       if ( sm ) {
286         isAllQuad = true;
287         SMDS_ElemIteratorPtr eIt = sm->GetElements();
288         while ( isAllQuad && eIt->more() ) {
289           const SMDS_MeshElement* elem =  eIt->next();
290           isAllQuad = ( elem->NbNodes()==4 ||(_quadraticMesh && elem->NbNodes()==8) );
291         }
292       }
293     }
294     if ( ! isAllQuad ) {
295       SMESH_ComputeErrorPtr err = ComputePentahedralMesh(aMesh, aShape);
296       return ClearAndReturn( aQuads, error(err));
297     }
298     StdMeshers_Quadrangle_2D *quadAlgo =
299       dynamic_cast < StdMeshers_Quadrangle_2D * >(algo);
300     ASSERT(quadAlgo);
301     try {
302       aQuads[i] = quadAlgo->CheckAnd2Dcompute(aMesh, aFace, _quadraticMesh);
303       if(!aQuads[i]) {
304         return error( quadAlgo->GetComputeError());
305       }
306     }
307     catch(SALOME_Exception & S_ex) {
308       return ClearAndReturn( aQuads, error(COMPERR_SLM_EXCEPTION,TComm(S_ex.what()) <<
309                                            " Raised by StdMeshers_Quadrangle_2D "
310                                            " on face #" << meshDS->ShapeToIndex( aFace )));
311     }
312
313     // 0.2.1 - number of points on the opposite edges must be the same
314     if (aQuads[i]->side[0]->NbPoints() != aQuads[i]->side[2]->NbPoints() ||
315         aQuads[i]->side[1]->NbPoints() != aQuads[i]->side[3]->NbPoints()
316         /*aQuads[i]->side[0]->NbEdges() != 1 ||
317         aQuads[i]->side[1]->NbEdges() != 1 ||
318         aQuads[i]->side[2]->NbEdges() != 1 ||
319         aQuads[i]->side[3]->NbEdges() != 1*/) {
320       MESSAGE("different number of points on the opposite edges of face " << i);
321       // Try to go into penta algorithm 'cause it has been improved.
322       SMESH_ComputeErrorPtr err = ComputePentahedralMesh(aMesh, aShape);
323       return ClearAndReturn( aQuads, error(err));
324     }
325   }
326
327   // 1.  - identify faces and vertices of the "cube"
328   // 1.1 - ancestor maps vertex->edges in the cube
329
330 //   TopTools_IndexedDataMapOfShapeListOfShape MS;
331 //   TopExp::MapShapesAndAncestors(aShape, TopAbs_VERTEX, TopAbs_EDGE, MS);
332
333   // 1.2 - first face is choosen as face Y=0 of the unit cube
334
335   const TopoDS_Shape & aFace = meshFaces[0]->GetSubShape();
336   //const TopoDS_Face & F = TopoDS::Face(aFace);
337
338   // 1.3 - identify the 4 vertices of the face Y=0: V000, V100, V101, V001
339
340   aCube.V000 = aQuads[0]->side[0]->FirstVertex(); // will be (0,0,0) on the unit cube
341   aCube.V100 = aQuads[0]->side[0]->LastVertex();  // will be (1,0,0) on the unit cube
342   aCube.V001 = aQuads[0]->side[2]->FirstVertex(); // will be (0,0,1) on the unit cube
343   aCube.V101 = aQuads[0]->side[2]->LastVertex();  // will be (1,0,1) on the unit cube
344
345   TopTools_IndexedMapOfShape MV0;
346   TopExp::MapShapes(aFace, TopAbs_VERTEX, MV0);
347
348   aCube.V010 = OppositeVertex( aCube.V000, MV0, aQuads);
349   aCube.V110 = OppositeVertex( aCube.V100, MV0, aQuads);
350   aCube.V011 = OppositeVertex( aCube.V001, MV0, aQuads);
351   aCube.V111 = OppositeVertex( aCube.V101, MV0, aQuads);
352
353   // 1.6 - find remaining faces given 4 vertices
354
355   int _indY0 = 0;
356   int _indY1 = GetFaceIndex(aMesh, aShape, meshFaces,
357                             aCube.V010, aCube.V011, aCube.V110, aCube.V111);
358   int _indZ0 = GetFaceIndex(aMesh, aShape, meshFaces,
359                             aCube.V000, aCube.V010, aCube.V100, aCube.V110);
360   int _indZ1 = GetFaceIndex(aMesh, aShape, meshFaces,
361                             aCube.V001, aCube.V011, aCube.V101, aCube.V111);
362   int _indX0 = GetFaceIndex(aMesh, aShape, meshFaces,
363                             aCube.V000, aCube.V001, aCube.V010, aCube.V011);
364   int _indX1 = GetFaceIndex(aMesh, aShape, meshFaces,
365                             aCube.V100, aCube.V101, aCube.V110, aCube.V111);
366
367   // IPAL21120: SIGSEGV on Meshing attached Compound with Automatic Hexadralization
368   if ( _indY1 < 1 || _indZ0 < 1 || _indZ1 < 1 || _indX0 < 1 || _indX1 < 1 )
369     return error(COMPERR_BAD_SHAPE);
370
371   aCube.quad_Y0 = aQuads[_indY0];
372   aCube.quad_Y1 = aQuads[_indY1];
373   aCube.quad_Z0 = aQuads[_indZ0];
374   aCube.quad_Z1 = aQuads[_indZ1];
375   aCube.quad_X0 = aQuads[_indX0];
376   aCube.quad_X1 = aQuads[_indX1];
377
378   // 1.7 - get convertion coefs from face 2D normalized to 3D normalized
379
380   Conv2DStruct cx0;                     // for face X=0
381   Conv2DStruct cx1;                     // for face X=1
382   Conv2DStruct cy0;
383   Conv2DStruct cy1;
384   Conv2DStruct cz0;
385   Conv2DStruct cz1;
386
387   GetConv2DCoefs(*aCube.quad_X0, meshFaces[_indX0]->GetSubShape(),
388                  aCube.V000, aCube.V010, aCube.V011, aCube.V001, cx0);
389   GetConv2DCoefs(*aCube.quad_X1, meshFaces[_indX1]->GetSubShape(),
390                  aCube.V100, aCube.V110, aCube.V111, aCube.V101, cx1);
391   GetConv2DCoefs(*aCube.quad_Y0, meshFaces[_indY0]->GetSubShape(),
392                  aCube.V000, aCube.V100, aCube.V101, aCube.V001, cy0);
393   GetConv2DCoefs(*aCube.quad_Y1, meshFaces[_indY1]->GetSubShape(),
394                  aCube.V010, aCube.V110, aCube.V111, aCube.V011, cy1);
395   GetConv2DCoefs(*aCube.quad_Z0, meshFaces[_indZ0]->GetSubShape(),
396                  aCube.V000, aCube.V100, aCube.V110, aCube.V010, cz0);
397   GetConv2DCoefs(*aCube.quad_Z1, meshFaces[_indZ1]->GetSubShape(),
398                  aCube.V001, aCube.V101, aCube.V111, aCube.V011, cz1);
399
400   // 1.8 - create a 3D structure for normalized values
401   
402   int nbx = aCube.quad_Z0->side[0]->NbPoints();
403   if (cz0.a1 == 0.) nbx = aCube.quad_Z0->side[1]->NbPoints();
404  
405   int nby = aCube.quad_X0->side[0]->NbPoints();
406   if (cx0.a1 == 0.) nby = aCube.quad_X0->side[1]->NbPoints();
407  
408   int nbz = aCube.quad_Y0->side[0]->NbPoints();
409   if (cy0.a1 != 0.) nbz = aCube.quad_Y0->side[1]->NbPoints();
410
411   int i1, j1, nbxyz = nbx * nby * nbz;
412   Point3DStruct *np = new Point3DStruct[nbxyz];
413
414   // 1.9 - store node indexes of faces
415
416   {
417     const TopoDS_Face & F = TopoDS::Face(meshFaces[_indX0]->GetSubShape());
418
419     faceQuadStruct *quad = aCube.quad_X0;
420     int i = 0;                          // j = x/face , k = y/face
421     int nbdown = quad->side[0]->NbPoints();
422     int nbright = quad->side[1]->NbPoints();
423
424     SMDS_NodeIteratorPtr itf= aMesh.GetSubMesh(F)->GetSubMeshDS()->GetNodes();
425                         
426     while(itf->more()) {
427       const SMDS_MeshNode * node = itf->next();
428       if(aTool.IsMedium(node))
429         continue;
430       if ( !findIJ( node, quad, i1, j1 ))
431         return ClearAndReturn( aQuads, false );
432       int ij1 = j1 * nbdown + i1;
433       quad->uv_grid[ij1].node = node;
434     }
435
436     for (int i1 = 0; i1 < nbdown; i1++)
437       for (int j1 = 0; j1 < nbright; j1++) {
438         int ij1 = j1 * nbdown + i1;
439         int j = cx0.ia * i1 + cx0.ib * j1 + cx0.ic;     // j = x/face
440         int k = cx0.ja * i1 + cx0.jb * j1 + cx0.jc;     // k = y/face
441         int ijk = k * nbx * nby + j * nbx + i;
442         //MESSAGE(" "<<ij1<<" "<<i<<" "<<j<<" "<<ijk);
443         np[ijk].node = quad->uv_grid[ij1].node;
444         //SCRUTE(np[ijk].nodeId);
445       }
446   }
447
448   {
449     const TopoDS_Face & F = TopoDS::Face(meshFaces[_indX1]->GetSubShape());
450
451     SMDS_NodeIteratorPtr itf= aMesh.GetSubMesh(F)->GetSubMeshDS()->GetNodes();
452
453     faceQuadStruct *quad = aCube.quad_X1;
454     int i = nbx - 1;            // j = x/face , k = y/face
455     int nbdown = quad->side[0]->NbPoints();
456     int nbright = quad->side[1]->NbPoints();
457
458     while(itf->more()) {
459       const SMDS_MeshNode * node = itf->next();
460       if(aTool.IsMedium(node))
461         continue;
462       if ( !findIJ( node, quad, i1, j1 ))
463         return ClearAndReturn( aQuads, false );
464       int ij1 = j1 * nbdown + i1;
465       quad->uv_grid[ij1].node = node;
466     }
467
468     for (int i1 = 0; i1 < nbdown; i1++)
469       for (int j1 = 0; j1 < nbright; j1++) {
470         int ij1 = j1 * nbdown + i1;
471         int j = cx1.ia * i1 + cx1.ib * j1 + cx1.ic;     // j = x/face
472         int k = cx1.ja * i1 + cx1.jb * j1 + cx1.jc;     // k = y/face
473         int ijk = k * nbx * nby + j * nbx + i;
474         //MESSAGE(" "<<ij1<<" "<<i<<" "<<j<<" "<<ijk);
475         np[ijk].node = quad->uv_grid[ij1].node;
476         //SCRUTE(np[ijk].nodeId);
477       }
478   }
479
480   {
481     const TopoDS_Face & F = TopoDS::Face(meshFaces[_indY0]->GetSubShape());
482
483     SMDS_NodeIteratorPtr itf= aMesh.GetSubMesh(F)->GetSubMeshDS()->GetNodes();
484
485     faceQuadStruct *quad = aCube.quad_Y0;
486     int j = 0;                          // i = x/face , k = y/face
487     int nbdown = quad->side[0]->NbPoints();
488     int nbright = quad->side[1]->NbPoints();
489
490     while(itf->more()) {
491       const SMDS_MeshNode * node = itf->next();
492       if(aTool.IsMedium(node))
493         continue;
494       if ( !findIJ( node, quad, i1, j1 ))
495         return ClearAndReturn( aQuads, false );
496       int ij1 = j1 * nbdown + i1;
497       quad->uv_grid[ij1].node = node;
498     }
499
500     for (int i1 = 0; i1 < nbdown; i1++)
501       for (int j1 = 0; j1 < nbright; j1++) {
502         int ij1 = j1 * nbdown + i1;
503         int i = cy0.ia * i1 + cy0.ib * j1 + cy0.ic;     // i = x/face
504         int k = cy0.ja * i1 + cy0.jb * j1 + cy0.jc;     // k = y/face
505         int ijk = k * nbx * nby + j * nbx + i;
506         //MESSAGE(" "<<ij1<<" "<<i<<" "<<j<<" "<<ijk);
507         np[ijk].node = quad->uv_grid[ij1].node;
508         //SCRUTE(np[ijk].nodeId);
509       }
510   }
511
512   {
513     const TopoDS_Face & F = TopoDS::Face(meshFaces[_indY1]->GetSubShape());
514
515     SMDS_NodeIteratorPtr itf= aMesh.GetSubMesh(F)->GetSubMeshDS()->GetNodes();
516
517     faceQuadStruct *quad = aCube.quad_Y1;
518     int j = nby - 1;            // i = x/face , k = y/face
519     int nbdown = quad->side[0]->NbPoints();
520     int nbright = quad->side[1]->NbPoints();
521
522     while(itf->more()) {
523       const SMDS_MeshNode * node = itf->next();
524       if(aTool.IsMedium(node))
525         continue;
526       if ( !findIJ( node, quad, i1, j1 ))
527         return ClearAndReturn( aQuads, false );
528       int ij1 = j1 * nbdown + i1;
529       quad->uv_grid[ij1].node = node;
530     }
531
532     for (int i1 = 0; i1 < nbdown; i1++)
533       for (int j1 = 0; j1 < nbright; j1++) {
534         int ij1 = j1 * nbdown + i1;
535         int i = cy1.ia * i1 + cy1.ib * j1 + cy1.ic;     // i = x/face
536         int k = cy1.ja * i1 + cy1.jb * j1 + cy1.jc;     // k = y/face
537         int ijk = k * nbx * nby + j * nbx + i;
538         //MESSAGE(" "<<ij1<<" "<<i<<" "<<j<<" "<<ijk);
539         np[ijk].node = quad->uv_grid[ij1].node;
540         //SCRUTE(np[ijk].nodeId);
541       }
542   }
543
544   {
545     const TopoDS_Face & F = TopoDS::Face(meshFaces[_indZ0]->GetSubShape());
546
547     SMDS_NodeIteratorPtr itf= aMesh.GetSubMesh(F)->GetSubMeshDS()->GetNodes();
548
549     faceQuadStruct *quad = aCube.quad_Z0;
550     int k = 0;                          // i = x/face , j = y/face
551     int nbdown = quad->side[0]->NbPoints();
552     int nbright = quad->side[1]->NbPoints();
553
554     while(itf->more()) {
555       const SMDS_MeshNode * node = itf->next();
556       if(aTool.IsMedium(node))
557         continue;
558       if ( !findIJ( node, quad, i1, j1 ))
559         return ClearAndReturn( aQuads, false );
560       int ij1 = j1 * nbdown + i1;
561       quad->uv_grid[ij1].node = node;
562     }
563
564     for (int i1 = 0; i1 < nbdown; i1++)
565       for (int j1 = 0; j1 < nbright; j1++) {
566         int ij1 = j1 * nbdown + i1;
567         int i = cz0.ia * i1 + cz0.ib * j1 + cz0.ic;     // i = x/face
568         int j = cz0.ja * i1 + cz0.jb * j1 + cz0.jc;     // j = y/face
569         int ijk = k * nbx * nby + j * nbx + i;
570         //MESSAGE(" "<<ij1<<" "<<i<<" "<<j<<" "<<ijk);
571         np[ijk].node = quad->uv_grid[ij1].node;
572         //SCRUTE(np[ijk].nodeId);
573       }
574   }
575
576   {
577     const TopoDS_Face & F = TopoDS::Face(meshFaces[_indZ1]->GetSubShape());
578
579     SMDS_NodeIteratorPtr itf= aMesh.GetSubMesh(F)->GetSubMeshDS()->GetNodes();
580
581     faceQuadStruct *quad = aCube.quad_Z1;
582     int k = nbz - 1;            // i = x/face , j = y/face
583     int nbdown = quad->side[0]->NbPoints();
584     int nbright = quad->side[1]->NbPoints();
585     
586     while(itf->more()) {
587       const SMDS_MeshNode * node = itf->next();
588       if(aTool.IsMedium(node))
589         continue;
590       if ( !findIJ( node, quad, i1, j1 ))
591         return ClearAndReturn( aQuads, false );
592       int ij1 = j1 * nbdown + i1;
593       quad->uv_grid[ij1].node = node;
594     }
595
596     for (int i1 = 0; i1 < nbdown; i1++)
597       for (int j1 = 0; j1 < nbright; j1++) {
598         int ij1 = j1 * nbdown + i1;
599         int i = cz1.ia * i1 + cz1.ib * j1 + cz1.ic;     // i = x/face
600         int j = cz1.ja * i1 + cz1.jb * j1 + cz1.jc;     // j = y/face
601         int ijk = k * nbx * nby + j * nbx + i;
602         //MESSAGE(" "<<ij1<<" "<<i<<" "<<j<<" "<<ijk);
603         np[ijk].node = quad->uv_grid[ij1].node;
604         //SCRUTE(np[ijk].nodeId);
605       }
606   }
607
608   // 2.0 - for each node of the cube:
609   //       - get the 8 points 3D = 8 vertices of the cube
610   //       - get the 12 points 3D on the 12 edges of the cube
611   //       - get the 6 points 3D on the 6 faces with their ID
612   //       - compute the point 3D
613   //       - store the point 3D in SMESHDS, store its ID in 3D structure
614
615   int shapeID = meshDS->ShapeToIndex( aShape );
616
617   Pt3 p000, p001, p010, p011, p100, p101, p110, p111;
618   Pt3 px00, px01, px10, px11;
619   Pt3 p0y0, p0y1, p1y0, p1y1;
620   Pt3 p00z, p01z, p10z, p11z;
621   Pt3 pxy0, pxy1, px0z, px1z, p0yz, p1yz;
622
623   GetPoint(p000, 0, 0, 0, nbx, nby, nbz, np, meshDS);
624   GetPoint(p001, 0, 0, nbz - 1, nbx, nby, nbz, np, meshDS);
625   GetPoint(p010, 0, nby - 1, 0, nbx, nby, nbz, np, meshDS);
626   GetPoint(p011, 0, nby - 1, nbz - 1, nbx, nby, nbz, np, meshDS);
627   GetPoint(p100, nbx - 1, 0, 0, nbx, nby, nbz, np, meshDS);
628   GetPoint(p101, nbx - 1, 0, nbz - 1, nbx, nby, nbz, np, meshDS);
629   GetPoint(p110, nbx - 1, nby - 1, 0, nbx, nby, nbz, np, meshDS);
630   GetPoint(p111, nbx - 1, nby - 1, nbz - 1, nbx, nby, nbz, np, meshDS);
631
632   for (int i = 1; i < nbx - 1; i++) {
633     for (int j = 1; j < nby - 1; j++) {
634       for (int k = 1; k < nbz - 1; k++) {
635         // *** seulement maillage regulier
636         // 12 points on edges
637         GetPoint(px00, i, 0, 0, nbx, nby, nbz, np, meshDS);
638         GetPoint(px01, i, 0, nbz - 1, nbx, nby, nbz, np, meshDS);
639         GetPoint(px10, i, nby - 1, 0, nbx, nby, nbz, np, meshDS);
640         GetPoint(px11, i, nby - 1, nbz - 1, nbx, nby, nbz, np, meshDS);
641
642         GetPoint(p0y0, 0, j, 0, nbx, nby, nbz, np, meshDS);
643         GetPoint(p0y1, 0, j, nbz - 1, nbx, nby, nbz, np, meshDS);
644         GetPoint(p1y0, nbx - 1, j, 0, nbx, nby, nbz, np, meshDS);
645         GetPoint(p1y1, nbx - 1, j, nbz - 1, nbx, nby, nbz, np, meshDS);
646
647         GetPoint(p00z, 0, 0, k, nbx, nby, nbz, np, meshDS);
648         GetPoint(p01z, 0, nby - 1, k, nbx, nby, nbz, np, meshDS);
649         GetPoint(p10z, nbx - 1, 0, k, nbx, nby, nbz, np, meshDS);
650         GetPoint(p11z, nbx - 1, nby - 1, k, nbx, nby, nbz, np, meshDS);
651
652         // 12 points on faces
653         GetPoint(pxy0, i, j, 0, nbx, nby, nbz, np, meshDS);
654         GetPoint(pxy1, i, j, nbz - 1, nbx, nby, nbz, np, meshDS);
655         GetPoint(px0z, i, 0, k, nbx, nby, nbz, np, meshDS);
656         GetPoint(px1z, i, nby - 1, k, nbx, nby, nbz, np, meshDS);
657         GetPoint(p0yz, 0, j, k, nbx, nby, nbz, np, meshDS);
658         GetPoint(p1yz, nbx - 1, j, k, nbx, nby, nbz, np, meshDS);
659
660         int ijk = k * nbx * nby + j * nbx + i;
661         double x = double (i) / double (nbx - 1);       // *** seulement
662         double y = double (j) / double (nby - 1);       // *** maillage
663         double z = double (k) / double (nbz - 1);       // *** regulier
664
665         Pt3 X;
666         for (int i = 0; i < 3; i++) {
667           X[i] = (1 - x) * p0yz[i] + x * p1yz[i]
668                  + (1 - y) * px0z[i] + y * px1z[i]
669                  + (1 - z) * pxy0[i] + z * pxy1[i]
670                  - (1 - x) * ((1 - y) * p00z[i] + y * p01z[i])
671                  - x * ((1 - y) * p10z[i] + y * p11z[i])
672                  - (1 - y) * ((1 - z) * px00[i] + z * px01[i])
673                  - y * ((1 - z) * px10[i] + z * px11[i])
674                  - (1 - z) * ((1 - x) * p0y0[i] + x * p1y0[i])
675                  - z * ((1 - x) * p0y1[i] + x * p1y1[i])
676                  + (1 - x) * ((1 - y) * ((1 - z) * p000[i] + z * p001[i])
677                  + y * ((1 - z) * p010[i] + z * p011[i]))
678                  + x * ((1 - y) * ((1 - z) * p100[i] + z * p101[i])
679                  + y * ((1 - z) * p110[i] + z * p111[i]));
680         }
681
682         SMDS_MeshNode * node = meshDS->AddNode(X[0], X[1], X[2]);
683         np[ijk].node = node;
684         meshDS->SetNodeInVolume(node, shapeID);
685       }
686     }
687   }
688
689   // find orientation of furute volumes according to MED convention
690   vector< bool > forward( nbx * nby );
691   SMDS_VolumeTool vTool;
692   for (int i = 0; i < nbx - 1; i++) {
693     for (int j = 0; j < nby - 1; j++) {
694       int n1 = j * nbx + i;
695       int n2 = j * nbx + i + 1;
696       int n3 = (j + 1) * nbx + i + 1;
697       int n4 = (j + 1) * nbx + i;
698       int n5 = nbx * nby + j * nbx + i;
699       int n6 = nbx * nby + j * nbx + i + 1;
700       int n7 = nbx * nby + (j + 1) * nbx + i + 1;
701       int n8 = nbx * nby + (j + 1) * nbx + i;
702
703       SMDS_VolumeOfNodes tmpVol (np[n1].node,np[n2].node,np[n3].node,np[n4].node,
704                                  np[n5].node,np[n6].node,np[n7].node,np[n8].node);
705       vTool.Set( &tmpVol );
706       forward[ n1 ] = vTool.IsForward();
707     }
708   }
709
710   //2.1 - for each node of the cube (less 3 *1 Faces):
711   //      - store hexahedron in SMESHDS
712   MESSAGE("Storing hexahedron into the DS");
713   for (int i = 0; i < nbx - 1; i++) {
714     for (int j = 0; j < nby - 1; j++) {
715       bool isForw = forward.at( j * nbx + i );
716       for (int k = 0; k < nbz - 1; k++) {
717         int n1 = k * nbx * nby + j * nbx + i;
718         int n2 = k * nbx * nby + j * nbx + i + 1;
719         int n3 = k * nbx * nby + (j + 1) * nbx + i + 1;
720         int n4 = k * nbx * nby + (j + 1) * nbx + i;
721         int n5 = (k + 1) * nbx * nby + j * nbx + i;
722         int n6 = (k + 1) * nbx * nby + j * nbx + i + 1;
723         int n7 = (k + 1) * nbx * nby + (j + 1) * nbx + i + 1;
724         int n8 = (k + 1) * nbx * nby + (j + 1) * nbx + i;
725
726         SMDS_MeshVolume * elt;
727         if ( isForw ) {
728           elt = aTool.AddVolume(np[n1].node, np[n2].node,
729                                 np[n3].node, np[n4].node,
730                                 np[n5].node, np[n6].node,
731                                 np[n7].node, np[n8].node);
732         }
733         else {
734           elt = aTool.AddVolume(np[n1].node, np[n4].node,
735                                 np[n3].node, np[n2].node,
736                                 np[n5].node, np[n8].node,
737                                 np[n7].node, np[n6].node);
738         }
739         
740         meshDS->SetMeshElementOnShape(elt, shapeID);
741       }
742     }
743   }
744   if ( np ) delete [] np;
745   return ClearAndReturn( aQuads, true );
746 }
747
748
749 //=============================================================================
750 /*!
751  *  Evaluate
752  */
753 //=============================================================================
754
755 bool StdMeshers_Hexa_3D::Evaluate(SMESH_Mesh & aMesh,
756                                   const TopoDS_Shape & aShape,
757                                   MapShapeNbElems& aResMap)
758 {
759   vector < SMESH_subMesh * >meshFaces;
760   TopTools_SequenceOfShape aFaces;
761   for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
762     aFaces.Append(exp.Current());
763     SMESH_subMesh *aSubMesh = aMesh.GetSubMeshContaining(exp.Current());
764     ASSERT(aSubMesh);
765     meshFaces.push_back(aSubMesh);
766   }
767   if (meshFaces.size() != 6) {
768     //return error(COMPERR_BAD_SHAPE, TComm(meshFaces.size())<<" instead of 6 faces in a block");
769     static StdMeshers_CompositeHexa_3D compositeHexa(-10, 0, aMesh.GetGen());
770     return compositeHexa.Evaluate(aMesh, aShape, aResMap);
771   }
772   
773   int i = 0;
774   for(; i<6; i++) {
775     //TopoDS_Shape aFace = meshFaces[i]->GetSubShape();
776     TopoDS_Shape aFace = aFaces.Value(i+1);
777     SMESH_Algo *algo = _gen->GetAlgo(aMesh, aFace);
778     if( !algo ) {
779       std::vector<int> aResVec(SMDSEntity_Last);
780       for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
781       SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
782       aResMap.insert(std::make_pair(sm,aResVec));
783       SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
784       smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
785       return false;
786     }
787     string algoName = algo->GetName();
788     bool isAllQuad = false;
789     if (algoName == "Quadrangle_2D") {
790       MapShapeNbElemsItr anIt = aResMap.find(meshFaces[i]);
791       if( anIt == aResMap.end() ) continue;
792       std::vector<int> aVec = (*anIt).second;
793       int nbtri = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
794       if( nbtri == 0 )
795         isAllQuad = true;
796     }
797     if ( ! isAllQuad ) {
798       return EvaluatePentahedralMesh(aMesh, aShape, aResMap);
799     }
800   }
801   
802   // find number of 1d elems for 1 face
803   int nb1d = 0;
804   TopTools_MapOfShape Edges1;
805   bool IsQuadratic = false;
806   bool IsFirst = true;
807   for (TopExp_Explorer exp(aFaces.Value(1), TopAbs_EDGE); exp.More(); exp.Next()) {
808     Edges1.Add(exp.Current());
809     SMESH_subMesh *sm = aMesh.GetSubMesh(exp.Current());
810     if( sm ) {
811       MapShapeNbElemsItr anIt = aResMap.find(sm);
812       if( anIt == aResMap.end() ) continue;
813       std::vector<int> aVec = (*anIt).second;
814       nb1d += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
815       if(IsFirst) {
816         IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
817         IsFirst = false;
818       }
819     }
820   }
821   // find face opposite to 1 face
822   int OppNum = 0;
823   for(i=2; i<=6; i++) {
824     bool IsOpposite = true;
825     for(TopExp_Explorer exp(aFaces.Value(i), TopAbs_EDGE); exp.More(); exp.Next()) {
826       if( Edges1.Contains(exp.Current()) ) {
827         IsOpposite = false;
828         break;
829       }
830     }
831     if(IsOpposite) {
832       OppNum = i;
833       break;
834     }
835   }
836   // find number of 2d elems on side faces
837   int nb2d = 0;
838   for(i=2; i<=6; i++) {
839     if( i == OppNum ) continue;
840     MapShapeNbElemsItr anIt = aResMap.find( meshFaces[i-1] );
841     if( anIt == aResMap.end() ) continue;
842     std::vector<int> aVec = (*anIt).second;
843     nb2d += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
844   }
845   
846   MapShapeNbElemsItr anIt = aResMap.find( meshFaces[0] );
847   std::vector<int> aVec = (*anIt).second;
848   int nb2d_face0 = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
849   int nb0d_face0 = aVec[SMDSEntity_Node];
850
851   std::vector<int> aResVec(SMDSEntity_Last);
852   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
853   if(IsQuadratic) {
854     aResVec[SMDSEntity_Quad_Hexa] = nb2d_face0 * ( nb2d/nb1d );
855     int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
856     aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
857   }
858   else {
859     aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
860     aResVec[SMDSEntity_Hexa] = nb2d_face0 * ( nb2d/nb1d );
861   }
862   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
863   aResMap.insert(std::make_pair(sm,aResVec));
864
865   return true;
866 }
867
868 //================================================================================
869 /*!
870  * \brief Computes hexahedral mesh from 2D mesh of block
871  */
872 //================================================================================
873
874 bool StdMeshers_Hexa_3D::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
875 {
876   static StdMeshers_HexaFromSkin_3D * algo = 0;
877   if ( !algo ) {
878     SMESH_Gen* gen = aMesh.GetGen();
879     algo = new StdMeshers_HexaFromSkin_3D( gen->GetANewId(), 0, gen );
880   }
881   algo->InitComputeError();
882   algo->Compute( aMesh, aHelper );
883   return error( algo->GetComputeError());
884 }
885
886 //=============================================================================
887 /*!
888  *
889  */
890 //=============================================================================
891
892 void StdMeshers_Hexa_3D::GetPoint(Pt3 p, int i, int j, int k, int nbx, int nby, int nbz,
893                                   Point3DStruct * np, const SMESHDS_Mesh * meshDS)
894 {
895         int ijk = k * nbx * nby + j * nbx + i;
896         const SMDS_MeshNode * node = np[ijk].node;
897         p[0] = node->X();
898         p[1] = node->Y();
899         p[2] = node->Z();
900         //MESSAGE(" "<<i<<" "<<j<<" "<<k<<" "<<p[0]<<" "<<p[1]<<" "<<p[2]);
901 }
902
903 //=============================================================================
904 /*!
905  *  
906  */
907 //=============================================================================
908
909 int StdMeshers_Hexa_3D::GetFaceIndex(SMESH_Mesh & aMesh,
910         const TopoDS_Shape & aShape,
911         const vector < SMESH_subMesh * >&meshFaces,
912         const TopoDS_Vertex & V0,
913         const TopoDS_Vertex & V1,
914         const TopoDS_Vertex & V2, const TopoDS_Vertex & V3)
915 {
916         //MESSAGE("StdMeshers_Hexa_3D::GetFaceIndex");
917         int faceIndex = -1;
918         for (int i = 1; i < 6; i++)
919         {
920                 const TopoDS_Shape & aFace = meshFaces[i]->GetSubShape();
921                 //const TopoDS_Face& F = TopoDS::Face(aFace);
922                 TopTools_IndexedMapOfShape M;
923                 TopExp::MapShapes(aFace, TopAbs_VERTEX, M);
924                 bool verticesInShape = false;
925                 if (M.Contains(V0))
926                         if (M.Contains(V1))
927                                 if (M.Contains(V2))
928                                         if (M.Contains(V3))
929                                                 verticesInShape = true;
930                 if (verticesInShape)
931                 {
932                         faceIndex = i;
933                         break;
934                 }
935         }
936         //IPAL21120 ASSERT(faceIndex > 0);
937         //SCRUTE(faceIndex);
938         return faceIndex;
939 }
940
941 //=============================================================================
942 /*!
943  *  
944  */
945 //=============================================================================
946
947 TopoDS_Edge
948         StdMeshers_Hexa_3D::EdgeNotInFace(SMESH_Mesh & aMesh,
949         const TopoDS_Shape & aShape,
950         const TopoDS_Face & aFace,
951         const TopoDS_Vertex & aVertex,
952         const TopTools_IndexedDataMapOfShapeListOfShape & MS)
953 {
954         //MESSAGE("StdMeshers_Hexa_3D::EdgeNotInFace");
955         TopTools_IndexedDataMapOfShapeListOfShape MF;
956         TopExp::MapShapesAndAncestors(aFace, TopAbs_VERTEX, TopAbs_EDGE, MF);
957         const TopTools_ListOfShape & ancestorsInSolid = MS.FindFromKey(aVertex);
958         const TopTools_ListOfShape & ancestorsInFace = MF.FindFromKey(aVertex);
959 //      SCRUTE(ancestorsInSolid.Extent());
960 //      SCRUTE(ancestorsInFace.Extent());
961         ASSERT(ancestorsInSolid.Extent() == 6); // 6 (edges doublees)
962         ASSERT(ancestorsInFace.Extent() == 2);
963
964         TopoDS_Edge E;
965         E.Nullify();
966         TopTools_ListIteratorOfListOfShape its(ancestorsInSolid);
967         for (; its.More(); its.Next())
968         {
969                 TopoDS_Shape ancestor = its.Value();
970                 TopTools_ListIteratorOfListOfShape itf(ancestorsInFace);
971                 bool isInFace = false;
972                 for (; itf.More(); itf.Next())
973                 {
974                         TopoDS_Shape ancestorInFace = itf.Value();
975                         if (ancestorInFace.IsSame(ancestor))
976                         {
977                                 isInFace = true;
978                                 break;
979                         }
980                 }
981                 if (!isInFace)
982                 {
983                         E = TopoDS::Edge(ancestor);
984                         break;
985                 }
986         }
987         return E;
988 }
989
990 //=============================================================================
991 /*!
992  *  
993  */
994 //=============================================================================
995
996 void StdMeshers_Hexa_3D::GetConv2DCoefs(const faceQuadStruct & quad,
997         const TopoDS_Shape & aShape,
998         const TopoDS_Vertex & V0,
999         const TopoDS_Vertex & V1,
1000         const TopoDS_Vertex & V2, const TopoDS_Vertex & V3, Conv2DStruct & conv)
1001 {
1002 //      MESSAGE("StdMeshers_Hexa_3D::GetConv2DCoefs");
1003 //      const TopoDS_Face & F = TopoDS::Face(aShape);
1004 //      TopoDS_Edge E = quad.edge[0];
1005 //      double f, l;
1006 //      Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
1007 //      TopoDS_Vertex VFirst, VLast;
1008 //      TopExp::Vertices(E, VFirst, VLast);     // corresponds to f and l
1009 //      bool isForward = (((l - f) * (quad.last[0] - quad.first[0])) > 0);
1010   TopoDS_Vertex VA, VB;
1011 //      if (isForward)
1012 //      {
1013 //              VA = VFirst;
1014 //              VB = VLast;
1015 //      }
1016 //      else
1017 //      {
1018 //              VA = VLast;
1019 //              VB = VFirst;
1020 //      }
1021   VA = quad.side[0]->FirstVertex();
1022   VB = quad.side[0]->LastVertex();
1023   
1024         int a1, b1, c1, a2, b2, c2;
1025         if (VA.IsSame(V0))
1026                 if (VB.IsSame(V1))
1027                 {
1028                         a1 = 1;
1029                         b1 = 0;
1030                         c1 = 0;                         // x
1031                         a2 = 0;
1032                         b2 = 1;
1033                         c2 = 0;                         // y
1034                 }
1035                 else
1036                 {
1037                         ASSERT(VB.IsSame(V3));
1038                         a1 = 0;
1039                         b1 = 1;
1040                         c1 = 0;                         // y
1041                         a2 = 1;
1042                         b2 = 0;
1043                         c2 = 0;                         // x
1044                 }
1045         if (VA.IsSame(V1))
1046                 if (VB.IsSame(V2))
1047                 {
1048                         a1 = 0;
1049                         b1 = -1;
1050                         c1 = 1;                         // 1-y
1051                         a2 = 1;
1052                         b2 = 0;
1053                         c2 = 0;                         // x
1054                 }
1055                 else
1056                 {
1057                         ASSERT(VB.IsSame(V0));
1058                         a1 = -1;
1059                         b1 = 0;
1060                         c1 = 1;                         // 1-x
1061                         a2 = 0;
1062                         b2 = 1;
1063                         c2 = 0;                         // y
1064                 }
1065         if (VA.IsSame(V2))
1066                 if (VB.IsSame(V3))
1067                 {
1068                         a1 = -1;
1069                         b1 = 0;
1070                         c1 = 1;                         // 1-x
1071                         a2 = 0;
1072                         b2 = -1;
1073                         c2 = 1;                         // 1-y
1074                 }
1075                 else
1076                 {
1077                         ASSERT(VB.IsSame(V1));
1078                         a1 = 0;
1079                         b1 = -1;
1080                         c1 = 1;                         // 1-y
1081                         a2 = -1;
1082                         b2 = 0;
1083                         c2 = 1;                         // 1-x
1084                 }
1085         if (VA.IsSame(V3))
1086                 if (VB.IsSame(V0))
1087                 {
1088                         a1 = 0;
1089                         b1 = 1;
1090                         c1 = 0;                         // y
1091                         a2 = -1;
1092                         b2 = 0;
1093                         c2 = 1;                         // 1-x
1094                 }
1095                 else
1096                 {
1097                         ASSERT(VB.IsSame(V2));
1098                         a1 = 1;
1099                         b1 = 0;
1100                         c1 = 0;                         // x
1101                         a2 = 0;
1102                         b2 = -1;
1103                         c2 = 1;                         // 1-y
1104                 }
1105 //      MESSAGE("X = " << c1 << "+ " << a1 << "*x + " << b1 << "*y");
1106 //      MESSAGE("Y = " << c2 << "+ " << a2 << "*x + " << b2 << "*y");
1107         conv.a1 = a1;
1108         conv.b1 = b1;
1109         conv.c1 = c1;
1110         conv.a2 = a2;
1111         conv.b2 = b2;
1112         conv.c2 = c2;
1113
1114         int nbdown = quad.side[0]->NbPoints();
1115         int nbright = quad.side[1]->NbPoints();
1116         conv.ia = int (a1);
1117         conv.ib = int (b1);
1118         conv.ic =
1119                 int (c1 * a1 * a1) * (nbdown - 1) + int (c1 * b1 * b1) * (nbright - 1);
1120         conv.ja = int (a2);
1121         conv.jb = int (b2);
1122         conv.jc =
1123                 int (c2 * a2 * a2) * (nbdown - 1) + int (c2 * b2 * b2) * (nbright - 1);
1124 //      MESSAGE("I " << conv.ia << " " << conv.ib << " " << conv.ic);
1125 //      MESSAGE("J " << conv.ja << " " << conv.jb << " " << conv.jc);
1126 }
1127
1128 //================================================================================
1129 /*!
1130  * \brief Find a vertex opposite to the given vertex of aQuads[0]
1131   * \param aVertex - the vertex
1132   * \param aFace - the face aVertex belongs to
1133   * \param aQuads - quads
1134   * \retval TopoDS_Vertex - found vertex
1135  */
1136 //================================================================================
1137
1138 TopoDS_Vertex StdMeshers_Hexa_3D::OppositeVertex(const TopoDS_Vertex& aVertex,
1139                                                  const TopTools_IndexedMapOfShape& aQuads0Vertices,
1140                                                  FaceQuadStruct* aQuads[6])
1141 {
1142   int i, j;
1143   for ( i = 1; i < 6; ++i )
1144   {
1145     TopoDS_Vertex VV[] = { aQuads[i]->side[0]->FirstVertex(),
1146                            aQuads[i]->side[0]->LastVertex() , 
1147                            aQuads[i]->side[2]->LastVertex() ,
1148                            aQuads[i]->side[2]->FirstVertex() };
1149     for ( j = 0; j < 4; ++j )
1150       if ( aVertex.IsSame( VV[ j ]))
1151         break;
1152     if ( j < 4 ) {
1153       int jPrev = j ? j - 1 : 3;
1154       int jNext = (j + 1) % 4;
1155       if ( aQuads0Vertices.Contains( VV[ jPrev ] ))
1156         return VV[ jNext ];
1157       else
1158         return VV[ jPrev ];
1159     }
1160   }
1161   return TopoDS_Vertex();
1162 }
1163
1164 //modified by NIZNHY-PKV Wed Nov 17 15:34:13 2004 f
1165 ///////////////////////////////////////////////////////////////////////////////
1166 //ZZ
1167 //#include <stdio.h>
1168
1169 //=======================================================================
1170 //function : ComputePentahedralMesh
1171 //purpose  : 
1172 //=======================================================================
1173
1174 SMESH_ComputeErrorPtr ComputePentahedralMesh(SMESH_Mesh &         aMesh,
1175                                              const TopoDS_Shape & aShape)
1176 {
1177   //printf(" ComputePentahedralMesh HERE\n");
1178   //
1179   bool bOK;
1180   SMESH_ComputeErrorPtr err = SMESH_ComputeError::New();
1181   //int iErr;
1182   StdMeshers_Penta_3D anAlgo;
1183   //
1184   bOK=anAlgo.Compute(aMesh, aShape);
1185   //
1186   err = anAlgo.GetComputeError();
1187   //
1188   if ( !bOK && anAlgo.ErrorStatus() == 5 )
1189   {
1190     static StdMeshers_Prism_3D * aPrism3D = 0;
1191     if ( !aPrism3D ) {
1192       SMESH_Gen* gen = aMesh.GetGen();
1193       aPrism3D = new StdMeshers_Prism_3D( gen->GetANewId(), 0, gen );
1194     }
1195     SMESH_Hypothesis::Hypothesis_Status aStatus;
1196     if ( aPrism3D->CheckHypothesis( aMesh, aShape, aStatus ) ) {
1197       aPrism3D->InitComputeError();
1198       bOK = aPrism3D->Compute( aMesh, aShape );
1199       err = aPrism3D->GetComputeError();
1200     }
1201   }
1202   return err;
1203 }
1204
1205
1206 //=======================================================================
1207 //function : EvaluatePentahedralMesh
1208 //purpose  : 
1209 //=======================================================================
1210
1211 bool EvaluatePentahedralMesh(SMESH_Mesh & aMesh,
1212                              const TopoDS_Shape & aShape,
1213                              MapShapeNbElems& aResMap)
1214 {
1215   StdMeshers_Penta_3D anAlgo;
1216   bool bOK = anAlgo.Evaluate(aMesh, aShape, aResMap);
1217
1218   //err = anAlgo.GetComputeError();
1219   //if ( !bOK && anAlgo.ErrorStatus() == 5 )
1220   if( !bOK ) {
1221     static StdMeshers_Prism_3D * aPrism3D = 0;
1222     if ( !aPrism3D ) {
1223       SMESH_Gen* gen = aMesh.GetGen();
1224       aPrism3D = new StdMeshers_Prism_3D( gen->GetANewId(), 0, gen );
1225     }
1226     SMESH_Hypothesis::Hypothesis_Status aStatus;
1227     if ( aPrism3D->CheckHypothesis( aMesh, aShape, aStatus ) ) {
1228       return aPrism3D->Evaluate(aMesh, aShape, aResMap);
1229     }
1230   }
1231
1232   return bOK;
1233 }
1234
1235