Salome HOME
merge from branch BR_SMDS_MEMIMP 29 nov 2010
[modules/smesh.git] / src / StdMeshers / StdMeshers_Prism_3D.cxx
1 //  Copyright (C) 2007-2010  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
23 //  SMESH SMESH : implementaion of SMESH idl descriptions
24 // File      : StdMeshers_Prism_3D.cxx
25 // Module    : SMESH
26 // Created   : Fri Oct 20 11:37:07 2006
27 // Author    : Edward AGAPOV (eap)
28 //
29 #include "StdMeshers_Prism_3D.hxx"
30
31 #include "StdMeshers_ProjectionUtils.hxx"
32 #include "SMESH_MesherHelper.hxx"
33 #include "SMDS_VolumeTool.hxx"
34 #include "SMDS_VolumeOfNodes.hxx"
35 #include "SMDS_EdgePosition.hxx"
36 #include "SMESH_Comment.hxx"
37
38 #include "utilities.h"
39
40 #include <BRep_Tool.hxx>
41 #include <Bnd_B3d.hxx>
42 #include <Geom2dAdaptor_Curve.hxx>
43 #include <Geom2d_Line.hxx>
44 #include <Geom_Curve.hxx>
45 #include <TopExp.hxx>
46 #include <TopExp_Explorer.hxx>
47 #include <TopTools_ListIteratorOfListOfShape.hxx>
48 #include <TopTools_MapOfShape.hxx>
49 #include <TopTools_SequenceOfShape.hxx>
50 #include <TopoDS.hxx>
51 #include <gp_Ax2.hxx>
52 #include <gp_Ax3.hxx>
53
54 using namespace std;
55
56 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
57 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
58 #define SHOWYXZ(msg, xyz) // {\
59 // gp_Pnt p (xyz); \
60 // cout << msg << " ("<< p.X() << "; " <<p.Y() << "; " <<p.Z() << ") " <<endl;\
61 // }
62
63 typedef StdMeshers_ProjectionUtils TAssocTool;
64 typedef SMESH_Comment              TCom;
65
66 enum { ID_BOT_FACE = SMESH_Block::ID_Fxy0,
67        ID_TOP_FACE = SMESH_Block::ID_Fxy1,
68        BOTTOM_EDGE = 0, TOP_EDGE, V0_EDGE, V1_EDGE, // edge IDs in face
69        NB_WALL_FACES = 4 }; //
70
71 namespace {
72
73   //================================================================================
74   /*!
75    * \brief Return iterator pointing to node column for the given parameter
76    * \param columnsMap - node column map
77    * \param parameter - parameter
78    * \retval TParam2ColumnMap::iterator - result
79    *
80    * it returns closest left column
81    */
82   //================================================================================
83
84   TParam2ColumnIt getColumn( const TParam2ColumnMap* columnsMap,
85                              const double            parameter )
86   {
87     TParam2ColumnIt u_col = columnsMap->upper_bound( parameter );
88     if ( u_col != columnsMap->begin() )
89       --u_col;
90     return u_col; // return left column
91   }
92
93   //================================================================================
94   /*!
95    * \brief Return nodes around given parameter and a ratio
96    * \param column - node column
97    * \param param - parameter
98    * \param node1 - lower node
99    * \param node2 - upper node
100    * \retval double - ratio
101    */
102   //================================================================================
103
104   double getRAndNodes( const TNodeColumn*     column,
105                        const double           param,
106                        const SMDS_MeshNode* & node1,
107                        const SMDS_MeshNode* & node2)
108   {
109     if ( param >= 1.0 || column->size() == 1) {
110       node1 = node2 = column->back();
111       return 0;
112     }
113
114     int i = int( param * ( column->size() - 1 ));
115     double u0 = double( i )/ double( column->size() - 1 );
116     double r = ( param - u0 ) * ( column->size() - 1 );
117
118     node1 = (*column)[ i ];
119     node2 = (*column)[ i + 1];
120     return r;
121   }
122
123   //================================================================================
124   /*!
125    * \brief Compute boundary parameters of face parts
126     * \param nbParts - nb of parts to split columns into
127     * \param columnsMap - node columns of the face to split
128     * \param params - computed parameters
129    */
130   //================================================================================
131
132   void splitParams( const int               nbParts,
133                     const TParam2ColumnMap* columnsMap,
134                     vector< double > &      params)
135   {
136     params.clear();
137     params.reserve( nbParts + 1 );
138     TParam2ColumnIt last_par_col = --columnsMap->end();
139     double par = columnsMap->begin()->first; // 0.
140     double parLast = last_par_col->first;
141     params.push_back( par );
142     for ( int i = 0; i < nbParts - 1; ++ i )
143     {
144       double partSize = ( parLast - par ) / double ( nbParts - i );
145       TParam2ColumnIt par_col = getColumn( columnsMap, par + partSize );
146       if ( par_col->first == par ) {
147         ++par_col;
148         if ( par_col == last_par_col ) {
149           while ( i < nbParts - 1 )
150             params.push_back( par + partSize * i++ );
151           break;
152         }
153       }
154       par = par_col->first;
155       params.push_back( par );
156     }
157     params.push_back( parLast ); // 1.
158   }
159
160   //================================================================================
161   /*!
162    * \brief Return coordinate system for z-th layer of nodes
163    */
164   //================================================================================
165
166   gp_Ax2 getLayerCoordSys(const int                           z,
167                           const vector< const TNodeColumn* >& columns,
168                           int&                                xColumn)
169   {
170     // gravity center of a layer
171     gp_XYZ O(0,0,0);
172     int vertexCol = -1;
173     for ( int i = 0; i < columns.size(); ++i )
174     {
175       O += gpXYZ( (*columns[ i ])[ z ]);
176       if ( vertexCol < 0 &&
177            columns[ i ]->front()->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX )
178         vertexCol = i;
179     }
180     O /= columns.size();
181
182     // Z axis
183     gp_Vec Z(0,0,0);
184     int iPrev = columns.size()-1;
185     for ( int i = 0; i < columns.size(); ++i )
186     {
187       gp_Vec v1( O, gpXYZ( (*columns[ iPrev ])[ z ]));
188       gp_Vec v2( O, gpXYZ( (*columns[ i ]    )[ z ]));
189       Z += v1 ^ v2;
190       iPrev = i;
191     }
192
193     if ( vertexCol >= 0 )
194     {
195       O = gpXYZ( (*columns[ vertexCol ])[ z ]);
196     }
197     if ( xColumn < 0 || xColumn >= columns.size() )
198     {
199       // select a column for X dir
200       double maxDist = 0;
201       for ( int i = 0; i < columns.size(); ++i )
202       {
203         double dist = ( O - gpXYZ((*columns[ i ])[ z ])).SquareModulus();
204         if ( dist > maxDist )
205         {
206           xColumn = i;
207           maxDist = dist;
208         }
209       }
210     }
211
212     // X axis
213     gp_Vec X( O, gpXYZ( (*columns[ xColumn ])[ z ]));
214
215     return gp_Ax2( O, Z, X);
216   }
217
218   //================================================================================
219   /*!
220    * \brief Removes submeshes meshed with regular grid from given list
221    *  \retval int - nb of removed submeshes
222    */
223   //================================================================================
224
225   int removeQuasiQuads(list< SMESH_subMesh* >& notQuadSubMesh)
226   {
227     int oldNbSM = notQuadSubMesh.size();
228     SMESHDS_Mesh* mesh = notQuadSubMesh.front()->GetFather()->GetMeshDS();
229     list< SMESH_subMesh* >::iterator smIt = notQuadSubMesh.begin();
230 #define __NEXT_SM { ++smIt; continue; }
231     while ( smIt != notQuadSubMesh.end() )
232     {
233       SMESH_subMesh* faceSm = *smIt;
234       SMESHDS_SubMesh* faceSmDS = faceSm->GetSubMeshDS();
235       int nbQuads = faceSmDS->NbElements();
236       if ( nbQuads == 0 ) __NEXT_SM;
237
238       // get oredered edges
239       list< TopoDS_Edge > orderedEdges;
240       list< int >         nbEdgesInWires;
241       TopoDS_Vertex       V000;
242       int nbWires = SMESH_Block::GetOrderedEdges( TopoDS::Face( faceSm->GetSubShape() ),
243                                                   V000, orderedEdges, nbEdgesInWires );
244       if ( nbWires != 1 || nbEdgesInWires.front() <= 4 )
245         __NEXT_SM;
246
247       // get nb of segements on edges
248       list<int> nbSegOnEdge;
249       list< TopoDS_Edge >::iterator edge = orderedEdges.begin();
250       for ( ; edge != orderedEdges.end(); ++edge )
251       {
252         if ( SMESHDS_SubMesh* edgeSmDS = mesh->MeshElements( *edge ))
253           nbSegOnEdge.push_back( edgeSmDS->NbElements() );
254         else
255           nbSegOnEdge.push_back(0);
256       }
257
258       // unite nbSegOnEdge of continues edges
259       int nbEdges = nbEdgesInWires.front();
260       list<int>::iterator nbSegIt = nbSegOnEdge.begin();
261       for ( edge = orderedEdges.begin(); edge != orderedEdges.end(); )
262       {
263         const TopoDS_Edge& e1 = *edge++;
264         const TopoDS_Edge& e2 = ( edge == orderedEdges.end() ? orderedEdges.front() : *edge );
265         if ( SMESH_Algo::IsContinuous( e1, e2 ))
266         {
267           // common vertex of continues edges must be shared by two 2D mesh elems of geom face
268           TopoDS_Vertex vCommon = TopExp::LastVertex( e1, true );
269           const SMDS_MeshNode* vNode = SMESH_Algo::VertexNode( vCommon, mesh );
270           int nbF = 0;
271           if ( vNode )
272           {
273             SMDS_ElemIteratorPtr fIt = vNode->GetInverseElementIterator(SMDSAbs_Face);
274             while ( fIt->more() )
275               nbF += faceSmDS->Contains( fIt->next() );
276           }
277           list<int>::iterator nbSegIt1 = nbSegIt++;
278           if ( !vNode || nbF == 2 ) // !vNode - two edges can be meshed as one
279           {
280             // unite
281             if ( nbSegIt == nbSegOnEdge.end() ) nbSegIt = nbSegOnEdge.begin();
282             *nbSegIt += *nbSegIt1;
283             nbSegOnEdge.erase( nbSegIt1 );
284             --nbEdges;
285           }
286         }
287         else
288         {
289           ++nbSegIt;
290         }
291       }
292       vector<int> nbSegVec( nbSegOnEdge.begin(), nbSegOnEdge.end());
293       if ( nbSegVec.size() == 4 &&
294            nbSegVec[0] == nbSegVec[2] &&
295            nbSegVec[1] == nbSegVec[3] &&
296            nbSegVec[0] * nbSegVec[1] == nbQuads
297            )
298         smIt = notQuadSubMesh.erase( smIt );
299       else
300         __NEXT_SM;
301     }
302
303     return oldNbSM - notQuadSubMesh.size();
304   }
305 }
306
307 //=======================================================================
308 //function : StdMeshers_Prism_3D
309 //purpose  : 
310 //=======================================================================
311
312 StdMeshers_Prism_3D::StdMeshers_Prism_3D(int hypId, int studyId, SMESH_Gen* gen)
313   :SMESH_3D_Algo(hypId, studyId, gen)
314 {
315   _name = "Prism_3D";
316   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);       // 1 bit per shape type
317   myProjectTriangles = false;
318 }
319
320 //================================================================================
321 /*!
322  * \brief Destructor
323  */
324 //================================================================================
325
326 StdMeshers_Prism_3D::~StdMeshers_Prism_3D()
327 {}
328
329 //=======================================================================
330 //function : CheckHypothesis
331 //purpose  : 
332 //=======================================================================
333
334 bool StdMeshers_Prism_3D::CheckHypothesis(SMESH_Mesh&                          aMesh,
335                                           const TopoDS_Shape&                  aShape,
336                                           SMESH_Hypothesis::Hypothesis_Status& aStatus)
337 {
338   // Check shape geometry
339 /*  PAL16229
340   aStatus = SMESH_Hypothesis::HYP_BAD_GEOMETRY;
341
342   // find not quadrangle faces
343   list< TopoDS_Shape > notQuadFaces;
344   int nbEdge, nbWire, nbFace = 0;
345   TopExp_Explorer exp( aShape, TopAbs_FACE );
346   for ( ; exp.More(); exp.Next() ) {
347     ++nbFace;
348     const TopoDS_Shape& face = exp.Current();
349     nbEdge = TAssocTool::Count( face, TopAbs_EDGE, 0 );
350     nbWire = TAssocTool::Count( face, TopAbs_WIRE, 0 );
351     if (  nbEdge!= 4 || nbWire!= 1 ) {
352       if ( !notQuadFaces.empty() ) {
353         if ( TAssocTool::Count( notQuadFaces.back(), TopAbs_EDGE, 0 ) != nbEdge ||
354              TAssocTool::Count( notQuadFaces.back(), TopAbs_WIRE, 0 ) != nbWire )
355           RETURN_BAD_RESULT("Different not quad faces");
356       }
357       notQuadFaces.push_back( face );
358     }
359   }
360   if ( !notQuadFaces.empty() )
361   {
362     if ( notQuadFaces.size() != 2 )
363       RETURN_BAD_RESULT("Bad nb not quad faces: " << notQuadFaces.size());
364
365     // check total nb faces
366     nbEdge = TAssocTool::Count( notQuadFaces.back(), TopAbs_EDGE, 0 );
367     if ( nbFace != nbEdge + 2 )
368       RETURN_BAD_RESULT("Bad nb of faces: " << nbFace << " but must be " << nbEdge + 2);
369   }
370 */
371   // no hypothesis
372   aStatus = SMESH_Hypothesis::HYP_OK;
373   return true;
374 }
375
376 //=======================================================================
377 //function : Compute
378 //purpose  : 
379 //=======================================================================
380
381 bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theShape)
382 {
383   SMESH_MesherHelper helper( theMesh );
384   myHelper = &helper;
385
386   myHelper->IsQuadraticSubMesh( theShape );
387
388   // Analyse mesh and geomerty to find block subshapes and submeshes
389   if ( !myBlock.Init( myHelper, theShape ))
390     return error( myBlock.GetError());
391
392   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
393
394   int volumeID = meshDS->ShapeToIndex( theShape );
395
396
397   // To compute coordinates of a node inside a block, it is necessary to know
398   // 1. normalized parameters of the node by which
399   // 2. coordinates of node projections on all block sub-shapes are computed
400
401   // So we fill projections on vertices at once as they are same for all nodes
402   myShapeXYZ.resize( myBlock.NbSubShapes() );
403   for ( int iV = SMESH_Block::ID_FirstV; iV < SMESH_Block::ID_FirstE; ++iV ) {
404     myBlock.VertexPoint( iV, myShapeXYZ[ iV ]);
405     SHOWYXZ("V point " <<iV << " ", myShapeXYZ[ iV ]);
406   }
407
408   // Projections on the top and bottom faces are taken from nodes existing
409   // on these faces; find correspondence between bottom and top nodes
410   myBotToColumnMap.clear();
411   if ( !assocOrProjBottom2Top() ) // it also fills myBotToColumnMap
412     return false;
413
414
415   // Create nodes inside the block
416
417   // try to use transformation (issue 0020680)
418   vector<gp_Trsf> trsf;
419   if ( myBlock.GetLayersTransformation(trsf))
420   {
421     // loop on nodes inside the bottom face
422     TNode2ColumnMap::iterator bot_column = myBotToColumnMap.begin();
423     for ( ; bot_column != myBotToColumnMap.end(); ++bot_column )
424     {
425       const TNode& tBotNode = bot_column->first; // bottom TNode
426       if ( tBotNode.GetPositionType() != SMDS_TOP_FACE )
427         continue; // node is not inside face 
428
429       // column nodes; middle part of the column are zero pointers
430       TNodeColumn& column = bot_column->second;
431       TNodeColumn::iterator columnNodes = column.begin();
432       for ( int z = 0; columnNodes != column.end(); ++columnNodes, ++z)
433       {
434         const SMDS_MeshNode* & node = *columnNodes;
435         if ( node ) continue; // skip bottom or top node
436
437         gp_XYZ coords = tBotNode.GetCoords();
438         trsf[z-1].Transforms( coords );
439         node = meshDS->AddNode( coords.X(), coords.Y(), coords.Z() );
440         meshDS->SetNodeInVolume( node, volumeID );
441       }
442     } // loop on bottom nodes
443   }
444   else // use block approach
445   {
446     // loop on nodes inside the bottom face
447     TNode prevBNode;
448     TNode2ColumnMap::iterator bot_column = myBotToColumnMap.begin();
449     for ( ; bot_column != myBotToColumnMap.end(); ++bot_column )
450     {
451       const TNode& tBotNode = bot_column->first; // bottom TNode
452       if ( tBotNode.GetPositionType() != SMDS_TOP_FACE )
453         continue; // node is not inside face 
454
455       // column nodes; middle part of the column are zero pointers
456       TNodeColumn& column = bot_column->second;
457
458       // compute bottom node parameters
459       gp_XYZ paramHint(-1,-1,-1);
460       if ( prevBNode.IsNeighbor( tBotNode ))
461         paramHint = prevBNode.GetParams();
462       if ( !myBlock.ComputeParameters( tBotNode.GetCoords(), tBotNode.ChangeParams(),
463                                        ID_BOT_FACE, paramHint ))
464         return error(TCom("Can't compute normalized parameters for node ")
465                      << tBotNode.myNode->GetID() << " on the face #"
466                      << myBlock.SubMesh( ID_BOT_FACE )->GetId() );
467       prevBNode = tBotNode;
468
469       myShapeXYZ[ ID_BOT_FACE ] = tBotNode.GetCoords();
470       gp_XYZ botParams          = tBotNode.GetParams();
471
472       // compute top node parameters
473       myShapeXYZ[ ID_TOP_FACE ] = gpXYZ( column.back() );
474       gp_XYZ topParams = botParams;
475       topParams.SetZ( 1 );
476       if ( column.size() > 2 ) {
477         gp_Pnt topCoords = myShapeXYZ[ ID_TOP_FACE ];
478         if ( !myBlock.ComputeParameters( topCoords, topParams, ID_TOP_FACE, topParams ))
479           return error(TCom("Can't compute normalized parameters ")
480                        << "for node " << column.back()->GetID()
481                        << " on the face #"<< column.back()->getshapeId() );
482       }
483
484       // vertical loop
485       TNodeColumn::iterator columnNodes = column.begin();
486       for ( int z = 0; columnNodes != column.end(); ++columnNodes, ++z)
487       {
488         const SMDS_MeshNode* & node = *columnNodes;
489         if ( node ) continue; // skip bottom or top node
490
491         // params of a node to create
492         double rz = (double) z / (double) ( column.size() - 1 );
493         gp_XYZ params = botParams * ( 1 - rz ) + topParams * rz;
494
495         // set coords on all faces and nodes
496         const int nbSideFaces = 4;
497         int sideFaceIDs[nbSideFaces] = { SMESH_Block::ID_Fx0z,
498                                          SMESH_Block::ID_Fx1z,
499                                          SMESH_Block::ID_F0yz,
500                                          SMESH_Block::ID_F1yz };
501         for ( int iF = 0; iF < nbSideFaces; ++iF )
502           if ( !setFaceAndEdgesXYZ( sideFaceIDs[ iF ], params, z ))
503             return false;
504
505         // compute coords for a new node
506         gp_XYZ coords;
507         if ( !SMESH_Block::ShellPoint( params, myShapeXYZ, coords ))
508           return error("Can't compute coordinates by normalized parameters");
509
510         SHOWYXZ("TOPFacePoint ",myShapeXYZ[ ID_TOP_FACE]);
511         SHOWYXZ("BOT Node "<< tBotNode.myNode->GetID(),gpXYZ(tBotNode.myNode));
512         SHOWYXZ("ShellPoint ",coords);
513
514         // create a node
515         node = meshDS->AddNode( coords.X(), coords.Y(), coords.Z() );
516         meshDS->SetNodeInVolume( node, volumeID );
517       }
518     } // loop on bottom nodes
519   }
520
521   // Create volumes
522
523   SMESHDS_SubMesh* smDS = myBlock.SubMeshDS( ID_BOT_FACE );
524   if ( !smDS ) return error(COMPERR_BAD_INPUT_MESH, "Null submesh");
525
526   // loop on bottom mesh faces
527   SMDS_ElemIteratorPtr faceIt = smDS->GetElements();
528   while ( faceIt->more() )
529   {
530     const SMDS_MeshElement* face = faceIt->next();
531     if ( !face || face->GetType() != SMDSAbs_Face )
532       continue;
533     int nbNodes = face->NbNodes();
534     if ( face->IsQuadratic() )
535       nbNodes /= 2;
536
537     // find node columns for each node
538     vector< const TNodeColumn* > columns( nbNodes );
539     for ( int i = 0; i < nbNodes; ++i )
540     {
541       const SMDS_MeshNode* n = face->GetNode( i );
542       if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ) {
543         TNode2ColumnMap::iterator bot_column = myBotToColumnMap.find( n );
544         if ( bot_column == myBotToColumnMap.end() )
545           return error(TCom("No nodes found above node ") << n->GetID() );
546         columns[ i ] = & bot_column->second;
547       }
548       else {
549         columns[ i ] = myBlock.GetNodeColumn( n );
550         if ( !columns[ i ] )
551           return error(TCom("No side nodes found above node ") << n->GetID() );
552       }
553     }
554     // create prisms
555     AddPrisms( columns, myHelper );
556
557   } // loop on bottom mesh faces
558         
559   return true;
560 }
561
562
563 //=======================================================================
564 //function : Evaluate
565 //purpose  : 
566 //=======================================================================
567
568 bool StdMeshers_Prism_3D::Evaluate(SMESH_Mesh& theMesh,
569                                    const TopoDS_Shape& theShape,
570                                    MapShapeNbElems& aResMap)
571 {
572   // find face contains only triangles
573   vector < SMESH_subMesh * >meshFaces;
574   TopTools_SequenceOfShape aFaces;
575   int NumBase = 0, i = 0, NbQFs = 0;
576   for (TopExp_Explorer exp(theShape, TopAbs_FACE); exp.More(); exp.Next()) {
577     i++;
578     aFaces.Append(exp.Current());
579     SMESH_subMesh *aSubMesh = theMesh.GetSubMesh(exp.Current());
580     meshFaces.push_back(aSubMesh);
581     MapShapeNbElemsItr anIt = aResMap.find(meshFaces[i-1]);
582     if( anIt==aResMap.end() ) {
583       SMESH_ComputeErrorPtr& smError = aSubMesh->GetComputeError();
584       smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
585       return false;
586     }
587     std::vector<int> aVec = (*anIt).second;
588     int nbtri = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
589     int nbqua = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
590     if( nbtri==0 && nbqua>0 ) {
591       NbQFs++;
592     }
593     if( nbtri>0 ) {
594       NumBase = i;
595     }
596   }
597
598   if(NbQFs<4) {
599     std::vector<int> aResVec(SMDSEntity_Last);
600     for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
601     SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
602     aResMap.insert(std::make_pair(sm,aResVec));
603     SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
604     smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
605     return false;
606   }
607
608   if(NumBase==0) NumBase = 1; // only quads => set 1 faces as base
609
610   // find number of 1d elems for base face
611   int nb1d = 0;
612   TopTools_MapOfShape Edges1;
613   for (TopExp_Explorer exp(aFaces.Value(NumBase), TopAbs_EDGE); exp.More(); exp.Next()) {
614     Edges1.Add(exp.Current());
615     SMESH_subMesh *sm = theMesh.GetSubMesh(exp.Current());
616     if( sm ) {
617       MapShapeNbElemsItr anIt = aResMap.find(sm);
618       if( anIt == aResMap.end() ) continue;
619       std::vector<int> aVec = (*anIt).second;
620       nb1d += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
621     }
622   }
623   // find face opposite to base face
624   int OppNum = 0;
625   for(i=1; i<=6; i++) {
626     if(i==NumBase) continue;
627     bool IsOpposite = true;
628     for(TopExp_Explorer exp(aFaces.Value(i), TopAbs_EDGE); exp.More(); exp.Next()) {
629       if( Edges1.Contains(exp.Current()) ) {
630         IsOpposite = false;
631         break;
632       }
633     }
634     if(IsOpposite) {
635       OppNum = i;
636       break;
637     }
638   }
639   // find number of 2d elems on side faces
640   int nb2d = 0;
641   for(i=1; i<=6; i++) {
642     if( i==OppNum || i==NumBase ) continue;
643     MapShapeNbElemsItr anIt = aResMap.find( meshFaces[i-1] );
644     if( anIt == aResMap.end() ) continue;
645     std::vector<int> aVec = (*anIt).second;
646     nb2d += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
647   }
648   
649   MapShapeNbElemsItr anIt = aResMap.find( meshFaces[NumBase-1] );
650   std::vector<int> aVec = (*anIt).second;
651   bool IsQuadratic = (aVec[SMDSEntity_Quad_Triangle]>aVec[SMDSEntity_Triangle]) ||
652                      (aVec[SMDSEntity_Quad_Quadrangle]>aVec[SMDSEntity_Quadrangle]);
653   int nb2d_face0_3 = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
654   int nb2d_face0_4 = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
655   int nb0d_face0 = aVec[SMDSEntity_Node];
656   int nb1d_face0_int = ( nb2d_face0_3*3 + nb2d_face0_4*4 - nb1d ) / 2;
657
658   std::vector<int> aResVec(SMDSEntity_Last);
659   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
660   if(IsQuadratic) {
661     aResVec[SMDSEntity_Quad_Penta] = nb2d_face0_3 * ( nb2d/nb1d );
662     aResVec[SMDSEntity_Quad_Hexa] = nb2d_face0_4 * ( nb2d/nb1d );
663     aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
664   }
665   else {
666     aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
667     aResVec[SMDSEntity_Penta] = nb2d_face0_3 * ( nb2d/nb1d );
668     aResVec[SMDSEntity_Hexa] = nb2d_face0_4 * ( nb2d/nb1d );
669   }
670   SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
671   aResMap.insert(std::make_pair(sm,aResVec));
672
673   return true;
674 }
675
676
677 //================================================================================
678 /*!
679  * \brief Create prisms
680  * \param columns - columns of nodes generated from nodes of a mesh face
681  * \param helper - helper initialized by mesh and shape to add prisms to
682  */
683 //================================================================================
684
685 void StdMeshers_Prism_3D::AddPrisms( vector<const TNodeColumn*> & columns,
686                                      SMESH_MesherHelper*          helper)
687 {
688   SMESHDS_Mesh * meshDS = helper->GetMeshDS();
689   int shapeID = helper->GetSubShapeID();
690
691   int nbNodes = columns.size();
692   int nbZ     = columns[0]->size();
693   if ( nbZ < 2 ) return;
694
695   // find out orientation
696   bool isForward = true;
697   SMDS_VolumeTool vTool;
698   int z = 1;
699   switch ( nbNodes ) {
700   case 3: {
701     const SMDS_MeshNode* botNodes[3] = { (*columns[0])[z-1],
702                                          (*columns[1])[z-1],
703                                          (*columns[2])[z-1] };
704     const SMDS_MeshNode* topNodes[3] = { (*columns[0])[z],
705                                          (*columns[1])[z],
706                                          (*columns[2])[z] };
707     SMDS_VolumeOfNodes tmpVol ( botNodes[0], botNodes[1], botNodes[2],
708                                 topNodes[0], topNodes[1], topNodes[2]);
709     vTool.Set( &tmpVol );
710     isForward  = vTool.IsForward();
711     break;
712   }
713   case 4: {
714     const SMDS_MeshNode* botNodes[4] = { (*columns[0])[z-1], (*columns[1])[z-1],
715                                          (*columns[2])[z-1], (*columns[3])[z-1] };
716     const SMDS_MeshNode* topNodes[4] = { (*columns[0])[z], (*columns[1])[z],
717                                          (*columns[2])[z], (*columns[3])[z] };
718     SMDS_VolumeOfNodes tmpVol ( botNodes[0], botNodes[1], botNodes[2], botNodes[3],
719                                 topNodes[0], topNodes[1], topNodes[2], topNodes[3]);
720     vTool.Set( &tmpVol );
721     isForward  = vTool.IsForward();
722     break;
723   }
724   }
725
726   // vertical loop on columns
727   for ( z = 1; z < nbZ; ++z )
728   {
729     SMDS_MeshElement* vol = 0;
730     switch ( nbNodes ) {
731
732     case 3: {
733       const SMDS_MeshNode* botNodes[3] = { (*columns[0])[z-1],
734                                            (*columns[1])[z-1],
735                                            (*columns[2])[z-1] };
736       const SMDS_MeshNode* topNodes[3] = { (*columns[0])[z],
737                                            (*columns[1])[z],
738                                            (*columns[2])[z] };
739       if ( isForward )
740         vol = helper->AddVolume( botNodes[0], botNodes[1], botNodes[2],
741                                  topNodes[0], topNodes[1], topNodes[2]);
742       else
743         vol = helper->AddVolume( topNodes[0], topNodes[1], topNodes[2],
744                                  botNodes[0], botNodes[1], botNodes[2]);
745       break;
746       }
747     case 4: {
748       const SMDS_MeshNode* botNodes[4] = { (*columns[0])[z-1], (*columns[1])[z-1],
749                                            (*columns[2])[z-1], (*columns[3])[z-1] };
750       const SMDS_MeshNode* topNodes[4] = { (*columns[0])[z], (*columns[1])[z],
751                                            (*columns[2])[z], (*columns[3])[z] };
752       if ( isForward )
753         vol = helper->AddVolume( botNodes[0], botNodes[1], botNodes[2], botNodes[3],
754                                  topNodes[0], topNodes[1], topNodes[2], topNodes[3]);
755       else
756         vol = helper->AddVolume( topNodes[0], topNodes[1], topNodes[2], topNodes[3],
757                                  botNodes[0], botNodes[1], botNodes[2], botNodes[3]);
758       break;
759       }
760     default:
761       // polyhedron
762       vector<const SMDS_MeshNode*> nodes( 2*nbNodes + 4*nbNodes);
763       vector<int> quantities( 2 + nbNodes, 4 );
764       quantities[0] = quantities[1] = nbNodes;
765       columns.resize( nbNodes + 1 );
766       columns[ nbNodes ] = columns[ 0 ];
767       for ( int i = 0; i < nbNodes; ++i ) {
768         nodes[ i         ] = (*columns[ i ])[z-1]; // bottom
769         nodes[ i+nbNodes ] = (*columns[ i ])[z  ]; // top
770         // side
771         int di = 2*nbNodes + 4*i - 1;
772         nodes[ di   ] = (*columns[i  ])[z-1];
773         nodes[ di+1 ] = (*columns[i+1])[z-1];
774         nodes[ di+2 ] = (*columns[i+1])[z  ];
775         nodes[ di+3 ] = (*columns[i  ])[z  ];
776       }
777       vol = meshDS->AddPolyhedralVolume( nodes, quantities );
778     }
779     if ( vol && shapeID > 0 )
780       meshDS->SetMeshElementOnShape( vol, shapeID );
781   }
782 }
783
784 //================================================================================
785 /*!
786  * \brief Find correspondence between bottom and top nodes
787  *  If elements on the bottom and top faces are topologically different,
788  *  and projection is possible and allowed, perform the projection
789  *  \retval bool - is a success or not
790  */
791 //================================================================================
792
793 bool StdMeshers_Prism_3D::assocOrProjBottom2Top()
794 {
795   SMESH_subMesh * botSM = myBlock.SubMesh( ID_BOT_FACE );
796   SMESH_subMesh * topSM = myBlock.SubMesh( ID_TOP_FACE );
797
798   SMESHDS_SubMesh * botSMDS = botSM->GetSubMeshDS();
799   SMESHDS_SubMesh * topSMDS = topSM->GetSubMeshDS();
800
801   if ( !botSMDS || botSMDS->NbElements() == 0 )
802     return error(TCom("No elememts on face #") << botSM->GetId());
803
804   bool needProject = false;
805   if ( !topSMDS || 
806        botSMDS->NbElements() != topSMDS->NbElements() ||
807        botSMDS->NbNodes()    != topSMDS->NbNodes())
808   {
809     MESSAGE("nb elem bot " << botSMDS->NbElements() << " top " << topSMDS->NbElements());
810     MESSAGE("nb node bot " << botSMDS->NbNodes() << " top " << topSMDS->NbNodes());
811     if ( myBlock.HasNotQuadElemOnTop() )
812       return error(TCom("Mesh on faces #") << botSM->GetId()
813                    <<" and #"<< topSM->GetId() << " seems different" );
814     needProject = true;
815   }
816
817   if ( 0/*needProject && !myProjectTriangles*/ )
818     return error(TCom("Mesh on faces #") << botSM->GetId()
819                  <<" and #"<< topSM->GetId() << " seems different" );
820   ///RETURN_BAD_RESULT("Need to project but not allowed");
821
822   if ( needProject )
823   {
824     return projectBottomToTop();
825   }
826
827   TopoDS_Face botFace = TopoDS::Face( myBlock.Shape( ID_BOT_FACE ));
828   TopoDS_Face topFace = TopoDS::Face( myBlock.Shape( ID_TOP_FACE ));
829   // associate top and bottom faces
830   TAssocTool::TShapeShapeMap shape2ShapeMap;
831   if ( !TAssocTool::FindSubShapeAssociation( botFace, myBlock.Mesh(),
832                                              topFace, myBlock.Mesh(),
833                                              shape2ShapeMap) )
834     return error(TCom("Topology of faces #") << botSM->GetId()
835                  <<" and #"<< topSM->GetId() << " seems different" );
836
837   // Find matching nodes of top and bottom faces
838   TNodeNodeMap n2nMap;
839   if ( ! TAssocTool::FindMatchingNodesOnFaces( botFace, myBlock.Mesh(),
840                                                topFace, myBlock.Mesh(),
841                                                shape2ShapeMap, n2nMap ))
842     return error(TCom("Mesh on faces #") << botSM->GetId()
843                  <<" and #"<< topSM->GetId() << " seems different" );
844
845   // Fill myBotToColumnMap
846
847   int zSize = myBlock.VerticalSize();
848   //TNode prevTNode;
849   TNodeNodeMap::iterator bN_tN = n2nMap.begin();
850   for ( ; bN_tN != n2nMap.end(); ++bN_tN )
851   {
852     const SMDS_MeshNode* botNode = bN_tN->first;
853     const SMDS_MeshNode* topNode = bN_tN->second;
854     if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
855       continue; // wall columns are contained in myBlock
856     // create node column
857     TNode bN( botNode );
858     TNode2ColumnMap::iterator bN_col = 
859       myBotToColumnMap.insert( make_pair ( bN, TNodeColumn() )).first;
860     TNodeColumn & column = bN_col->second;
861     column.resize( zSize );
862     column.front() = botNode;
863     column.back()  = topNode;
864   }
865   return true;
866 }
867
868 //================================================================================
869 /*!
870  * \brief Remove quadrangles from the top face and
871  * create triangles there by projection from the bottom
872  * \retval bool - a success or not
873  */
874 //================================================================================
875
876 bool StdMeshers_Prism_3D::projectBottomToTop()
877 {
878   SMESH_subMesh * botSM = myBlock.SubMesh( ID_BOT_FACE );
879   SMESH_subMesh * topSM = myBlock.SubMesh( ID_TOP_FACE );
880
881   SMESHDS_SubMesh * botSMDS = botSM->GetSubMeshDS();
882   SMESHDS_SubMesh * topSMDS = topSM->GetSubMeshDS();
883
884   if ( topSMDS )
885     topSM->ComputeStateEngine( SMESH_subMesh::CLEAN );
886
887   SMESHDS_Mesh* meshDS = myBlock.MeshDS();
888   int shapeID = myHelper->GetSubShapeID();
889   int topFaceID = meshDS->ShapeToIndex( topSM->GetSubShape() );
890
891   // Fill myBotToColumnMap
892
893   int zSize = myBlock.VerticalSize();
894   TNode prevTNode;
895   SMDS_NodeIteratorPtr nIt = botSMDS->GetNodes();
896   while ( nIt->more() )
897   {
898     const SMDS_MeshNode* botNode = nIt->next();
899     if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
900       continue; // strange
901     // compute bottom node params
902     TNode bN( botNode );
903     gp_XYZ paramHint(-1,-1,-1);
904     if ( prevTNode.IsNeighbor( bN ))
905       paramHint = prevTNode.GetParams();
906     if ( !myBlock.ComputeParameters( bN.GetCoords(), bN.ChangeParams(),
907                                      ID_BOT_FACE, paramHint ))
908       return error(TCom("Can't compute normalized parameters for node ")
909                    << botNode->GetID() << " on the face #"<< botSM->GetId() );
910     prevTNode = bN;
911     // compute top node coords
912     gp_XYZ topXYZ; gp_XY topUV;
913     if ( !myBlock.FacePoint( ID_TOP_FACE, bN.GetParams(), topXYZ ) ||
914          !myBlock.FaceUV   ( ID_TOP_FACE, bN.GetParams(), topUV ))
915       return error(TCom("Can't compute coordinates "
916                         "by normalized parameters on the face #")<< topSM->GetId() );
917     SMDS_MeshNode * topNode = meshDS->AddNode( topXYZ.X(),topXYZ.Y(),topXYZ.Z() );
918     meshDS->SetNodeOnFace( topNode, topFaceID, topUV.X(), topUV.Y() );
919     // create node column
920     TNode2ColumnMap::iterator bN_col = 
921       myBotToColumnMap.insert( make_pair ( bN, TNodeColumn() )).first;
922     TNodeColumn & column = bN_col->second;
923     column.resize( zSize );
924     column.front() = botNode;
925     column.back()  = topNode;
926   }
927
928   // Create top faces
929
930   // loop on bottom mesh faces
931   SMDS_ElemIteratorPtr faceIt = botSMDS->GetElements();
932   while ( faceIt->more() )
933   {
934     const SMDS_MeshElement* face = faceIt->next();
935     if ( !face || face->GetType() != SMDSAbs_Face )
936       continue;
937     int nbNodes = face->NbNodes();
938     if ( face->IsQuadratic() )
939       nbNodes /= 2;
940
941     // find top node in columns for each bottom node
942     vector< const SMDS_MeshNode* > nodes( nbNodes );
943     for ( int i = 0; i < nbNodes; ++i )
944     {
945       const SMDS_MeshNode* n = face->GetNode( nbNodes - i - 1 );
946       if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ) {
947         TNode2ColumnMap::iterator bot_column = myBotToColumnMap.find( n );
948         if ( bot_column == myBotToColumnMap.end() )
949           return error(TCom("No nodes found above node ") << n->GetID() );
950         nodes[ i ] = bot_column->second.back();
951       }
952       else {
953         const TNodeColumn* column = myBlock.GetNodeColumn( n );
954         if ( !column )
955           return error(TCom("No side nodes found above node ") << n->GetID() );
956         nodes[ i ] = column->back();
957       }
958     }
959     // create a face, with reversed orientation
960     SMDS_MeshElement* newFace = 0;
961     switch ( nbNodes ) {
962
963     case 3: {
964       newFace = myHelper->AddFace(nodes[0], nodes[1], nodes[2]);
965       break;
966       }
967     case 4: {
968       newFace = myHelper->AddFace( nodes[0], nodes[1], nodes[2], nodes[3] );
969       break;
970       }
971     default:
972       newFace = meshDS->AddPolygonalFace( nodes );
973     }
974     if ( newFace && shapeID > 0 )
975       meshDS->SetMeshElementOnShape( newFace, shapeID );
976   }
977
978   return true;
979 }
980
981 //================================================================================
982 /*!
983  * \brief Set projection coordinates of a node to a face and it's subshapes
984  * \param faceID - the face given by in-block ID
985  * \param params - node normalized parameters
986  * \retval bool - is a success
987  */
988 //================================================================================
989
990 bool StdMeshers_Prism_3D::setFaceAndEdgesXYZ( const int faceID, const gp_XYZ& params, int z )
991 {
992   // find base and top edges of the face
993   enum { BASE = 0, TOP, LEFT, RIGHT };
994   vector< int > edgeVec; // 0-base, 1-top
995   SMESH_Block::GetFaceEdgesIDs( faceID, edgeVec );
996
997   myBlock.EdgePoint( edgeVec[ BASE ], params, myShapeXYZ[ edgeVec[ BASE ]]);
998   myBlock.EdgePoint( edgeVec[ TOP  ], params, myShapeXYZ[ edgeVec[ TOP ]]);
999
1000   SHOWYXZ("\nparams ", params);
1001   SHOWYXZ("TOP is " <<edgeVec[ TOP ], myShapeXYZ[ edgeVec[ TOP]]);
1002   SHOWYXZ("BASE is "<<edgeVec[ BASE], myShapeXYZ[ edgeVec[ BASE]]);
1003
1004   if ( faceID == SMESH_Block::ID_Fx0z || faceID == SMESH_Block::ID_Fx1z )
1005   {
1006     myBlock.EdgePoint( edgeVec[ LEFT ], params, myShapeXYZ[ edgeVec[ LEFT ]]);
1007     myBlock.EdgePoint( edgeVec[ RIGHT ], params, myShapeXYZ[ edgeVec[ RIGHT ]]);
1008
1009     SHOWYXZ("VER "<<edgeVec[ LEFT], myShapeXYZ[ edgeVec[ LEFT]]);
1010     SHOWYXZ("VER "<<edgeVec[ RIGHT], myShapeXYZ[ edgeVec[ RIGHT]]);
1011   }
1012   myBlock.FacePoint( faceID, params, myShapeXYZ[ faceID ]);
1013   SHOWYXZ("FacePoint "<<faceID, myShapeXYZ[ faceID]);
1014
1015   return true;
1016 }
1017
1018 //================================================================================
1019 /*!
1020  * \brief Return true if this node and other one belong to one face
1021  */
1022 //================================================================================
1023
1024 bool TNode::IsNeighbor( const TNode& other ) const
1025 {
1026   if ( !other.myNode || !myNode ) return false;
1027
1028   SMDS_ElemIteratorPtr fIt = other.myNode->GetInverseElementIterator(SMDSAbs_Face);
1029   while ( fIt->more() )
1030     if ( fIt->next()->GetNodeIndex( myNode ) >= 0 )
1031       return true;
1032   return false;
1033 }
1034
1035 //================================================================================
1036 /*!
1037  * \brief Constructor. Initialization is needed
1038  */
1039 //================================================================================
1040
1041 StdMeshers_PrismAsBlock::StdMeshers_PrismAsBlock()
1042 {
1043   mySide = 0;
1044 }
1045
1046 StdMeshers_PrismAsBlock::~StdMeshers_PrismAsBlock()
1047 {
1048   if ( mySide ) {
1049     delete mySide; mySide = 0;
1050   }
1051 }
1052
1053 //================================================================================
1054 /*!
1055  * \brief Initialization.
1056  * \param helper - helper loaded with mesh and 3D shape
1057  * \param shape3D - a closed shell or solid
1058  * \retval bool - false if a mesh or a shape are KO
1059  */
1060 //================================================================================
1061
1062 bool StdMeshers_PrismAsBlock::Init(SMESH_MesherHelper* helper,
1063                                    const TopoDS_Shape& shape3D)
1064 {
1065   if ( mySide ) {
1066     delete mySide; mySide = 0;
1067   }
1068   vector< TSideFace* > sideFaces( NB_WALL_FACES, 0 );
1069   vector< pair< double, double> > params ( NB_WALL_FACES );
1070   mySide = new TSideFace( sideFaces, params );
1071
1072   myHelper = helper;
1073   SMESHDS_Mesh* meshDS = myHelper->GetMeshDS();
1074
1075   SMESH_Block::init();
1076   myShapeIDMap.Clear();
1077   myShapeIndex2ColumnMap.clear();
1078   
1079   int wallFaceIds[ NB_WALL_FACES ] = { // to walk around a block
1080     SMESH_Block::ID_Fx0z, SMESH_Block::ID_F1yz,
1081     SMESH_Block::ID_Fx1z, SMESH_Block::ID_F0yz
1082   };
1083
1084   myError = SMESH_ComputeError::New();
1085
1086   // -------------------------------------------------------------
1087   // Look for top and bottom faces: not quadrangle ones or meshed
1088   // with not quadrangle elements
1089   // -------------------------------------------------------------
1090
1091   list< SMESH_subMesh* > notQuadGeomSubMesh;
1092   list< SMESH_subMesh* > notQuadElemSubMesh;
1093   int nbFaces = 0;
1094   //
1095   SMESH_subMesh* mainSubMesh = myHelper->GetMesh()->GetSubMeshContaining( shape3D );
1096   if ( !mainSubMesh ) return error(COMPERR_BAD_INPUT_MESH,"Null submesh of shape3D");
1097
1098   // analyse face submeshes
1099   SMESH_subMeshIteratorPtr smIt = mainSubMesh->getDependsOnIterator(false,false);
1100   while ( smIt->more() )
1101   {
1102     SMESH_subMesh* sm = smIt->next();
1103     const TopoDS_Shape& face = sm->GetSubShape();
1104     if ( face.ShapeType() != TopAbs_FACE )
1105       continue;
1106     nbFaces++;
1107
1108     // is quadrangle face?
1109     list< TopoDS_Edge > orderedEdges;
1110     list< int >         nbEdgesInWires;
1111     TopoDS_Vertex       V000;
1112     int nbWires = GetOrderedEdges( TopoDS::Face( face ),
1113                                    V000, orderedEdges, nbEdgesInWires );
1114     if ( nbWires != 1 || nbEdgesInWires.front() != 4 )
1115       notQuadGeomSubMesh.push_back( sm );
1116
1117     // look for not quadrangle mesh elements
1118     if ( SMESHDS_SubMesh* smDS = sm->GetSubMeshDS() ) {
1119       bool hasNotQuad = false;
1120       SMDS_ElemIteratorPtr eIt = smDS->GetElements();
1121       while ( eIt->more() && !hasNotQuad ) {
1122         const SMDS_MeshElement* elem = eIt->next();
1123         if ( elem->GetType() == SMDSAbs_Face ) {
1124           int nbNodes = elem->NbNodes();
1125           if ( elem->IsQuadratic() )
1126             nbNodes /= 2;
1127           hasNotQuad = ( nbNodes != 4 );
1128         }
1129       }
1130       if ( hasNotQuad )
1131         notQuadElemSubMesh.push_back( sm );
1132     }
1133     else {
1134       return error(COMPERR_BAD_INPUT_MESH,TCom("Not meshed face #")<<sm->GetId());
1135     }
1136     // check if a quadrangle face is meshed with a quadranglar grid
1137     if ( notQuadGeomSubMesh.back() != sm &&
1138          notQuadElemSubMesh.back() != sm )
1139     {
1140       // count nb edges on face sides
1141       vector< int > nbEdges;
1142       nbEdges.reserve( nbEdgesInWires.front() );
1143       for ( list< TopoDS_Edge >::iterator edge = orderedEdges.begin();
1144             edge != orderedEdges.end(); ++edge )
1145       {
1146         if ( SMESHDS_SubMesh* smDS = meshDS->MeshElements( *edge ))
1147           nbEdges.push_back ( smDS->NbElements() );
1148         else
1149           nbEdges.push_back ( 0 );
1150       }
1151       int nbQuads = sm->GetSubMeshDS()->NbElements();
1152       if ( nbEdges[0] *  nbEdges[1] != nbQuads ||
1153            nbEdges[0] != nbEdges[2] ||
1154            nbEdges[1] != nbEdges[3] )
1155         notQuadElemSubMesh.push_back( sm );
1156     }
1157   }
1158
1159   // ----------------------------------------------------------------------
1160   // Analyse mesh and topology of faces: choose the bottom submesh.
1161   // If there are not quadrangle geom faces, they are top and bottom ones.
1162   // Not quadrangle geom faces must be only on top and bottom.
1163   // ----------------------------------------------------------------------
1164
1165   SMESH_subMesh * botSM = 0;
1166   SMESH_subMesh * topSM = 0;
1167
1168   int nbNotQuad       = notQuadGeomSubMesh.size();
1169   int nbNotQuadMeshed = notQuadElemSubMesh.size();
1170   bool hasNotQuad = ( nbNotQuad || nbNotQuadMeshed );
1171
1172   // detect bad cases
1173   if ( nbNotQuadMeshed > 2 )
1174   {
1175     return error(COMPERR_BAD_INPUT_MESH,
1176                  TCom("More than 2 faces with not quadrangle elements: ")
1177                  <<nbNotQuadMeshed);
1178   }
1179   int nbQuasiQuads = 0;
1180   if ( nbNotQuad > 0 && nbNotQuad != 2 )
1181   {
1182     // Issue 0020843 - one of side faces is quasi-quadrilateral.
1183     // Remove from notQuadGeomSubMesh faces meshed with regular grid
1184     nbQuasiQuads = removeQuasiQuads( notQuadGeomSubMesh );
1185     nbNotQuad -= nbQuasiQuads;
1186     if ( nbNotQuad > 0 && nbNotQuad != 2 )
1187       return error(COMPERR_BAD_SHAPE,
1188                    TCom("More than 2 not quadrilateral faces: ")
1189                    <<nbNotQuad);
1190   }
1191
1192   // get found submeshes
1193   if ( hasNotQuad )
1194   {
1195     if ( nbNotQuadMeshed > 0 ) botSM = notQuadElemSubMesh.front();
1196     else                       botSM = notQuadGeomSubMesh.front();
1197     if ( nbNotQuadMeshed > 1 ) topSM = notQuadElemSubMesh.back();
1198     else if ( nbNotQuad  > 1 ) topSM = notQuadGeomSubMesh.back();
1199   }
1200   // detect other bad cases
1201   if ( nbNotQuad == 2 && nbNotQuadMeshed > 0 ) {
1202     bool ok = false;
1203     if ( nbNotQuadMeshed == 1 )
1204       ok = ( find( notQuadGeomSubMesh.begin(),
1205                    notQuadGeomSubMesh.end(), botSM ) != notQuadGeomSubMesh.end() );
1206     else
1207       ok = ( notQuadGeomSubMesh == notQuadElemSubMesh );
1208     if ( !ok )
1209       return error(COMPERR_BAD_INPUT_MESH, "Side face meshed with not quadrangle elements");
1210   }
1211
1212   myNotQuadOnTop = ( nbNotQuadMeshed > 1 );
1213   MESSAGE("myNotQuadOnTop " << myNotQuadOnTop << " nbNotQuadMeshed " << nbNotQuadMeshed);
1214  
1215   // ----------------------------------------------------------
1216
1217   if ( nbNotQuad == 0 ) // Standard block of 6 quadrangle faces ?
1218   {
1219     // SMESH_Block will perform geometry analysis, we need just to find 2
1220     // connected vertices on top and bottom
1221
1222     TopoDS_Vertex Vbot, Vtop;
1223     if ( nbNotQuadMeshed > 0 ) // Look for vertices
1224     {
1225       TopTools_IndexedMapOfShape edgeMap;
1226       TopExp::MapShapes( botSM->GetSubShape(), TopAbs_EDGE, edgeMap );
1227       // vertex 1 is any vertex of the bottom face
1228       Vbot = TopExp::FirstVertex( TopoDS::Edge( edgeMap( 1 )));
1229       // vertex 2 is end vertex of edge sharing Vbot and not belonging to the bottom face
1230       TopTools_ListIteratorOfListOfShape ancestIt = Mesh()->GetAncestors( Vbot );
1231       for ( ; Vtop.IsNull() && ancestIt.More(); ancestIt.Next() )
1232       {
1233         const TopoDS_Shape & ancestor = ancestIt.Value();
1234         if ( ancestor.ShapeType() == TopAbs_EDGE && !edgeMap.FindIndex( ancestor ))
1235         {
1236           TopoDS_Vertex V1, V2;
1237           TopExp::Vertices( TopoDS::Edge( ancestor ), V1, V2);
1238           if      ( Vbot.IsSame ( V1 )) Vtop = V2;
1239           else if ( Vbot.IsSame ( V2 )) Vtop = V1;
1240           // check that Vtop belongs to shape3D
1241           TopExp_Explorer exp( shape3D, TopAbs_VERTEX );
1242           for ( ; exp.More(); exp.Next() )
1243             if ( Vtop.IsSame( exp.Current() ))
1244               break;
1245           if ( !exp.More() )
1246             Vtop.Nullify();
1247         }
1248       }
1249     }
1250     // get shell from shape3D
1251     TopoDS_Shell shell;
1252     TopExp_Explorer exp( shape3D, TopAbs_SHELL );
1253     int nbShell = 0;
1254     for ( ; exp.More(); exp.Next(), ++nbShell )
1255       shell = TopoDS::Shell( exp.Current() );
1256 //     if ( nbShell != 1 )
1257 //       RETURN_BAD_RESULT("There must be 1 shell in the block");
1258
1259     // Load geometry in SMESH_Block
1260     if ( !SMESH_Block::FindBlockShapes( shell, Vbot, Vtop, myShapeIDMap )) {
1261       if ( !hasNotQuad )
1262         return error(COMPERR_BAD_SHAPE, "Can't detect top and bottom of a prism");
1263     }
1264     else {
1265       if ( !botSM ) botSM = Mesh()->GetSubMeshContaining( myShapeIDMap( ID_BOT_FACE ));
1266       if ( !topSM ) topSM = Mesh()->GetSubMeshContaining( myShapeIDMap( ID_TOP_FACE ));
1267     }
1268
1269   } // end  Standard block of 6 quadrangle faces
1270   // --------------------------------------------------------
1271
1272   // Here the top and bottom faces are found
1273   if ( nbNotQuadMeshed == 2 ) // roughly check correspondence of horiz meshes
1274   {
1275 //     SMESHDS_SubMesh* topSMDS = topSM->GetSubMeshDS();
1276 //     SMESHDS_SubMesh* botSMDS = botSM->GetSubMeshDS();
1277 //     if ( topSMDS->NbNodes() != botSMDS->NbNodes() ||
1278 //          topSMDS->NbElements() != botSMDS->NbElements() )
1279 //       RETURN_BAD_RESULT("Top mesh doesn't correspond to bottom one");
1280   }
1281
1282   // ---------------------------------------------------------
1283   // If there are not quadrangle geom faces, we emulate
1284   // a block of 6 quadrangle faces.
1285   // Load SMESH_Block with faces and edges geometry
1286   // ---------------------------------------------------------
1287
1288   
1289   // find vertex 000 - the one with smallest coordinates (for easy DEBUG :-)
1290   TopoDS_Vertex V000;
1291   double minVal = DBL_MAX, minX, val;
1292   for ( TopExp_Explorer exp( botSM->GetSubShape(), TopAbs_VERTEX );
1293         exp.More(); exp.Next() )
1294   {
1295     const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
1296     gp_Pnt P = BRep_Tool::Pnt( v );
1297     val = P.X() + P.Y() + P.Z();
1298     if ( val < minVal || ( val == minVal && P.X() < minX )) {
1299       V000 = v;
1300       minVal = val;
1301       minX = P.X();
1302     }
1303   }
1304
1305   // Get ordered bottom edges
1306   list< TopoDS_Edge > orderedEdges;
1307   list< int >         nbVertexInWires;
1308   SMESH_Block::GetOrderedEdges( TopoDS::Face( botSM->GetSubShape().Reversed() ),
1309                                 V000, orderedEdges, nbVertexInWires );
1310 //   if ( nbVertexInWires.size() != 1 )
1311 //     RETURN_BAD_RESULT("Wrong prism geometry");
1312
1313   // Get Wall faces corresponding to the ordered bottom edges
1314   list< TopoDS_Face > wallFaces;
1315   if ( !GetWallFaces( Mesh(), shape3D, botSM->GetSubShape(), orderedEdges, wallFaces))
1316     return error(COMPERR_BAD_SHAPE, "Can't find side faces");
1317
1318   // Find columns of wall nodes and calculate edges' lengths
1319   // --------------------------------------------------------
1320
1321   myParam2ColumnMaps.clear();
1322   myParam2ColumnMaps.resize( orderedEdges.size() ); // total nb edges
1323
1324   int iE, nbEdges = nbVertexInWires.front(); // nb outer edges
1325   vector< double > edgeLength( nbEdges );
1326   map< double, int > len2edgeMap;
1327
1328   list< TopoDS_Edge >::iterator edgeIt = orderedEdges.begin();
1329   list< TopoDS_Face >::iterator faceIt = wallFaces.begin();
1330   for ( iE = 0; iE < nbEdges; ++edgeIt, ++faceIt )
1331   {
1332     TParam2ColumnMap & faceColumns = myParam2ColumnMaps[ iE ];
1333     if ( !myHelper->LoadNodeColumns( faceColumns, *faceIt, *edgeIt, meshDS ))
1334       return error(COMPERR_BAD_INPUT_MESH, TCom("Can't find regular quadrangle mesh ")
1335                    << "on a side face #" << MeshDS()->ShapeToIndex( *faceIt ));
1336
1337     SHOWYXZ("\np1 F "<<iE, gpXYZ(faceColumns.begin()->second.front() ));
1338     SHOWYXZ("p2 F "<<iE, gpXYZ(faceColumns.rbegin()->second.front() ));
1339     SHOWYXZ("V First "<<iE, BRep_Tool::Pnt( TopExp::FirstVertex(*edgeIt,true )));
1340
1341     edgeLength[ iE ] = SMESH_Algo::EdgeLength( *edgeIt );
1342
1343     if ( nbEdges < NB_WALL_FACES ) // fill map used to split faces
1344     {
1345       SMESHDS_SubMesh* smDS = meshDS->MeshElements( *edgeIt);
1346       if ( !smDS )
1347         return error(COMPERR_BAD_INPUT_MESH, TCom("Null submesh on the edge #")
1348                      << MeshDS()->ShapeToIndex( *edgeIt ));
1349       // assure length uniqueness
1350       edgeLength[ iE ] *= smDS->NbNodes() + edgeLength[ iE ] / ( 1000 + iE );
1351       len2edgeMap[ edgeLength[ iE ]] = iE;
1352     }
1353     ++iE;
1354   }
1355   // Load columns of internal edges (forming holes)
1356   // and fill map ShapeIndex to TParam2ColumnMap for them
1357   for ( ; edgeIt != orderedEdges.end() ; ++edgeIt, ++faceIt )
1358   {
1359     TParam2ColumnMap & faceColumns = myParam2ColumnMaps[ iE ];
1360     if ( !myHelper->LoadNodeColumns( faceColumns, *faceIt, *edgeIt, meshDS ))
1361       return error(COMPERR_BAD_INPUT_MESH, TCom("Can't find regular quadrangle mesh ")
1362                    << "on a side face #" << MeshDS()->ShapeToIndex( *faceIt ));
1363     // edge columns
1364     int id = MeshDS()->ShapeToIndex( *edgeIt );
1365     bool isForward = true; // meaningless for intenal wires
1366     myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1367     // columns for vertices
1368     // 1
1369     const SMDS_MeshNode* n0 = faceColumns.begin()->second.front();
1370     id = n0->getshapeId();
1371     myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1372     // 2
1373     const SMDS_MeshNode* n1 = faceColumns.rbegin()->second.front();
1374     id = n1->getshapeId();
1375     myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1376 //     SHOWYXZ("\np1 F "<<iE, gpXYZ(faceColumns.begin()->second.front() ));
1377 //     SHOWYXZ("p2 F "<<iE, gpXYZ(faceColumns.rbegin()->second.front() ));
1378 //     SHOWYXZ("V First "<<iE, BRep_Tool::Pnt( TopExp::FirstVertex(*edgeIt,true )));
1379     ++iE;
1380   }
1381
1382   // Create 4 wall faces of a block
1383   // -------------------------------
1384
1385   if ( nbEdges <= NB_WALL_FACES ) // ************* Split faces if necessary
1386   {
1387     map< int, int > iE2nbSplit;
1388     if ( nbEdges != NB_WALL_FACES ) // define how to split
1389     {
1390       if ( len2edgeMap.size() != nbEdges )
1391         RETURN_BAD_RESULT("Uniqueness of edge lengths not assured");
1392       map< double, int >::reverse_iterator maxLen_i = len2edgeMap.rbegin();
1393       map< double, int >::reverse_iterator midLen_i = ++len2edgeMap.rbegin();
1394       double maxLen = maxLen_i->first;
1395       double midLen = ( len2edgeMap.size() == 1 ) ? 0 : midLen_i->first;
1396       switch ( nbEdges ) {
1397       case 1: // 0-th edge is split into 4 parts
1398         iE2nbSplit.insert( make_pair( 0, 4 )); break;
1399       case 2: // either the longest edge is split into 3 parts, or both edges into halves
1400         if ( maxLen / 3 > midLen / 2 ) {
1401           iE2nbSplit.insert( make_pair( maxLen_i->second, 3 ));
1402         }
1403         else {
1404           iE2nbSplit.insert( make_pair( maxLen_i->second, 2 ));
1405           iE2nbSplit.insert( make_pair( midLen_i->second, 2 ));
1406         }
1407         break;
1408       case 3:
1409         // split longest into halves
1410         iE2nbSplit.insert( make_pair( maxLen_i->second, 2 ));
1411       }
1412     }
1413     // Create TSideFace's
1414     faceIt = wallFaces.begin();
1415     edgeIt = orderedEdges.begin();
1416     int iSide = 0;
1417     for ( iE = 0; iE < nbEdges; ++edgeIt, ++faceIt )
1418     {
1419      // split?
1420       map< int, int >::iterator i_nb = iE2nbSplit.find( iE );
1421       if ( i_nb != iE2nbSplit.end() ) {
1422         // split!
1423         int nbSplit = i_nb->second;
1424         vector< double > params;
1425         splitParams( nbSplit, &myParam2ColumnMaps[ iE ], params );
1426         bool isForward = ( edgeIt->Orientation() == TopAbs_FORWARD );
1427         for ( int i = 0; i < nbSplit; ++i ) {
1428           double f = ( isForward ? params[ i ] : params[ nbSplit - i-1 ]);
1429           double l = ( isForward ? params[ i+1 ] : params[ nbSplit - i ]);
1430           TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1431                                            *faceIt, *edgeIt,
1432                                            &myParam2ColumnMaps[ iE ], f, l );
1433           mySide->SetComponent( iSide++, comp );
1434         }
1435       }
1436       else {
1437         TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1438                                          *faceIt, *edgeIt,
1439                                          &myParam2ColumnMaps[ iE ]);
1440         mySide->SetComponent( iSide++, comp );
1441       }
1442       ++iE;
1443     }
1444   }
1445   else { // **************************** Unite faces
1446
1447     // unite first faces
1448     int nbExraFaces = nbEdges - 3;
1449     int iSide = 0, iE;
1450     double u0 = 0, sumLen = 0;
1451     for ( iE = 0; iE < nbExraFaces; ++iE )
1452       sumLen += edgeLength[ iE ];
1453
1454     vector< TSideFace* > components( nbExraFaces );
1455     vector< pair< double, double> > params( nbExraFaces );
1456     faceIt = wallFaces.begin();
1457     edgeIt = orderedEdges.begin();
1458     for ( iE = 0; iE < nbExraFaces; ++edgeIt, ++faceIt )
1459     {
1460       components[ iE ] = new TSideFace( myHelper, wallFaceIds[ iSide ],
1461                                         *faceIt, *edgeIt,
1462                                         &myParam2ColumnMaps[ iE ]);
1463       double u1 = u0 + edgeLength[ iE ] / sumLen;
1464       params[ iE ] = make_pair( u0 , u1 );
1465       u0 = u1;
1466       ++iE;
1467     }
1468     mySide->SetComponent( iSide++, new TSideFace( components, params ));
1469
1470     // fill the rest faces
1471     for ( ; iE < nbEdges; ++faceIt, ++edgeIt )
1472     {
1473       TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1474                                        *faceIt, *edgeIt,
1475                                        &myParam2ColumnMaps[ iE ]);
1476       mySide->SetComponent( iSide++, comp );
1477       ++iE;
1478     }
1479   }
1480
1481
1482   // Fill geometry fields of SMESH_Block
1483   // ------------------------------------
1484
1485   TopoDS_Face botF = TopoDS::Face( botSM->GetSubShape() );
1486   TopoDS_Face topF = TopoDS::Face( topSM->GetSubShape() );
1487
1488   vector< int > botEdgeIdVec;
1489   SMESH_Block::GetFaceEdgesIDs( ID_BOT_FACE, botEdgeIdVec );
1490
1491   bool isForward[NB_WALL_FACES] = { true, true, true, true };
1492   Adaptor2d_Curve2d* botPcurves[NB_WALL_FACES];
1493   Adaptor2d_Curve2d* topPcurves[NB_WALL_FACES];
1494
1495   for ( int iF = 0; iF < NB_WALL_FACES; ++iF )
1496   {
1497     TSideFace * sideFace = mySide->GetComponent( iF );
1498     if ( !sideFace )
1499       RETURN_BAD_RESULT("NULL TSideFace");
1500     int fID = sideFace->FaceID();
1501
1502     // fill myShapeIDMap
1503     if ( sideFace->InsertSubShapes( myShapeIDMap ) != 8 &&
1504          !sideFace->IsComplex())
1505       MESSAGE( ": Warning : InsertSubShapes() < 8 on side " << iF );
1506
1507     // side faces geometry
1508     Adaptor2d_Curve2d* pcurves[NB_WALL_FACES];
1509     if ( !sideFace->GetPCurves( pcurves ))
1510       RETURN_BAD_RESULT("TSideFace::GetPCurves() failed");
1511
1512     SMESH_Block::TFace& tFace = myFace[ fID - ID_FirstF ];
1513     tFace.Set( fID, sideFace->Surface(), pcurves, isForward );
1514
1515     SHOWYXZ( endl<<"F "<< iF << " id " << fID << " FRW " << sideFace->IsForward(), sideFace->Value(0,0));
1516     // edges 3D geometry
1517     vector< int > edgeIdVec;
1518     SMESH_Block::GetFaceEdgesIDs( fID, edgeIdVec );
1519     for ( int isMax = 0; isMax < 2; ++isMax ) {
1520       {
1521         int eID = edgeIdVec[ isMax ];
1522         SMESH_Block::TEdge& tEdge = myEdge[ eID - ID_FirstE ];
1523         tEdge.Set( eID, sideFace->HorizCurve(isMax), true);
1524         SHOWYXZ(eID<<" HOR"<<isMax<<"(0)", sideFace->HorizCurve(isMax)->Value(0));
1525         SHOWYXZ(eID<<" HOR"<<isMax<<"(1)", sideFace->HorizCurve(isMax)->Value(1));
1526       }
1527       {
1528         int eID = edgeIdVec[ isMax+2 ];
1529         SMESH_Block::TEdge& tEdge = myEdge[ eID - ID_FirstE  ];
1530         tEdge.Set( eID, sideFace->VertiCurve(isMax), true);
1531         SHOWYXZ(eID<<" VER"<<isMax<<"(0)", sideFace->VertiCurve(isMax)->Value(0));
1532         SHOWYXZ(eID<<" VER"<<isMax<<"(1)", sideFace->VertiCurve(isMax)->Value(1));
1533
1534         // corner points
1535         vector< int > vertexIdVec;
1536         SMESH_Block::GetEdgeVertexIDs( eID, vertexIdVec );
1537         myPnt[ vertexIdVec[0] - ID_FirstV ] = tEdge.GetCurve()->Value(0).XYZ();
1538         myPnt[ vertexIdVec[1] - ID_FirstV ] = tEdge.GetCurve()->Value(1).XYZ();
1539       }
1540     }
1541     // pcurves on horizontal faces
1542     for ( iE = 0; iE < NB_WALL_FACES; ++iE ) {
1543       if ( edgeIdVec[ BOTTOM_EDGE ] == botEdgeIdVec[ iE ] ) {
1544         botPcurves[ iE ] = sideFace->HorizPCurve( false, botF );
1545         topPcurves[ iE ] = sideFace->HorizPCurve( true,  topF );
1546         break;
1547       }
1548     }
1549     //sideFace->dumpNodes( 4 ); // debug
1550   }
1551   // horizontal faces geometry
1552   {
1553     SMESH_Block::TFace& tFace = myFace[ ID_BOT_FACE - ID_FirstF ];
1554     tFace.Set( ID_BOT_FACE, new BRepAdaptor_Surface( botF ), botPcurves, isForward );
1555     SMESH_Block::Insert( botF, ID_BOT_FACE, myShapeIDMap );
1556   }
1557   {
1558     SMESH_Block::TFace& tFace = myFace[ ID_TOP_FACE - ID_FirstF ];
1559     tFace.Set( ID_TOP_FACE, new BRepAdaptor_Surface( topF ), topPcurves, isForward );
1560     SMESH_Block::Insert( topF, ID_TOP_FACE, myShapeIDMap );
1561   }
1562
1563   // Fill map ShapeIndex to TParam2ColumnMap
1564   // ----------------------------------------
1565
1566   list< TSideFace* > fList;
1567   list< TSideFace* >::iterator fListIt;
1568   fList.push_back( mySide );
1569   for ( fListIt = fList.begin(); fListIt != fList.end(); ++fListIt)
1570   {
1571     int nb = (*fListIt)->NbComponents();
1572     for ( int i = 0; i < nb; ++i ) {
1573       if ( TSideFace* comp = (*fListIt)->GetComponent( i ))
1574         fList.push_back( comp );
1575     }
1576     if ( TParam2ColumnMap* cols = (*fListIt)->GetColumns()) {
1577       // columns for a base edge
1578       int id = MeshDS()->ShapeToIndex( (*fListIt)->BaseEdge() );
1579       bool isForward = (*fListIt)->IsForward();
1580       myShapeIndex2ColumnMap[ id ] = make_pair( cols, isForward );
1581
1582       // columns for vertices
1583       const SMDS_MeshNode* n0 = cols->begin()->second.front();
1584       id = n0->getshapeId();
1585       myShapeIndex2ColumnMap[ id ] = make_pair( cols, isForward );
1586
1587       const SMDS_MeshNode* n1 = cols->rbegin()->second.front();
1588       id = n1->getshapeId();
1589       myShapeIndex2ColumnMap[ id ] = make_pair( cols, !isForward );
1590     }
1591   }
1592
1593 //   gp_XYZ testPar(0.25, 0.25, 0), testCoord;
1594 //   if ( !FacePoint( ID_BOT_FACE, testPar, testCoord ))
1595 //     RETURN_BAD_RESULT("TEST FacePoint() FAILED");
1596 //   SHOWYXZ("IN TEST PARAM" , testPar);
1597 //   SHOWYXZ("OUT TEST CORD" , testCoord);
1598 //   if ( !ComputeParameters( testCoord, testPar , ID_BOT_FACE))
1599 //     RETURN_BAD_RESULT("TEST ComputeParameters() FAILED");
1600 //   SHOWYXZ("OUT TEST PARAM" , testPar);
1601
1602   return true;
1603 }
1604
1605 //================================================================================
1606 /*!
1607  * \brief Return pointer to column of nodes
1608  * \param node - bottom node from which the returned column goes up
1609  * \retval const TNodeColumn* - the found column
1610  */
1611 //================================================================================
1612
1613 const TNodeColumn* StdMeshers_PrismAsBlock::GetNodeColumn(const SMDS_MeshNode* node) const
1614 {
1615   int sID = node->getshapeId();
1616
1617   map<int, pair< TParam2ColumnMap*, bool > >::const_iterator col_frw =
1618     myShapeIndex2ColumnMap.find( sID );
1619   if ( col_frw != myShapeIndex2ColumnMap.end() ) {
1620     const TParam2ColumnMap* cols = col_frw->second.first;
1621     TParam2ColumnIt u_col = cols->begin();
1622     for ( ; u_col != cols->end(); ++u_col )
1623       if ( u_col->second[ 0 ] == node )
1624         return & u_col->second;
1625   }
1626   return 0;
1627 }
1628
1629 //=======================================================================
1630 //function : GetLayersTransformation
1631 //purpose  : Return transformations to get coordinates of nodes of each layer
1632 //           by nodes of the bottom. Layer is a set of nodes at a certain step
1633 //           from bottom to top.
1634 //=======================================================================
1635
1636 bool StdMeshers_PrismAsBlock::GetLayersTransformation(vector<gp_Trsf> & trsf) const
1637 {
1638   const int zSize = VerticalSize();
1639   if ( zSize < 3 ) return true;
1640   trsf.resize( zSize - 2 );
1641
1642   // Select some node columns by which we will define coordinate system of layers
1643
1644   vector< const TNodeColumn* > columns;
1645   {
1646     const TopoDS_Shape& baseFace = Shape(ID_BOT_FACE);
1647     list< TopoDS_Edge > orderedEdges;
1648     list< int >         nbEdgesInWires;
1649     GetOrderedEdges( TopoDS::Face( baseFace ), TopoDS_Vertex(), orderedEdges, nbEdgesInWires );
1650     bool isReverse;
1651     list< TopoDS_Edge >::iterator edgeIt = orderedEdges.begin();
1652     for ( int iE = 0; iE < nbEdgesInWires.front(); ++iE, ++edgeIt )
1653     {
1654       const TParam2ColumnMap& u2colMap =
1655         GetParam2ColumnMap( myHelper->GetMeshDS()->ShapeToIndex( *edgeIt ), isReverse );
1656       isReverse = ( edgeIt->Orientation() == TopAbs_REVERSED );
1657       double f = u2colMap.begin()->first, l = u2colMap.rbegin()->first;
1658       if ( isReverse ) swap ( f, l );
1659       const int nbCol = 5;
1660       for ( int i = 0; i < nbCol; ++i )
1661       {
1662         double u = f + i/double(nbCol) * ( l - f );
1663         const TNodeColumn* col = & getColumn( & u2colMap, u )->second;
1664         if ( columns.empty() || col != columns.back() )
1665           columns.push_back( col );
1666       }
1667     }
1668   }
1669
1670   // Find tolerance to check transformations
1671
1672   double tol2;
1673   {
1674     Bnd_B3d bndBox;
1675     for ( int i = 0; i < columns.size(); ++i )
1676       bndBox.Add( gpXYZ( columns[i]->front() ));
1677     tol2 = bndBox.SquareExtent() * 4 * 1e-4;
1678   }
1679
1680   // Compute transformations
1681
1682   int xCol = -1;
1683   gp_Trsf fromCsZ, toCs0;
1684   gp_Ax3 cs0 = getLayerCoordSys(0, columns, xCol );
1685   //double dist0 = cs0.Location().Distance( gpXYZ( (*columns[0])[0]));
1686   toCs0.SetTransformation( cs0 );
1687   for ( int z = 1; z < zSize-1; ++z )
1688   {
1689     gp_Ax3 csZ = getLayerCoordSys(z, columns, xCol );
1690     //double distZ = csZ.Location().Distance( gpXYZ( (*columns[0])[z]));
1691     fromCsZ.SetTransformation( csZ );
1692     fromCsZ.Invert();
1693     gp_Trsf& t = trsf[ z-1 ];
1694     t = fromCsZ * toCs0;
1695     //t.SetScaleFactor( distZ/dist0 ); - it does not work properly, wrong base point
1696
1697     // check a transformation
1698     for ( int i = 0; i < columns.size(); ++i )
1699     {
1700       gp_Pnt p0 = gpXYZ( (*columns[i])[0] );
1701       gp_Pnt pz = gpXYZ( (*columns[i])[z] );
1702       t.Transforms( p0.ChangeCoord() );
1703       if ( p0.SquareDistance( pz ) > tol2 )
1704         return false;
1705     }
1706   }
1707   return true;
1708 }
1709
1710 //================================================================================
1711 /*!
1712  * \brief Check curve orientation of a bootom edge
1713   * \param meshDS - mesh DS
1714   * \param columnsMap - node columns map of side face
1715   * \param bottomEdge - the bootom edge
1716   * \param sideFaceID - side face in-block ID
1717   * \retval bool - true if orientation coinside with in-block froward orientation
1718  */
1719 //================================================================================
1720
1721 bool StdMeshers_PrismAsBlock::IsForwardEdge(SMESHDS_Mesh*           meshDS,
1722                                             const TParam2ColumnMap& columnsMap,
1723                                             const TopoDS_Edge &     bottomEdge,
1724                                             const int               sideFaceID)
1725 {
1726   bool isForward = false;
1727   if ( TAssocTool::IsClosedEdge( bottomEdge ))
1728   {
1729     isForward = ( bottomEdge.Orientation() == TopAbs_FORWARD );
1730   }
1731   else
1732   {
1733     const TNodeColumn& firstCol = columnsMap.begin()->second;
1734     const SMDS_MeshNode* bottomNode = firstCol[0];
1735     TopoDS_Shape firstVertex = SMESH_MesherHelper::GetSubShapeByNode( bottomNode, meshDS );
1736     isForward = ( firstVertex.IsSame( TopExp::FirstVertex( bottomEdge, true )));
1737   }
1738   // on 2 of 4 sides first vertex is end
1739   if ( sideFaceID == ID_Fx1z || sideFaceID == ID_F0yz )
1740     isForward = !isForward;
1741   return isForward;
1742 }
1743
1744 //================================================================================
1745 /*!
1746  * \brief Find wall faces by bottom edges
1747  * \param mesh - the mesh
1748  * \param mainShape - the prism
1749  * \param bottomFace - the bottom face
1750  * \param bottomEdges - edges bounding the bottom face
1751  * \param wallFaces - faces list to fill in
1752  */
1753 //================================================================================
1754
1755 bool StdMeshers_PrismAsBlock::GetWallFaces( SMESH_Mesh*                     mesh,
1756                                             const TopoDS_Shape &            mainShape,
1757                                             const TopoDS_Shape &            bottomFace,
1758                                             const std::list< TopoDS_Edge >& bottomEdges,
1759                                             std::list< TopoDS_Face >&       wallFaces)
1760 {
1761   wallFaces.clear();
1762
1763   TopTools_IndexedMapOfShape faceMap;
1764   TopExp::MapShapes( mainShape, TopAbs_FACE, faceMap );
1765
1766   list< TopoDS_Edge >::const_iterator edge = bottomEdges.begin();
1767   for ( ; edge != bottomEdges.end(); ++edge )
1768   {
1769     TopTools_ListIteratorOfListOfShape ancestIt = mesh->GetAncestors( *edge );
1770     for ( ; ancestIt.More(); ancestIt.Next() )
1771     {
1772       const TopoDS_Shape& ancestor = ancestIt.Value();
1773       if ( ancestor.ShapeType() == TopAbs_FACE && // face
1774            !bottomFace.IsSame( ancestor ) &&      // not bottom
1775            faceMap.FindIndex( ancestor ))         // belongs to the prism
1776       {
1777         wallFaces.push_back( TopoDS::Face( ancestor ));
1778         break;
1779       }
1780     }
1781   }
1782   return ( wallFaces.size() == bottomEdges.size() );
1783 }
1784
1785 //================================================================================
1786 /*!
1787  * \brief Constructor
1788   * \param faceID - in-block ID
1789   * \param face - geom face
1790   * \param columnsMap - map of node columns
1791   * \param first - first normalized param
1792   * \param last - last normalized param
1793  */
1794 //================================================================================
1795
1796 StdMeshers_PrismAsBlock::TSideFace::TSideFace(SMESH_MesherHelper* helper,
1797                                               const int           faceID,
1798                                               const TopoDS_Face&  face,
1799                                               const TopoDS_Edge&  baseEdge,
1800                                               TParam2ColumnMap*   columnsMap,
1801                                               const double        first,
1802                                               const double        last):
1803   myID( faceID ),
1804   myParamToColumnMap( columnsMap ),
1805   myBaseEdge( baseEdge ),
1806   myHelper( helper )
1807 {
1808   mySurface.Initialize( face );
1809   myParams.resize( 1 );
1810   myParams[ 0 ] = make_pair( first, last );
1811   myIsForward = StdMeshers_PrismAsBlock::IsForwardEdge( myHelper->GetMeshDS(),
1812                                                         *myParamToColumnMap,
1813                                                         myBaseEdge, myID );
1814 }
1815
1816 //================================================================================
1817 /*!
1818  * \brief Constructor of complex side face
1819  */
1820 //================================================================================
1821
1822 StdMeshers_PrismAsBlock::TSideFace::
1823 TSideFace(const vector< TSideFace* >&             components,
1824           const vector< pair< double, double> > & params)
1825   :myID( components[0] ? components[0]->myID : 0 ),
1826    myParamToColumnMap( 0 ),
1827    myParams( params ),
1828    myIsForward( true ),
1829    myComponents( components ),
1830    myHelper( components[0] ? components[0]->myHelper : 0 )
1831 {}
1832 //================================================================================
1833 /*!
1834  * \brief Copy constructor
1835   * \param other - other side
1836  */
1837 //================================================================================
1838
1839 StdMeshers_PrismAsBlock::TSideFace::TSideFace( const TSideFace& other )
1840 {
1841   myID               = other.myID;
1842   mySurface          = other.mySurface;
1843   myBaseEdge         = other.myBaseEdge;
1844   myParams           = other.myParams;
1845   myIsForward        = other.myIsForward;
1846   myHelper           = other.myHelper;
1847   myParamToColumnMap = other.myParamToColumnMap;
1848
1849   myComponents.resize( other.myComponents.size());
1850   for (int i = 0 ; i < myComponents.size(); ++i )
1851     myComponents[ i ] = new TSideFace( *other.myComponents[ i ]);
1852 }
1853
1854 //================================================================================
1855 /*!
1856  * \brief Deletes myComponents
1857  */
1858 //================================================================================
1859
1860 StdMeshers_PrismAsBlock::TSideFace::~TSideFace()
1861 {
1862   for (int i = 0 ; i < myComponents.size(); ++i )
1863     if ( myComponents[ i ] )
1864       delete myComponents[ i ];
1865 }
1866
1867 //================================================================================
1868 /*!
1869  * \brief Return geometry of the vertical curve
1870   * \param isMax - true means curve located closer to (1,1,1) block point
1871   * \retval Adaptor3d_Curve* - curve adaptor
1872  */
1873 //================================================================================
1874
1875 Adaptor3d_Curve* StdMeshers_PrismAsBlock::TSideFace::VertiCurve(const bool isMax) const
1876 {
1877   if ( !myComponents.empty() ) {
1878     if ( isMax )
1879       return myComponents.back()->VertiCurve(isMax);
1880     else
1881       return myComponents.front()->VertiCurve(isMax);
1882   }
1883   double f = myParams[0].first, l = myParams[0].second;
1884   if ( !myIsForward ) std::swap( f, l );
1885   return new TVerticalEdgeAdaptor( myParamToColumnMap, isMax ? l : f );
1886 }
1887
1888 //================================================================================
1889 /*!
1890  * \brief Return geometry of the top or bottom curve
1891   * \param isTop - 
1892   * \retval Adaptor3d_Curve* - 
1893  */
1894 //================================================================================
1895
1896 Adaptor3d_Curve* StdMeshers_PrismAsBlock::TSideFace::HorizCurve(const bool isTop) const
1897 {
1898   return new THorizontalEdgeAdaptor( this, isTop );
1899 }
1900
1901 //================================================================================
1902 /*!
1903  * \brief Return pcurves
1904   * \param pcurv - array of 4 pcurves
1905   * \retval bool - is a success
1906  */
1907 //================================================================================
1908
1909 bool StdMeshers_PrismAsBlock::TSideFace::GetPCurves(Adaptor2d_Curve2d* pcurv[4]) const
1910 {
1911   int iEdge[ 4 ] = { BOTTOM_EDGE, TOP_EDGE, V0_EDGE, V1_EDGE };
1912
1913   for ( int i = 0 ; i < 4 ; ++i ) {
1914     Handle(Geom2d_Line) line;
1915     switch ( iEdge[ i ] ) {
1916     case TOP_EDGE:
1917       line = new Geom2d_Line( gp_Pnt2d( 0, 1 ), gp::DX2d() ); break;
1918     case BOTTOM_EDGE:
1919       line = new Geom2d_Line( gp::Origin2d(), gp::DX2d() ); break;
1920     case V0_EDGE:
1921       line = new Geom2d_Line( gp::Origin2d(), gp::DY2d() ); break;
1922     case V1_EDGE:
1923       line = new Geom2d_Line( gp_Pnt2d( 1, 0 ), gp::DY2d() ); break;
1924     }
1925     pcurv[ i ] = new Geom2dAdaptor_Curve( line, 0, 1 );
1926   }
1927   return true;
1928 }
1929
1930 //================================================================================
1931 /*!
1932  * \brief Returns geometry of pcurve on a horizontal face
1933   * \param isTop - is top or bottom face
1934   * \param horFace - a horizontal face
1935   * \retval Adaptor2d_Curve2d* - curve adaptor
1936  */
1937 //================================================================================
1938
1939 Adaptor2d_Curve2d*
1940 StdMeshers_PrismAsBlock::TSideFace::HorizPCurve(const bool         isTop,
1941                                                 const TopoDS_Face& horFace) const
1942 {
1943   return new TPCurveOnHorFaceAdaptor( this, isTop, horFace );
1944 }
1945
1946 //================================================================================
1947 /*!
1948  * \brief Return a component corresponding to parameter
1949   * \param U - parameter along a horizontal size
1950   * \param localU - parameter along a horizontal size of a component
1951   * \retval TSideFace* - found component
1952  */
1953 //================================================================================
1954
1955 StdMeshers_PrismAsBlock::TSideFace*
1956 StdMeshers_PrismAsBlock::TSideFace::GetComponent(const double U,double & localU) const
1957 {
1958   localU = U;
1959   if ( myComponents.empty() )
1960     return const_cast<TSideFace*>( this );
1961
1962   int i;
1963   for ( i = 0; i < myComponents.size(); ++i )
1964     if ( U < myParams[ i ].second )
1965       break;
1966   if ( i >= myComponents.size() )
1967     i = myComponents.size() - 1;
1968
1969   double f = myParams[ i ].first, l = myParams[ i ].second;
1970   localU = ( U - f ) / ( l - f );
1971   return myComponents[ i ];
1972 }
1973
1974 //================================================================================
1975 /*!
1976  * \brief Find node columns for a parameter
1977   * \param U - parameter along a horizontal edge
1978   * \param col1 - the 1st found column
1979   * \param col2 - the 2nd found column
1980   * \retval r - normalized position of U between the found columns
1981  */
1982 //================================================================================
1983
1984 double StdMeshers_PrismAsBlock::TSideFace::GetColumns(const double      U,
1985                                                       TParam2ColumnIt & col1,
1986                                                       TParam2ColumnIt & col2) const
1987 {
1988   double u = U, r = 0;
1989   if ( !myComponents.empty() ) {
1990     TSideFace * comp = GetComponent(U,u);
1991     return comp->GetColumns( u, col1, col2 );
1992   }
1993
1994   if ( !myIsForward )
1995     u = 1 - u;
1996   double f = myParams[0].first, l = myParams[0].second;
1997   u = f + u * ( l - f );
1998
1999   col1 = col2 = getColumn( myParamToColumnMap, u );
2000   if ( ++col2 == myParamToColumnMap->end() ) {
2001     --col2;
2002     r = 0.5;
2003   }
2004   else {
2005     double uf = col1->first;
2006     double ul = col2->first;
2007     r = ( u - uf ) / ( ul - uf );
2008   }
2009   return r;
2010 }
2011
2012 //================================================================================
2013 /*!
2014  * \brief Return coordinates by normalized params
2015   * \param U - horizontal param
2016   * \param V - vertical param
2017   * \retval gp_Pnt - result point
2018  */
2019 //================================================================================
2020
2021 gp_Pnt StdMeshers_PrismAsBlock::TSideFace::Value(const Standard_Real U,
2022                                                  const Standard_Real V) const
2023 {
2024   if ( !myComponents.empty() ) {
2025     double u;
2026     TSideFace * comp = GetComponent(U,u);
2027     return comp->Value( u, V );
2028   }
2029
2030   TParam2ColumnIt u_col1, u_col2;
2031   double vR, hR = GetColumns( U, u_col1, u_col2 );
2032
2033   const SMDS_MeshNode* n1 = 0;
2034   const SMDS_MeshNode* n2 = 0;
2035   const SMDS_MeshNode* n3 = 0;
2036   const SMDS_MeshNode* n4 = 0;
2037
2038   // BEGIN issue 0020680: EDF 1252 SMESH: Bad cell created by Radial prism in center of torus
2039   // Workaround for a wrongly located point returned by mySurface.Value() for
2040   // UV located near boundary of BSpline surface.
2041   // To bypass the problem, we take point from 3D curve of edge.
2042   // It solves pb of the bloc_fiss_new.py
2043   const double tol = 1e-3;
2044   if ( V < tol || V+tol >= 1. )
2045   {
2046     n1 = V < tol ? u_col1->second.front() : u_col1->second.back();
2047     n3 = V < tol ? u_col2->second.front() : u_col2->second.back();
2048     TopoDS_Edge edge;
2049     if ( V < tol )
2050     {
2051       edge = myBaseEdge;
2052     }
2053     else
2054     {
2055       TopoDS_Shape s = myHelper->GetSubShapeByNode( n1, myHelper->GetMeshDS() );
2056       if ( s.ShapeType() != TopAbs_EDGE )
2057         s = myHelper->GetSubShapeByNode( n3, myHelper->GetMeshDS() );
2058       if ( s.ShapeType() == TopAbs_EDGE )
2059         edge = TopoDS::Edge( s );
2060     }
2061     if ( !edge.IsNull() )
2062     {
2063       double u1 = myHelper->GetNodeU( edge, n1 );
2064       double u3 = myHelper->GetNodeU( edge, n3 );
2065       double u = u1 * ( 1 - hR ) + u3 * hR;
2066       TopLoc_Location loc; double f,l;
2067       Handle(Geom_Curve) curve = BRep_Tool::Curve( edge,loc,f,l );
2068       return curve->Value( u ).Transformed( loc );
2069     }
2070   }
2071   // END issue 0020680: EDF 1252 SMESH: Bad cell created by Radial prism in center of torus
2072
2073   vR = getRAndNodes( & u_col1->second, V, n1, n2 );
2074   vR = getRAndNodes( & u_col2->second, V, n3, n4 );
2075   
2076   gp_XY uv1 = myHelper->GetNodeUV( mySurface.Face(), n1, n4);
2077   gp_XY uv2 = myHelper->GetNodeUV( mySurface.Face(), n2, n3);
2078   gp_XY uv12 = uv1 * ( 1 - vR ) + uv2 * vR;
2079
2080   gp_XY uv3 = myHelper->GetNodeUV( mySurface.Face(), n3, n2);
2081   gp_XY uv4 = myHelper->GetNodeUV( mySurface.Face(), n4, n1);
2082   gp_XY uv34 = uv3 * ( 1 - vR ) + uv4 * vR;
2083
2084   gp_XY uv = uv12 * ( 1 - hR ) + uv34 * hR;
2085
2086   gp_Pnt p = mySurface.Value( uv.X(), uv.Y() );
2087   return p;
2088 }
2089
2090
2091 //================================================================================
2092 /*!
2093  * \brief Return boundary edge
2094   * \param edge - edge index
2095   * \retval TopoDS_Edge - found edge
2096  */
2097 //================================================================================
2098
2099 TopoDS_Edge StdMeshers_PrismAsBlock::TSideFace::GetEdge(const int iEdge) const
2100 {
2101   if ( !myComponents.empty() ) {
2102     switch ( iEdge ) {
2103     case V0_EDGE : return myComponents.front()->GetEdge( iEdge );
2104     case V1_EDGE : return myComponents.back() ->GetEdge( iEdge );
2105     default: return TopoDS_Edge();
2106     }
2107   }
2108   TopoDS_Shape edge;
2109   const SMDS_MeshNode* node = 0;
2110   SMESHDS_Mesh * meshDS = myHelper->GetMesh()->GetMeshDS();
2111   TNodeColumn* column;
2112
2113   switch ( iEdge ) {
2114   case TOP_EDGE:
2115   case BOTTOM_EDGE:
2116     column = & (( ++myParamToColumnMap->begin())->second );
2117     node = ( iEdge == TOP_EDGE ) ? column->back() : column->front();
2118     edge = myHelper->GetSubShapeByNode ( node, meshDS );
2119     if ( edge.ShapeType() == TopAbs_VERTEX ) {
2120       column = & ( myParamToColumnMap->begin()->second );
2121       node = ( iEdge == TOP_EDGE ) ? column->back() : column->front();
2122     }
2123     break;
2124   case V0_EDGE:
2125   case V1_EDGE: {
2126     bool back = ( iEdge == V1_EDGE );
2127     if ( !myIsForward ) back = !back;
2128     if ( back )
2129       column = & ( myParamToColumnMap->rbegin()->second );
2130     else
2131       column = & ( myParamToColumnMap->begin()->second );
2132     if ( column->size() > 0 )
2133       edge = myHelper->GetSubShapeByNode( (*column)[ 1 ], meshDS );
2134     if ( edge.IsNull() || edge.ShapeType() == TopAbs_VERTEX )
2135       node = column->front();
2136     break;
2137   }
2138   default:;
2139   }
2140   if ( !edge.IsNull() && edge.ShapeType() == TopAbs_EDGE )
2141     return TopoDS::Edge( edge );
2142
2143   // find edge by 2 vertices
2144   TopoDS_Shape V1 = edge;
2145   TopoDS_Shape V2 = myHelper->GetSubShapeByNode( node, meshDS );
2146   if ( V2.ShapeType() == TopAbs_VERTEX && !V2.IsSame( V1 ))
2147   {
2148     TopTools_ListIteratorOfListOfShape ancestIt =
2149       myHelper->GetMesh()->GetAncestors( V1 );
2150     for ( ; ancestIt.More(); ancestIt.Next() )
2151     {
2152       const TopoDS_Shape & ancestor = ancestIt.Value();
2153       if ( ancestor.ShapeType() == TopAbs_EDGE )
2154         for ( TopExp_Explorer e( ancestor, TopAbs_VERTEX ); e.More(); e.Next() )
2155           if ( V2.IsSame( e.Current() ))
2156             return TopoDS::Edge( ancestor );
2157     }
2158   }
2159   return TopoDS_Edge();
2160 }
2161
2162 //================================================================================
2163 /*!
2164  * \brief Fill block subshapes
2165   * \param shapeMap - map to fill in
2166   * \retval int - nb inserted subshapes
2167  */
2168 //================================================================================
2169
2170 int StdMeshers_PrismAsBlock::TSideFace::InsertSubShapes(TBlockShapes& shapeMap) const
2171 {
2172   int nbInserted = 0;
2173
2174   // Insert edges
2175   vector< int > edgeIdVec;
2176   SMESH_Block::GetFaceEdgesIDs( myID, edgeIdVec );
2177
2178   for ( int i = BOTTOM_EDGE; i <=V1_EDGE ; ++i ) {
2179     TopoDS_Edge e = GetEdge( i );
2180     if ( !e.IsNull() ) {
2181       nbInserted += SMESH_Block::Insert( e, edgeIdVec[ i ], shapeMap);
2182     }
2183   }
2184
2185   // Insert corner vertices
2186
2187   TParam2ColumnIt col1, col2 ;
2188   vector< int > vertIdVec;
2189
2190   // from V0 column
2191   SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ V0_EDGE ], vertIdVec);
2192   GetColumns(0, col1, col2 );
2193   const SMDS_MeshNode* node0 = col1->second.front();
2194   const SMDS_MeshNode* node1 = col1->second.back();
2195   TopoDS_Shape v0 = myHelper->GetSubShapeByNode( node0, myHelper->GetMeshDS());
2196   TopoDS_Shape v1 = myHelper->GetSubShapeByNode( node1, myHelper->GetMeshDS());
2197   if ( v0.ShapeType() == TopAbs_VERTEX ) {
2198     nbInserted += SMESH_Block::Insert( v0, vertIdVec[ 0 ], shapeMap);
2199   }
2200   if ( v1.ShapeType() == TopAbs_VERTEX ) {
2201     nbInserted += SMESH_Block::Insert( v1, vertIdVec[ 1 ], shapeMap);
2202   }
2203   
2204   // from V1 column
2205   SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ V1_EDGE ], vertIdVec);
2206   GetColumns(1, col1, col2 );
2207   node0 = col2->second.front();
2208   node1 = col2->second.back();
2209   v0 = myHelper->GetSubShapeByNode( node0, myHelper->GetMeshDS());
2210   v1 = myHelper->GetSubShapeByNode( node1, myHelper->GetMeshDS());
2211   if ( v0.ShapeType() == TopAbs_VERTEX ) {
2212     nbInserted += SMESH_Block::Insert( v0, vertIdVec[ 0 ], shapeMap);
2213   }
2214   if ( v1.ShapeType() == TopAbs_VERTEX ) {
2215     nbInserted += SMESH_Block::Insert( v1, vertIdVec[ 1 ], shapeMap);
2216   }
2217
2218 //   TopoDS_Vertex V0, V1, Vcom;
2219 //   TopExp::Vertices( myBaseEdge, V0, V1, true );
2220 //   if ( !myIsForward ) std::swap( V0, V1 );
2221
2222 //   // bottom vertex IDs
2223 //   SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ _u0 ], vertIdVec);
2224 //   SMESH_Block::Insert( V0, vertIdVec[ 0 ], shapeMap);
2225 //   SMESH_Block::Insert( V1, vertIdVec[ 1 ], shapeMap);
2226
2227 //   TopoDS_Edge sideEdge = GetEdge( V0_EDGE );
2228 //   if ( sideEdge.IsNull() || !TopExp::CommonVertex( botEdge, sideEdge, Vcom ))
2229 //     return false;
2230
2231 //   // insert one side edge
2232 //   int edgeID;
2233 //   if ( Vcom.IsSame( V0 )) edgeID = edgeIdVec[ _v0 ];
2234 //   else                    edgeID = edgeIdVec[ _v1 ];
2235 //   SMESH_Block::Insert( sideEdge, edgeID, shapeMap);
2236
2237 //   // top vertex of the side edge
2238 //   SMESH_Block::GetEdgeVertexIDs( edgeID, vertIdVec);
2239 //   TopoDS_Vertex Vtop = TopExp::FirstVertex( sideEdge );
2240 //   if ( Vcom.IsSame( Vtop ))
2241 //     Vtop = TopExp::LastVertex( sideEdge );
2242 //   SMESH_Block::Insert( Vtop, vertIdVec[ 1 ], shapeMap);
2243
2244 //   // other side edge
2245 //   sideEdge = GetEdge( V1_EDGE );
2246 //   if ( sideEdge.IsNull() )
2247 //     return false;
2248 //   if ( edgeID = edgeIdVec[ _v1 ]) edgeID = edgeIdVec[ _v0 ];
2249 //   else                            edgeID = edgeIdVec[ _v1 ];
2250 //   SMESH_Block::Insert( sideEdge, edgeID, shapeMap);
2251   
2252 //   // top edge
2253 //   TopoDS_Edge topEdge = GetEdge( TOP_EDGE );
2254 //   SMESH_Block::Insert( topEdge, edgeIdVec[ _u1 ], shapeMap);
2255
2256 //   // top vertex of the other side edge
2257 //   if ( !TopExp::CommonVertex( topEdge, sideEdge, Vcom ))
2258 //     return false;
2259 //   SMESH_Block::GetEdgeVertexIDs( edgeID, vertIdVec );
2260 //   SMESH_Block::Insert( Vcom, vertIdVec[ 1 ], shapeMap);
2261
2262   return nbInserted;
2263 }
2264
2265 //================================================================================
2266 /*!
2267  * \brief Dump ids of nodes of sides
2268  */
2269 //================================================================================
2270
2271 void StdMeshers_PrismAsBlock::TSideFace::dumpNodes(int nbNodes) const
2272 {
2273 #ifdef _DEBUG_
2274   cout << endl << "NODES OF FACE "; SMESH_Block::DumpShapeID( myID, cout ) << endl;
2275   THorizontalEdgeAdaptor* hSize0 = (THorizontalEdgeAdaptor*) HorizCurve(0);
2276   cout << "Horiz side 0: "; hSize0->dumpNodes(nbNodes); cout << endl;
2277   THorizontalEdgeAdaptor* hSize1 = (THorizontalEdgeAdaptor*) HorizCurve(1);
2278   cout << "Horiz side 1: "; hSize1->dumpNodes(nbNodes); cout << endl;
2279   TVerticalEdgeAdaptor* vSide0 = (TVerticalEdgeAdaptor*) VertiCurve(0);
2280   cout << "Verti side 0: "; vSide0->dumpNodes(nbNodes); cout << endl;
2281   TVerticalEdgeAdaptor* vSide1 = (TVerticalEdgeAdaptor*) VertiCurve(1);
2282   cout << "Verti side 1: "; vSide1->dumpNodes(nbNodes); cout << endl;
2283   delete hSize0; delete hSize1; delete vSide0; delete vSide1;
2284 #endif
2285 }
2286
2287 //================================================================================
2288 /*!
2289  * \brief Creates TVerticalEdgeAdaptor 
2290   * \param columnsMap - node column map
2291   * \param parameter - normalized parameter
2292  */
2293 //================================================================================
2294
2295 StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::
2296 TVerticalEdgeAdaptor( const TParam2ColumnMap* columnsMap, const double parameter)
2297 {
2298   myNodeColumn = & getColumn( columnsMap, parameter )->second;
2299 }
2300
2301 //================================================================================
2302 /*!
2303  * \brief Return coordinates for the given normalized parameter
2304   * \param U - normalized parameter
2305   * \retval gp_Pnt - coordinates
2306  */
2307 //================================================================================
2308
2309 gp_Pnt StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::Value(const Standard_Real U) const
2310 {
2311   const SMDS_MeshNode* n1;
2312   const SMDS_MeshNode* n2;
2313   double r = getRAndNodes( myNodeColumn, U, n1, n2 );
2314   return gpXYZ(n1) * ( 1 - r ) + gpXYZ(n2) * r;
2315 }
2316
2317 //================================================================================
2318 /*!
2319  * \brief Dump ids of nodes
2320  */
2321 //================================================================================
2322
2323 void StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::dumpNodes(int nbNodes) const
2324 {
2325 #ifdef _DEBUG_
2326   for ( int i = 0; i < nbNodes && i < myNodeColumn->size(); ++i )
2327     cout << (*myNodeColumn)[i]->GetID() << " ";
2328   if ( nbNodes < myNodeColumn->size() )
2329     cout << myNodeColumn->back()->GetID();
2330 #endif
2331 }
2332
2333 //================================================================================
2334 /*!
2335  * \brief Return coordinates for the given normalized parameter
2336   * \param U - normalized parameter
2337   * \retval gp_Pnt - coordinates
2338  */
2339 //================================================================================
2340
2341 gp_Pnt StdMeshers_PrismAsBlock::THorizontalEdgeAdaptor::Value(const Standard_Real U) const
2342 {
2343   return mySide->TSideFace::Value( U, myV );
2344 }
2345
2346 //================================================================================
2347 /*!
2348  * \brief Dump ids of <nbNodes> first nodes and the last one
2349  */
2350 //================================================================================
2351
2352 void StdMeshers_PrismAsBlock::THorizontalEdgeAdaptor::dumpNodes(int nbNodes) const
2353 {
2354 #ifdef _DEBUG_
2355   // Not bedugged code. Last node is sometimes incorrect
2356   const TSideFace* side = mySide;
2357   double u = 0;
2358   if ( mySide->IsComplex() )
2359     side = mySide->GetComponent(0,u);
2360
2361   TParam2ColumnIt col, col2;
2362   TParam2ColumnMap* u2cols = side->GetColumns();
2363   side->GetColumns( u , col, col2 );
2364   
2365   int j, i = myV ? mySide->ColumnHeight()-1 : 0;
2366
2367   const SMDS_MeshNode* n = 0;
2368   const SMDS_MeshNode* lastN
2369     = side->IsForward() ? u2cols->rbegin()->second[ i ] : u2cols->begin()->second[ i ];
2370   for ( j = 0; j < nbNodes && n != lastN; ++j )
2371   {
2372     n = col->second[ i ];
2373     cout << n->GetID() << " ";
2374     if ( side->IsForward() )
2375       ++col;
2376     else
2377       --col;
2378   }
2379
2380   // last node
2381   u = 1;
2382   if ( mySide->IsComplex() )
2383     side = mySide->GetComponent(1,u);
2384
2385   side->GetColumns( u , col, col2 );
2386   if ( n != col->second[ i ] )
2387     cout << col->second[ i ]->GetID();
2388 #endif
2389 }
2390 //================================================================================
2391 /*!
2392  * \brief Return UV on pcurve for the given normalized parameter
2393   * \param U - normalized parameter
2394   * \retval gp_Pnt - coordinates
2395  */
2396 //================================================================================
2397
2398 gp_Pnt2d StdMeshers_PrismAsBlock::TPCurveOnHorFaceAdaptor::Value(const Standard_Real U) const
2399 {
2400   TParam2ColumnIt u_col1, u_col2;
2401   double r = mySide->GetColumns( U, u_col1, u_col2 );
2402   gp_XY uv1 = mySide->GetNodeUV( myFace, u_col1->second[ myZ ]);
2403   gp_XY uv2 = mySide->GetNodeUV( myFace, u_col2->second[ myZ ]);
2404   return uv1 * ( 1 - r ) + uv2 * r;
2405 }