Salome HOME
Merge from V5_1_main 14/05/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()->GetPosition()->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     if ( myBlock.HasNotQuadElemOnTop() )
810       return error(TCom("Mesh on faces #") << botSM->GetId()
811                    <<" and #"<< topSM->GetId() << " seems different" );
812     needProject = true;
813   }
814
815   if ( 0/*needProject && !myProjectTriangles*/ )
816     return error(TCom("Mesh on faces #") << botSM->GetId()
817                  <<" and #"<< topSM->GetId() << " seems different" );
818   ///RETURN_BAD_RESULT("Need to project but not allowed");
819
820   if ( needProject )
821   {
822     return projectBottomToTop();
823   }
824
825   TopoDS_Face botFace = TopoDS::Face( myBlock.Shape( ID_BOT_FACE ));
826   TopoDS_Face topFace = TopoDS::Face( myBlock.Shape( ID_TOP_FACE ));
827   // associate top and bottom faces
828   TAssocTool::TShapeShapeMap shape2ShapeMap;
829   if ( !TAssocTool::FindSubShapeAssociation( botFace, myBlock.Mesh(),
830                                              topFace, myBlock.Mesh(),
831                                              shape2ShapeMap) )
832     return error(TCom("Topology of faces #") << botSM->GetId()
833                  <<" and #"<< topSM->GetId() << " seems different" );
834
835   // Find matching nodes of top and bottom faces
836   TNodeNodeMap n2nMap;
837   if ( ! TAssocTool::FindMatchingNodesOnFaces( botFace, myBlock.Mesh(),
838                                                topFace, myBlock.Mesh(),
839                                                shape2ShapeMap, n2nMap ))
840     return error(TCom("Mesh on faces #") << botSM->GetId()
841                  <<" and #"<< topSM->GetId() << " seems different" );
842
843   // Fill myBotToColumnMap
844
845   int zSize = myBlock.VerticalSize();
846   //TNode prevTNode;
847   TNodeNodeMap::iterator bN_tN = n2nMap.begin();
848   for ( ; bN_tN != n2nMap.end(); ++bN_tN )
849   {
850     const SMDS_MeshNode* botNode = bN_tN->first;
851     const SMDS_MeshNode* topNode = bN_tN->second;
852     if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
853       continue; // wall columns are contained in myBlock
854     // create node column
855     TNode bN( botNode );
856     TNode2ColumnMap::iterator bN_col = 
857       myBotToColumnMap.insert( make_pair ( bN, TNodeColumn() )).first;
858     TNodeColumn & column = bN_col->second;
859     column.resize( zSize );
860     column.front() = botNode;
861     column.back()  = topNode;
862   }
863   return true;
864 }
865
866 //================================================================================
867 /*!
868  * \brief Remove quadrangles from the top face and
869  * create triangles there by projection from the bottom
870  * \retval bool - a success or not
871  */
872 //================================================================================
873
874 bool StdMeshers_Prism_3D::projectBottomToTop()
875 {
876   SMESH_subMesh * botSM = myBlock.SubMesh( ID_BOT_FACE );
877   SMESH_subMesh * topSM = myBlock.SubMesh( ID_TOP_FACE );
878
879   SMESHDS_SubMesh * botSMDS = botSM->GetSubMeshDS();
880   SMESHDS_SubMesh * topSMDS = topSM->GetSubMeshDS();
881
882   if ( topSMDS )
883     topSM->ComputeStateEngine( SMESH_subMesh::CLEAN );
884
885   SMESHDS_Mesh* meshDS = myBlock.MeshDS();
886   int shapeID = myHelper->GetSubShapeID();
887   int topFaceID = meshDS->ShapeToIndex( topSM->GetSubShape() );
888
889   // Fill myBotToColumnMap
890
891   int zSize = myBlock.VerticalSize();
892   TNode prevTNode;
893   SMDS_NodeIteratorPtr nIt = botSMDS->GetNodes();
894   while ( nIt->more() )
895   {
896     const SMDS_MeshNode* botNode = nIt->next();
897     if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
898       continue; // strange
899     // compute bottom node params
900     TNode bN( botNode );
901     gp_XYZ paramHint(-1,-1,-1);
902     if ( prevTNode.IsNeighbor( bN ))
903       paramHint = prevTNode.GetParams();
904     if ( !myBlock.ComputeParameters( bN.GetCoords(), bN.ChangeParams(),
905                                      ID_BOT_FACE, paramHint ))
906       return error(TCom("Can't compute normalized parameters for node ")
907                    << botNode->GetID() << " on the face #"<< botSM->GetId() );
908     prevTNode = bN;
909     // compute top node coords
910     gp_XYZ topXYZ; gp_XY topUV;
911     if ( !myBlock.FacePoint( ID_TOP_FACE, bN.GetParams(), topXYZ ) ||
912          !myBlock.FaceUV   ( ID_TOP_FACE, bN.GetParams(), topUV ))
913       return error(TCom("Can't compute coordinates "
914                         "by normalized parameters on the face #")<< topSM->GetId() );
915     SMDS_MeshNode * topNode = meshDS->AddNode( topXYZ.X(),topXYZ.Y(),topXYZ.Z() );
916     meshDS->SetNodeOnFace( topNode, topFaceID, topUV.X(), topUV.Y() );
917     // create node column
918     TNode2ColumnMap::iterator bN_col = 
919       myBotToColumnMap.insert( make_pair ( bN, TNodeColumn() )).first;
920     TNodeColumn & column = bN_col->second;
921     column.resize( zSize );
922     column.front() = botNode;
923     column.back()  = topNode;
924   }
925
926   // Create top faces
927
928   // loop on bottom mesh faces
929   SMDS_ElemIteratorPtr faceIt = botSMDS->GetElements();
930   while ( faceIt->more() )
931   {
932     const SMDS_MeshElement* face = faceIt->next();
933     if ( !face || face->GetType() != SMDSAbs_Face )
934       continue;
935     int nbNodes = face->NbNodes();
936     if ( face->IsQuadratic() )
937       nbNodes /= 2;
938
939     // find top node in columns for each bottom node
940     vector< const SMDS_MeshNode* > nodes( nbNodes );
941     for ( int i = 0; i < nbNodes; ++i )
942     {
943       const SMDS_MeshNode* n = face->GetNode( nbNodes - i - 1 );
944       if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ) {
945         TNode2ColumnMap::iterator bot_column = myBotToColumnMap.find( n );
946         if ( bot_column == myBotToColumnMap.end() )
947           return error(TCom("No nodes found above node ") << n->GetID() );
948         nodes[ i ] = bot_column->second.back();
949       }
950       else {
951         const TNodeColumn* column = myBlock.GetNodeColumn( n );
952         if ( !column )
953           return error(TCom("No side nodes found above node ") << n->GetID() );
954         nodes[ i ] = column->back();
955       }
956     }
957     // create a face, with reversed orientation
958     SMDS_MeshElement* newFace = 0;
959     switch ( nbNodes ) {
960
961     case 3: {
962       newFace = myHelper->AddFace(nodes[0], nodes[1], nodes[2]);
963       break;
964       }
965     case 4: {
966       newFace = myHelper->AddFace( nodes[0], nodes[1], nodes[2], nodes[3] );
967       break;
968       }
969     default:
970       newFace = meshDS->AddPolygonalFace( nodes );
971     }
972     if ( newFace && shapeID > 0 )
973       meshDS->SetMeshElementOnShape( newFace, shapeID );
974   }
975
976   return true;
977 }
978
979 //================================================================================
980 /*!
981  * \brief Set projection coordinates of a node to a face and it's subshapes
982  * \param faceID - the face given by in-block ID
983  * \param params - node normalized parameters
984  * \retval bool - is a success
985  */
986 //================================================================================
987
988 bool StdMeshers_Prism_3D::setFaceAndEdgesXYZ( const int faceID, const gp_XYZ& params, int z )
989 {
990   // find base and top edges of the face
991   enum { BASE = 0, TOP, LEFT, RIGHT };
992   vector< int > edgeVec; // 0-base, 1-top
993   SMESH_Block::GetFaceEdgesIDs( faceID, edgeVec );
994
995   myBlock.EdgePoint( edgeVec[ BASE ], params, myShapeXYZ[ edgeVec[ BASE ]]);
996   myBlock.EdgePoint( edgeVec[ TOP  ], params, myShapeXYZ[ edgeVec[ TOP ]]);
997
998   SHOWYXZ("\nparams ", params);
999   SHOWYXZ("TOP is " <<edgeVec[ TOP ], myShapeXYZ[ edgeVec[ TOP]]);
1000   SHOWYXZ("BASE is "<<edgeVec[ BASE], myShapeXYZ[ edgeVec[ BASE]]);
1001
1002   if ( faceID == SMESH_Block::ID_Fx0z || faceID == SMESH_Block::ID_Fx1z )
1003   {
1004     myBlock.EdgePoint( edgeVec[ LEFT ], params, myShapeXYZ[ edgeVec[ LEFT ]]);
1005     myBlock.EdgePoint( edgeVec[ RIGHT ], params, myShapeXYZ[ edgeVec[ RIGHT ]]);
1006
1007     SHOWYXZ("VER "<<edgeVec[ LEFT], myShapeXYZ[ edgeVec[ LEFT]]);
1008     SHOWYXZ("VER "<<edgeVec[ RIGHT], myShapeXYZ[ edgeVec[ RIGHT]]);
1009   }
1010   myBlock.FacePoint( faceID, params, myShapeXYZ[ faceID ]);
1011   SHOWYXZ("FacePoint "<<faceID, myShapeXYZ[ faceID]);
1012
1013   return true;
1014 }
1015
1016 //================================================================================
1017 /*!
1018  * \brief Return true if this node and other one belong to one face
1019  */
1020 //================================================================================
1021
1022 bool TNode::IsNeighbor( const TNode& other ) const
1023 {
1024   if ( !other.myNode || !myNode ) return false;
1025
1026   SMDS_ElemIteratorPtr fIt = other.myNode->GetInverseElementIterator(SMDSAbs_Face);
1027   while ( fIt->more() )
1028     if ( fIt->next()->GetNodeIndex( myNode ) >= 0 )
1029       return true;
1030   return false;
1031 }
1032
1033 //================================================================================
1034 /*!
1035  * \brief Constructor. Initialization is needed
1036  */
1037 //================================================================================
1038
1039 StdMeshers_PrismAsBlock::StdMeshers_PrismAsBlock()
1040 {
1041   mySide = 0;
1042 }
1043
1044 StdMeshers_PrismAsBlock::~StdMeshers_PrismAsBlock()
1045 {
1046   if ( mySide ) {
1047     delete mySide; mySide = 0;
1048   }
1049 }
1050
1051 //================================================================================
1052 /*!
1053  * \brief Initialization.
1054  * \param helper - helper loaded with mesh and 3D shape
1055  * \param shape3D - a closed shell or solid
1056  * \retval bool - false if a mesh or a shape are KO
1057  */
1058 //================================================================================
1059
1060 bool StdMeshers_PrismAsBlock::Init(SMESH_MesherHelper* helper,
1061                                    const TopoDS_Shape& shape3D)
1062 {
1063   if ( mySide ) {
1064     delete mySide; mySide = 0;
1065   }
1066   vector< TSideFace* > sideFaces( NB_WALL_FACES, 0 );
1067   vector< pair< double, double> > params ( NB_WALL_FACES );
1068   mySide = new TSideFace( sideFaces, params );
1069
1070   myHelper = helper;
1071   SMESHDS_Mesh* meshDS = myHelper->GetMeshDS();
1072
1073   SMESH_Block::init();
1074   myShapeIDMap.Clear();
1075   myShapeIndex2ColumnMap.clear();
1076   
1077   int wallFaceIds[ NB_WALL_FACES ] = { // to walk around a block
1078     SMESH_Block::ID_Fx0z, SMESH_Block::ID_F1yz,
1079     SMESH_Block::ID_Fx1z, SMESH_Block::ID_F0yz
1080   };
1081
1082   myError = SMESH_ComputeError::New();
1083
1084   // -------------------------------------------------------------
1085   // Look for top and bottom faces: not quadrangle ones or meshed
1086   // with not quadrangle elements
1087   // -------------------------------------------------------------
1088
1089   list< SMESH_subMesh* > notQuadGeomSubMesh;
1090   list< SMESH_subMesh* > notQuadElemSubMesh;
1091   int nbFaces = 0;
1092   //
1093   SMESH_subMesh* mainSubMesh = myHelper->GetMesh()->GetSubMeshContaining( shape3D );
1094   if ( !mainSubMesh ) return error(COMPERR_BAD_INPUT_MESH,"Null submesh of shape3D");
1095
1096   // analyse face submeshes
1097   SMESH_subMeshIteratorPtr smIt = mainSubMesh->getDependsOnIterator(false,false);
1098   while ( smIt->more() )
1099   {
1100     SMESH_subMesh* sm = smIt->next();
1101     const TopoDS_Shape& face = sm->GetSubShape();
1102     if ( face.ShapeType() != TopAbs_FACE )
1103       continue;
1104     nbFaces++;
1105
1106     // is quadrangle face?
1107     list< TopoDS_Edge > orderedEdges;
1108     list< int >         nbEdgesInWires;
1109     TopoDS_Vertex       V000;
1110     int nbWires = GetOrderedEdges( TopoDS::Face( face ),
1111                                    V000, orderedEdges, nbEdgesInWires );
1112     if ( nbWires != 1 || nbEdgesInWires.front() != 4 )
1113       notQuadGeomSubMesh.push_back( sm );
1114
1115     // look for not quadrangle mesh elements
1116     if ( SMESHDS_SubMesh* smDS = sm->GetSubMeshDS() ) {
1117       bool hasNotQuad = false;
1118       SMDS_ElemIteratorPtr eIt = smDS->GetElements();
1119       while ( eIt->more() && !hasNotQuad ) {
1120         const SMDS_MeshElement* elem = eIt->next();
1121         if ( elem->GetType() == SMDSAbs_Face ) {
1122           int nbNodes = elem->NbNodes();
1123           if ( elem->IsQuadratic() )
1124             nbNodes /= 2;
1125           hasNotQuad = ( nbNodes != 4 );
1126         }
1127       }
1128       if ( hasNotQuad )
1129         notQuadElemSubMesh.push_back( sm );
1130     }
1131     else {
1132       return error(COMPERR_BAD_INPUT_MESH,TCom("Not meshed face #")<<sm->GetId());
1133     }
1134     // check if a quadrangle face is meshed with a quadranglar grid
1135     if ( notQuadGeomSubMesh.back() != sm &&
1136          notQuadElemSubMesh.back() != sm )
1137     {
1138       // count nb edges on face sides
1139       vector< int > nbEdges;
1140       nbEdges.reserve( nbEdgesInWires.front() );
1141       for ( list< TopoDS_Edge >::iterator edge = orderedEdges.begin();
1142             edge != orderedEdges.end(); ++edge )
1143       {
1144         if ( SMESHDS_SubMesh* smDS = meshDS->MeshElements( *edge ))
1145           nbEdges.push_back ( smDS->NbElements() );
1146         else
1147           nbEdges.push_back ( 0 );
1148       }
1149       int nbQuads = sm->GetSubMeshDS()->NbElements();
1150       if ( nbEdges[0] *  nbEdges[1] != nbQuads ||
1151            nbEdges[0] != nbEdges[2] ||
1152            nbEdges[1] != nbEdges[3] )
1153         notQuadElemSubMesh.push_back( sm );
1154     }
1155   }
1156
1157   // ----------------------------------------------------------------------
1158   // Analyse mesh and topology of faces: choose the bottom submesh.
1159   // If there are not quadrangle geom faces, they are top and bottom ones.
1160   // Not quadrangle geom faces must be only on top and bottom.
1161   // ----------------------------------------------------------------------
1162
1163   SMESH_subMesh * botSM = 0;
1164   SMESH_subMesh * topSM = 0;
1165
1166   int nbNotQuad       = notQuadGeomSubMesh.size();
1167   int nbNotQuadMeshed = notQuadElemSubMesh.size();
1168   bool hasNotQuad = ( nbNotQuad || nbNotQuadMeshed );
1169
1170   // detect bad cases
1171   if ( nbNotQuadMeshed > 2 )
1172   {
1173     return error(COMPERR_BAD_INPUT_MESH,
1174                  TCom("More than 2 faces with not quadrangle elements: ")
1175                  <<nbNotQuadMeshed);
1176   }
1177   int nbQuasiQuads = 0;
1178   if ( nbNotQuad > 0 && nbNotQuad != 2 )
1179   {
1180     // Issue 0020843 - one of side faces is quasi-quadrilateral.
1181     // Remove from notQuadGeomSubMesh faces meshed with regular grid
1182     nbQuasiQuads = removeQuasiQuads( notQuadGeomSubMesh );
1183     nbNotQuad -= nbQuasiQuads;
1184     if ( nbNotQuad > 0 && nbNotQuad != 2 )
1185       return error(COMPERR_BAD_SHAPE,
1186                    TCom("More than 2 not quadrilateral faces: ")
1187                    <<nbNotQuad);
1188   }
1189
1190   // get found submeshes
1191   if ( hasNotQuad )
1192   {
1193     if ( nbNotQuadMeshed > 0 ) botSM = notQuadElemSubMesh.front();
1194     else                       botSM = notQuadGeomSubMesh.front();
1195     if ( nbNotQuadMeshed > 1 ) topSM = notQuadElemSubMesh.back();
1196     else if ( nbNotQuad  > 1 ) topSM = notQuadGeomSubMesh.back();
1197   }
1198   // detect other bad cases
1199   if ( nbNotQuad == 2 && nbNotQuadMeshed > 0 ) {
1200     bool ok = false;
1201     if ( nbNotQuadMeshed == 1 )
1202       ok = ( find( notQuadGeomSubMesh.begin(),
1203                    notQuadGeomSubMesh.end(), botSM ) != notQuadGeomSubMesh.end() );
1204     else
1205       ok = ( notQuadGeomSubMesh == notQuadElemSubMesh );
1206     if ( !ok )
1207       return error(COMPERR_BAD_INPUT_MESH, "Side face meshed with not quadrangle elements");
1208   }
1209
1210   myNotQuadOnTop = ( nbNotQuadMeshed > 1 );
1211  
1212   // ----------------------------------------------------------
1213
1214   if ( nbNotQuad == 0 ) // Standard block of 6 quadrangle faces ?
1215   {
1216     // SMESH_Block will perform geometry analysis, we need just to find 2
1217     // connected vertices on top and bottom
1218
1219     TopoDS_Vertex Vbot, Vtop;
1220     if ( nbNotQuadMeshed > 0 ) // Look for vertices
1221     {
1222       TopTools_IndexedMapOfShape edgeMap;
1223       TopExp::MapShapes( botSM->GetSubShape(), TopAbs_EDGE, edgeMap );
1224       // vertex 1 is any vertex of the bottom face
1225       Vbot = TopExp::FirstVertex( TopoDS::Edge( edgeMap( 1 )));
1226       // vertex 2 is end vertex of edge sharing Vbot and not belonging to the bottom face
1227       TopTools_ListIteratorOfListOfShape ancestIt = Mesh()->GetAncestors( Vbot );
1228       for ( ; Vtop.IsNull() && ancestIt.More(); ancestIt.Next() )
1229       {
1230         const TopoDS_Shape & ancestor = ancestIt.Value();
1231         if ( ancestor.ShapeType() == TopAbs_EDGE && !edgeMap.FindIndex( ancestor ))
1232         {
1233           TopoDS_Vertex V1, V2;
1234           TopExp::Vertices( TopoDS::Edge( ancestor ), V1, V2);
1235           if      ( Vbot.IsSame ( V1 )) Vtop = V2;
1236           else if ( Vbot.IsSame ( V2 )) Vtop = V1;
1237           // check that Vtop belongs to shape3D
1238           TopExp_Explorer exp( shape3D, TopAbs_VERTEX );
1239           for ( ; exp.More(); exp.Next() )
1240             if ( Vtop.IsSame( exp.Current() ))
1241               break;
1242           if ( !exp.More() )
1243             Vtop.Nullify();
1244         }
1245       }
1246     }
1247     // get shell from shape3D
1248     TopoDS_Shell shell;
1249     TopExp_Explorer exp( shape3D, TopAbs_SHELL );
1250     int nbShell = 0;
1251     for ( ; exp.More(); exp.Next(), ++nbShell )
1252       shell = TopoDS::Shell( exp.Current() );
1253 //     if ( nbShell != 1 )
1254 //       RETURN_BAD_RESULT("There must be 1 shell in the block");
1255
1256     // Load geometry in SMESH_Block
1257     if ( !SMESH_Block::FindBlockShapes( shell, Vbot, Vtop, myShapeIDMap )) {
1258       if ( !hasNotQuad )
1259         return error(COMPERR_BAD_SHAPE, "Can't detect top and bottom of a prism");
1260     }
1261     else {
1262       if ( !botSM ) botSM = Mesh()->GetSubMeshContaining( myShapeIDMap( ID_BOT_FACE ));
1263       if ( !topSM ) topSM = Mesh()->GetSubMeshContaining( myShapeIDMap( ID_TOP_FACE ));
1264     }
1265
1266   } // end  Standard block of 6 quadrangle faces
1267   // --------------------------------------------------------
1268
1269   // Here the top and bottom faces are found
1270   if ( nbNotQuadMeshed == 2 ) // roughly check correspondence of horiz meshes
1271   {
1272 //     SMESHDS_SubMesh* topSMDS = topSM->GetSubMeshDS();
1273 //     SMESHDS_SubMesh* botSMDS = botSM->GetSubMeshDS();
1274 //     if ( topSMDS->NbNodes() != botSMDS->NbNodes() ||
1275 //          topSMDS->NbElements() != botSMDS->NbElements() )
1276 //       RETURN_BAD_RESULT("Top mesh doesn't correspond to bottom one");
1277   }
1278
1279   // ---------------------------------------------------------
1280   // If there are not quadrangle geom faces, we emulate
1281   // a block of 6 quadrangle faces.
1282   // Load SMESH_Block with faces and edges geometry
1283   // ---------------------------------------------------------
1284
1285   
1286   // find vertex 000 - the one with smallest coordinates (for easy DEBUG :-)
1287   TopoDS_Vertex V000;
1288   double minVal = DBL_MAX, minX, val;
1289   for ( TopExp_Explorer exp( botSM->GetSubShape(), TopAbs_VERTEX );
1290         exp.More(); exp.Next() )
1291   {
1292     const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
1293     gp_Pnt P = BRep_Tool::Pnt( v );
1294     val = P.X() + P.Y() + P.Z();
1295     if ( val < minVal || ( val == minVal && P.X() < minX )) {
1296       V000 = v;
1297       minVal = val;
1298       minX = P.X();
1299     }
1300   }
1301
1302   // Get ordered bottom edges
1303   list< TopoDS_Edge > orderedEdges;
1304   list< int >         nbVertexInWires;
1305   SMESH_Block::GetOrderedEdges( TopoDS::Face( botSM->GetSubShape().Reversed() ),
1306                                 V000, orderedEdges, nbVertexInWires );
1307 //   if ( nbVertexInWires.size() != 1 )
1308 //     RETURN_BAD_RESULT("Wrong prism geometry");
1309
1310   // Get Wall faces corresponding to the ordered bottom edges
1311   list< TopoDS_Face > wallFaces;
1312   if ( !GetWallFaces( Mesh(), shape3D, botSM->GetSubShape(), orderedEdges, wallFaces))
1313     return error(COMPERR_BAD_SHAPE, "Can't find side faces");
1314
1315   // Find columns of wall nodes and calculate edges' lengths
1316   // --------------------------------------------------------
1317
1318   myParam2ColumnMaps.clear();
1319   myParam2ColumnMaps.resize( orderedEdges.size() ); // total nb edges
1320
1321   int iE, nbEdges = nbVertexInWires.front(); // nb outer edges
1322   vector< double > edgeLength( nbEdges );
1323   map< double, int > len2edgeMap;
1324
1325   list< TopoDS_Edge >::iterator edgeIt = orderedEdges.begin();
1326   list< TopoDS_Face >::iterator faceIt = wallFaces.begin();
1327   for ( iE = 0; iE < nbEdges; ++edgeIt, ++faceIt )
1328   {
1329     TParam2ColumnMap & faceColumns = myParam2ColumnMaps[ iE ];
1330     if ( !myHelper->LoadNodeColumns( faceColumns, *faceIt, *edgeIt, meshDS ))
1331       return error(COMPERR_BAD_INPUT_MESH, TCom("Can't find regular quadrangle mesh ")
1332                    << "on a side face #" << MeshDS()->ShapeToIndex( *faceIt ));
1333
1334     SHOWYXZ("\np1 F "<<iE, gpXYZ(faceColumns.begin()->second.front() ));
1335     SHOWYXZ("p2 F "<<iE, gpXYZ(faceColumns.rbegin()->second.front() ));
1336     SHOWYXZ("V First "<<iE, BRep_Tool::Pnt( TopExp::FirstVertex(*edgeIt,true )));
1337
1338     edgeLength[ iE ] = SMESH_Algo::EdgeLength( *edgeIt );
1339
1340     if ( nbEdges < NB_WALL_FACES ) // fill map used to split faces
1341     {
1342       SMESHDS_SubMesh* smDS = meshDS->MeshElements( *edgeIt);
1343       if ( !smDS )
1344         return error(COMPERR_BAD_INPUT_MESH, TCom("Null submesh on the edge #")
1345                      << MeshDS()->ShapeToIndex( *edgeIt ));
1346       // assure length uniqueness
1347       edgeLength[ iE ] *= smDS->NbNodes() + edgeLength[ iE ] / ( 1000 + iE );
1348       len2edgeMap[ edgeLength[ iE ]] = iE;
1349     }
1350     ++iE;
1351   }
1352   // Load columns of internal edges (forming holes)
1353   // and fill map ShapeIndex to TParam2ColumnMap for them
1354   for ( ; edgeIt != orderedEdges.end() ; ++edgeIt, ++faceIt )
1355   {
1356     TParam2ColumnMap & faceColumns = myParam2ColumnMaps[ iE ];
1357     if ( !myHelper->LoadNodeColumns( faceColumns, *faceIt, *edgeIt, meshDS ))
1358       return error(COMPERR_BAD_INPUT_MESH, TCom("Can't find regular quadrangle mesh ")
1359                    << "on a side face #" << MeshDS()->ShapeToIndex( *faceIt ));
1360     // edge columns
1361     int id = MeshDS()->ShapeToIndex( *edgeIt );
1362     bool isForward = true; // meaningless for intenal wires
1363     myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1364     // columns for vertices
1365     // 1
1366     const SMDS_MeshNode* n0 = faceColumns.begin()->second.front();
1367     id = n0->GetPosition()->GetShapeId();
1368     myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1369     // 2
1370     const SMDS_MeshNode* n1 = faceColumns.rbegin()->second.front();
1371     id = n1->GetPosition()->GetShapeId();
1372     myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1373 //     SHOWYXZ("\np1 F "<<iE, gpXYZ(faceColumns.begin()->second.front() ));
1374 //     SHOWYXZ("p2 F "<<iE, gpXYZ(faceColumns.rbegin()->second.front() ));
1375 //     SHOWYXZ("V First "<<iE, BRep_Tool::Pnt( TopExp::FirstVertex(*edgeIt,true )));
1376     ++iE;
1377   }
1378
1379   // Create 4 wall faces of a block
1380   // -------------------------------
1381
1382   if ( nbEdges <= NB_WALL_FACES ) // ************* Split faces if necessary
1383   {
1384     map< int, int > iE2nbSplit;
1385     if ( nbEdges != NB_WALL_FACES ) // define how to split
1386     {
1387       if ( len2edgeMap.size() != nbEdges )
1388         RETURN_BAD_RESULT("Uniqueness of edge lengths not assured");
1389       map< double, int >::reverse_iterator maxLen_i = len2edgeMap.rbegin();
1390       map< double, int >::reverse_iterator midLen_i = ++len2edgeMap.rbegin();
1391       double maxLen = maxLen_i->first;
1392       double midLen = ( len2edgeMap.size() == 1 ) ? 0 : midLen_i->first;
1393       switch ( nbEdges ) {
1394       case 1: // 0-th edge is split into 4 parts
1395         iE2nbSplit.insert( make_pair( 0, 4 )); break;
1396       case 2: // either the longest edge is split into 3 parts, or both edges into halves
1397         if ( maxLen / 3 > midLen / 2 ) {
1398           iE2nbSplit.insert( make_pair( maxLen_i->second, 3 ));
1399         }
1400         else {
1401           iE2nbSplit.insert( make_pair( maxLen_i->second, 2 ));
1402           iE2nbSplit.insert( make_pair( midLen_i->second, 2 ));
1403         }
1404         break;
1405       case 3:
1406         // split longest into halves
1407         iE2nbSplit.insert( make_pair( maxLen_i->second, 2 ));
1408       }
1409     }
1410     // Create TSideFace's
1411     faceIt = wallFaces.begin();
1412     edgeIt = orderedEdges.begin();
1413     int iSide = 0;
1414     for ( iE = 0; iE < nbEdges; ++edgeIt, ++faceIt )
1415     {
1416      // split?
1417       map< int, int >::iterator i_nb = iE2nbSplit.find( iE );
1418       if ( i_nb != iE2nbSplit.end() ) {
1419         // split!
1420         int nbSplit = i_nb->second;
1421         vector< double > params;
1422         splitParams( nbSplit, &myParam2ColumnMaps[ iE ], params );
1423         bool isForward = ( edgeIt->Orientation() == TopAbs_FORWARD );
1424         for ( int i = 0; i < nbSplit; ++i ) {
1425           double f = ( isForward ? params[ i ] : params[ nbSplit - i-1 ]);
1426           double l = ( isForward ? params[ i+1 ] : params[ nbSplit - i ]);
1427           TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1428                                            *faceIt, *edgeIt,
1429                                            &myParam2ColumnMaps[ iE ], f, l );
1430           mySide->SetComponent( iSide++, comp );
1431         }
1432       }
1433       else {
1434         TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1435                                          *faceIt, *edgeIt,
1436                                          &myParam2ColumnMaps[ iE ]);
1437         mySide->SetComponent( iSide++, comp );
1438       }
1439       ++iE;
1440     }
1441   }
1442   else { // **************************** Unite faces
1443
1444     // unite first faces
1445     int nbExraFaces = nbEdges - 3;
1446     int iSide = 0, iE;
1447     double u0 = 0, sumLen = 0;
1448     for ( iE = 0; iE < nbExraFaces; ++iE )
1449       sumLen += edgeLength[ iE ];
1450
1451     vector< TSideFace* > components( nbExraFaces );
1452     vector< pair< double, double> > params( nbExraFaces );
1453     faceIt = wallFaces.begin();
1454     edgeIt = orderedEdges.begin();
1455     for ( iE = 0; iE < nbExraFaces; ++edgeIt, ++faceIt )
1456     {
1457       components[ iE ] = new TSideFace( myHelper, wallFaceIds[ iSide ],
1458                                         *faceIt, *edgeIt,
1459                                         &myParam2ColumnMaps[ iE ]);
1460       double u1 = u0 + edgeLength[ iE ] / sumLen;
1461       params[ iE ] = make_pair( u0 , u1 );
1462       u0 = u1;
1463       ++iE;
1464     }
1465     mySide->SetComponent( iSide++, new TSideFace( components, params ));
1466
1467     // fill the rest faces
1468     for ( ; iE < nbEdges; ++faceIt, ++edgeIt )
1469     {
1470       TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1471                                        *faceIt, *edgeIt,
1472                                        &myParam2ColumnMaps[ iE ]);
1473       mySide->SetComponent( iSide++, comp );
1474       ++iE;
1475     }
1476   }
1477
1478
1479   // Fill geometry fields of SMESH_Block
1480   // ------------------------------------
1481
1482   TopoDS_Face botF = TopoDS::Face( botSM->GetSubShape() );
1483   TopoDS_Face topF = TopoDS::Face( topSM->GetSubShape() );
1484
1485   vector< int > botEdgeIdVec;
1486   SMESH_Block::GetFaceEdgesIDs( ID_BOT_FACE, botEdgeIdVec );
1487
1488   bool isForward[NB_WALL_FACES] = { true, true, true, true };
1489   Adaptor2d_Curve2d* botPcurves[NB_WALL_FACES];
1490   Adaptor2d_Curve2d* topPcurves[NB_WALL_FACES];
1491
1492   for ( int iF = 0; iF < NB_WALL_FACES; ++iF )
1493   {
1494     TSideFace * sideFace = mySide->GetComponent( iF );
1495     if ( !sideFace )
1496       RETURN_BAD_RESULT("NULL TSideFace");
1497     int fID = sideFace->FaceID();
1498
1499     // fill myShapeIDMap
1500     if ( sideFace->InsertSubShapes( myShapeIDMap ) != 8 &&
1501          !sideFace->IsComplex())
1502       MESSAGE( ": Warning : InsertSubShapes() < 8 on side " << iF );
1503
1504     // side faces geometry
1505     Adaptor2d_Curve2d* pcurves[NB_WALL_FACES];
1506     if ( !sideFace->GetPCurves( pcurves ))
1507       RETURN_BAD_RESULT("TSideFace::GetPCurves() failed");
1508
1509     SMESH_Block::TFace& tFace = myFace[ fID - ID_FirstF ];
1510     tFace.Set( fID, sideFace->Surface(), pcurves, isForward );
1511
1512     SHOWYXZ( endl<<"F "<< iF << " id " << fID << " FRW " << sideFace->IsForward(), sideFace->Value(0,0));
1513     // edges 3D geometry
1514     vector< int > edgeIdVec;
1515     SMESH_Block::GetFaceEdgesIDs( fID, edgeIdVec );
1516     for ( int isMax = 0; isMax < 2; ++isMax ) {
1517       {
1518         int eID = edgeIdVec[ isMax ];
1519         SMESH_Block::TEdge& tEdge = myEdge[ eID - ID_FirstE ];
1520         tEdge.Set( eID, sideFace->HorizCurve(isMax), true);
1521         SHOWYXZ(eID<<" HOR"<<isMax<<"(0)", sideFace->HorizCurve(isMax)->Value(0));
1522         SHOWYXZ(eID<<" HOR"<<isMax<<"(1)", sideFace->HorizCurve(isMax)->Value(1));
1523       }
1524       {
1525         int eID = edgeIdVec[ isMax+2 ];
1526         SMESH_Block::TEdge& tEdge = myEdge[ eID - ID_FirstE  ];
1527         tEdge.Set( eID, sideFace->VertiCurve(isMax), true);
1528         SHOWYXZ(eID<<" VER"<<isMax<<"(0)", sideFace->VertiCurve(isMax)->Value(0));
1529         SHOWYXZ(eID<<" VER"<<isMax<<"(1)", sideFace->VertiCurve(isMax)->Value(1));
1530
1531         // corner points
1532         vector< int > vertexIdVec;
1533         SMESH_Block::GetEdgeVertexIDs( eID, vertexIdVec );
1534         myPnt[ vertexIdVec[0] - ID_FirstV ] = tEdge.GetCurve()->Value(0).XYZ();
1535         myPnt[ vertexIdVec[1] - ID_FirstV ] = tEdge.GetCurve()->Value(1).XYZ();
1536       }
1537     }
1538     // pcurves on horizontal faces
1539     for ( iE = 0; iE < NB_WALL_FACES; ++iE ) {
1540       if ( edgeIdVec[ BOTTOM_EDGE ] == botEdgeIdVec[ iE ] ) {
1541         botPcurves[ iE ] = sideFace->HorizPCurve( false, botF );
1542         topPcurves[ iE ] = sideFace->HorizPCurve( true,  topF );
1543         break;
1544       }
1545     }
1546     //sideFace->dumpNodes( 4 ); // debug
1547   }
1548   // horizontal faces geometry
1549   {
1550     SMESH_Block::TFace& tFace = myFace[ ID_BOT_FACE - ID_FirstF ];
1551     tFace.Set( ID_BOT_FACE, new BRepAdaptor_Surface( botF ), botPcurves, isForward );
1552     SMESH_Block::Insert( botF, ID_BOT_FACE, myShapeIDMap );
1553   }
1554   {
1555     SMESH_Block::TFace& tFace = myFace[ ID_TOP_FACE - ID_FirstF ];
1556     tFace.Set( ID_TOP_FACE, new BRepAdaptor_Surface( topF ), topPcurves, isForward );
1557     SMESH_Block::Insert( topF, ID_TOP_FACE, myShapeIDMap );
1558   }
1559
1560   // Fill map ShapeIndex to TParam2ColumnMap
1561   // ----------------------------------------
1562
1563   list< TSideFace* > fList;
1564   list< TSideFace* >::iterator fListIt;
1565   fList.push_back( mySide );
1566   for ( fListIt = fList.begin(); fListIt != fList.end(); ++fListIt)
1567   {
1568     int nb = (*fListIt)->NbComponents();
1569     for ( int i = 0; i < nb; ++i ) {
1570       if ( TSideFace* comp = (*fListIt)->GetComponent( i ))
1571         fList.push_back( comp );
1572     }
1573     if ( TParam2ColumnMap* cols = (*fListIt)->GetColumns()) {
1574       // columns for a base edge
1575       int id = MeshDS()->ShapeToIndex( (*fListIt)->BaseEdge() );
1576       bool isForward = (*fListIt)->IsForward();
1577       myShapeIndex2ColumnMap[ id ] = make_pair( cols, isForward );
1578
1579       // columns for vertices
1580       const SMDS_MeshNode* n0 = cols->begin()->second.front();
1581       id = n0->GetPosition()->GetShapeId();
1582       myShapeIndex2ColumnMap[ id ] = make_pair( cols, isForward );
1583
1584       const SMDS_MeshNode* n1 = cols->rbegin()->second.front();
1585       id = n1->GetPosition()->GetShapeId();
1586       myShapeIndex2ColumnMap[ id ] = make_pair( cols, !isForward );
1587     }
1588   }
1589
1590 //   gp_XYZ testPar(0.25, 0.25, 0), testCoord;
1591 //   if ( !FacePoint( ID_BOT_FACE, testPar, testCoord ))
1592 //     RETURN_BAD_RESULT("TEST FacePoint() FAILED");
1593 //   SHOWYXZ("IN TEST PARAM" , testPar);
1594 //   SHOWYXZ("OUT TEST CORD" , testCoord);
1595 //   if ( !ComputeParameters( testCoord, testPar , ID_BOT_FACE))
1596 //     RETURN_BAD_RESULT("TEST ComputeParameters() FAILED");
1597 //   SHOWYXZ("OUT TEST PARAM" , testPar);
1598
1599   return true;
1600 }
1601
1602 //================================================================================
1603 /*!
1604  * \brief Return pointer to column of nodes
1605  * \param node - bottom node from which the returned column goes up
1606  * \retval const TNodeColumn* - the found column
1607  */
1608 //================================================================================
1609
1610 const TNodeColumn* StdMeshers_PrismAsBlock::GetNodeColumn(const SMDS_MeshNode* node) const
1611 {
1612   int sID = node->GetPosition()->GetShapeId();
1613
1614   map<int, pair< TParam2ColumnMap*, bool > >::const_iterator col_frw =
1615     myShapeIndex2ColumnMap.find( sID );
1616   if ( col_frw != myShapeIndex2ColumnMap.end() ) {
1617     const TParam2ColumnMap* cols = col_frw->second.first;
1618     TParam2ColumnIt u_col = cols->begin();
1619     for ( ; u_col != cols->end(); ++u_col )
1620       if ( u_col->second[ 0 ] == node )
1621         return & u_col->second;
1622   }
1623   return 0;
1624 }
1625
1626 //=======================================================================
1627 //function : GetLayersTransformation
1628 //purpose  : Return transformations to get coordinates of nodes of each layer
1629 //           by nodes of the bottom. Layer is a set of nodes at a certain step
1630 //           from bottom to top.
1631 //=======================================================================
1632
1633 bool StdMeshers_PrismAsBlock::GetLayersTransformation(vector<gp_Trsf> & trsf) const
1634 {
1635   const int zSize = VerticalSize();
1636   if ( zSize < 3 ) return true;
1637   trsf.resize( zSize - 2 );
1638
1639   // Select some node columns by which we will define coordinate system of layers
1640
1641   vector< const TNodeColumn* > columns;
1642   {
1643     const TopoDS_Shape& baseFace = Shape(ID_BOT_FACE);
1644     list< TopoDS_Edge > orderedEdges;
1645     list< int >         nbEdgesInWires;
1646     GetOrderedEdges( TopoDS::Face( baseFace ), TopoDS_Vertex(), orderedEdges, nbEdgesInWires );
1647     bool isReverse;
1648     list< TopoDS_Edge >::iterator edgeIt = orderedEdges.begin();
1649     for ( int iE = 0; iE < nbEdgesInWires.front(); ++iE, ++edgeIt )
1650     {
1651       const TParam2ColumnMap& u2colMap =
1652         GetParam2ColumnMap( myHelper->GetMeshDS()->ShapeToIndex( *edgeIt ), isReverse );
1653       isReverse = ( edgeIt->Orientation() == TopAbs_REVERSED );
1654       double f = u2colMap.begin()->first, l = u2colMap.rbegin()->first;
1655       if ( isReverse ) swap ( f, l );
1656       const int nbCol = 5;
1657       for ( int i = 0; i < nbCol; ++i )
1658       {
1659         double u = f + i/double(nbCol) * ( l - f );
1660         const TNodeColumn* col = & getColumn( & u2colMap, u )->second;
1661         if ( columns.empty() || col != columns.back() )
1662           columns.push_back( col );
1663       }
1664     }
1665   }
1666
1667   // Find tolerance to check transformations
1668
1669   double tol2;
1670   {
1671     Bnd_B3d bndBox;
1672     for ( int i = 0; i < columns.size(); ++i )
1673       bndBox.Add( gpXYZ( columns[i]->front() ));
1674     tol2 = bndBox.SquareExtent() * 4 * 1e-4;
1675   }
1676
1677   // Compute transformations
1678
1679   int xCol = -1;
1680   gp_Trsf fromCsZ, toCs0;
1681   gp_Ax3 cs0 = getLayerCoordSys(0, columns, xCol );
1682   //double dist0 = cs0.Location().Distance( gpXYZ( (*columns[0])[0]));
1683   toCs0.SetTransformation( cs0 );
1684   for ( int z = 1; z < zSize-1; ++z )
1685   {
1686     gp_Ax3 csZ = getLayerCoordSys(z, columns, xCol );
1687     //double distZ = csZ.Location().Distance( gpXYZ( (*columns[0])[z]));
1688     fromCsZ.SetTransformation( csZ );
1689     fromCsZ.Invert();
1690     gp_Trsf& t = trsf[ z-1 ];
1691     t = fromCsZ * toCs0;
1692     //t.SetScaleFactor( distZ/dist0 ); - it does not work properly, wrong base point
1693
1694     // check a transformation
1695     for ( int i = 0; i < columns.size(); ++i )
1696     {
1697       gp_Pnt p0 = gpXYZ( (*columns[i])[0] );
1698       gp_Pnt pz = gpXYZ( (*columns[i])[z] );
1699       t.Transforms( p0.ChangeCoord() );
1700       if ( p0.SquareDistance( pz ) > tol2 )
1701         return false;
1702     }
1703   }
1704   return true;
1705 }
1706
1707 //================================================================================
1708 /*!
1709  * \brief Check curve orientation of a bootom edge
1710   * \param meshDS - mesh DS
1711   * \param columnsMap - node columns map of side face
1712   * \param bottomEdge - the bootom edge
1713   * \param sideFaceID - side face in-block ID
1714   * \retval bool - true if orientation coinside with in-block froward orientation
1715  */
1716 //================================================================================
1717
1718 bool StdMeshers_PrismAsBlock::IsForwardEdge(SMESHDS_Mesh*           meshDS,
1719                                             const TParam2ColumnMap& columnsMap,
1720                                             const TopoDS_Edge &     bottomEdge,
1721                                             const int               sideFaceID)
1722 {
1723   bool isForward = false;
1724   if ( TAssocTool::IsClosedEdge( bottomEdge ))
1725   {
1726     isForward = ( bottomEdge.Orientation() == TopAbs_FORWARD );
1727   }
1728   else
1729   {
1730     const TNodeColumn& firstCol = columnsMap.begin()->second;
1731     const SMDS_MeshNode* bottomNode = firstCol[0];
1732     TopoDS_Shape firstVertex = SMESH_MesherHelper::GetSubShapeByNode( bottomNode, meshDS );
1733     isForward = ( firstVertex.IsSame( TopExp::FirstVertex( bottomEdge, true )));
1734   }
1735   // on 2 of 4 sides first vertex is end
1736   if ( sideFaceID == ID_Fx1z || sideFaceID == ID_F0yz )
1737     isForward = !isForward;
1738   return isForward;
1739 }
1740
1741 //================================================================================
1742 /*!
1743  * \brief Find wall faces by bottom edges
1744  * \param mesh - the mesh
1745  * \param mainShape - the prism
1746  * \param bottomFace - the bottom face
1747  * \param bottomEdges - edges bounding the bottom face
1748  * \param wallFaces - faces list to fill in
1749  */
1750 //================================================================================
1751
1752 bool StdMeshers_PrismAsBlock::GetWallFaces( SMESH_Mesh*                     mesh,
1753                                             const TopoDS_Shape &            mainShape,
1754                                             const TopoDS_Shape &            bottomFace,
1755                                             const std::list< TopoDS_Edge >& bottomEdges,
1756                                             std::list< TopoDS_Face >&       wallFaces)
1757 {
1758   wallFaces.clear();
1759
1760   TopTools_IndexedMapOfShape faceMap;
1761   TopExp::MapShapes( mainShape, TopAbs_FACE, faceMap );
1762
1763   list< TopoDS_Edge >::const_iterator edge = bottomEdges.begin();
1764   for ( ; edge != bottomEdges.end(); ++edge )
1765   {
1766     TopTools_ListIteratorOfListOfShape ancestIt = mesh->GetAncestors( *edge );
1767     for ( ; ancestIt.More(); ancestIt.Next() )
1768     {
1769       const TopoDS_Shape& ancestor = ancestIt.Value();
1770       if ( ancestor.ShapeType() == TopAbs_FACE && // face
1771            !bottomFace.IsSame( ancestor ) &&      // not bottom
1772            faceMap.FindIndex( ancestor ))         // belongs to the prism
1773       {
1774         wallFaces.push_back( TopoDS::Face( ancestor ));
1775         break;
1776       }
1777     }
1778   }
1779   return ( wallFaces.size() == bottomEdges.size() );
1780 }
1781
1782 //================================================================================
1783 /*!
1784  * \brief Constructor
1785   * \param faceID - in-block ID
1786   * \param face - geom face
1787   * \param columnsMap - map of node columns
1788   * \param first - first normalized param
1789   * \param last - last normalized param
1790  */
1791 //================================================================================
1792
1793 StdMeshers_PrismAsBlock::TSideFace::TSideFace(SMESH_MesherHelper* helper,
1794                                               const int           faceID,
1795                                               const TopoDS_Face&  face,
1796                                               const TopoDS_Edge&  baseEdge,
1797                                               TParam2ColumnMap*   columnsMap,
1798                                               const double        first,
1799                                               const double        last):
1800   myID( faceID ),
1801   myParamToColumnMap( columnsMap ),
1802   myBaseEdge( baseEdge ),
1803   myHelper( helper )
1804 {
1805   mySurface.Initialize( face );
1806   myParams.resize( 1 );
1807   myParams[ 0 ] = make_pair( first, last );
1808   myIsForward = StdMeshers_PrismAsBlock::IsForwardEdge( myHelper->GetMeshDS(),
1809                                                         *myParamToColumnMap,
1810                                                         myBaseEdge, myID );
1811 }
1812
1813 //================================================================================
1814 /*!
1815  * \brief Constructor of complex side face
1816  */
1817 //================================================================================
1818
1819 StdMeshers_PrismAsBlock::TSideFace::
1820 TSideFace(const vector< TSideFace* >&             components,
1821           const vector< pair< double, double> > & params)
1822   :myID( components[0] ? components[0]->myID : 0 ),
1823    myParamToColumnMap( 0 ),
1824    myParams( params ),
1825    myIsForward( true ),
1826    myComponents( components ),
1827    myHelper( components[0] ? components[0]->myHelper : 0 )
1828 {}
1829 //================================================================================
1830 /*!
1831  * \brief Copy constructor
1832   * \param other - other side
1833  */
1834 //================================================================================
1835
1836 StdMeshers_PrismAsBlock::TSideFace::TSideFace( const TSideFace& other )
1837 {
1838   myID               = other.myID;
1839   mySurface          = other.mySurface;
1840   myBaseEdge         = other.myBaseEdge;
1841   myParams           = other.myParams;
1842   myIsForward        = other.myIsForward;
1843   myHelper           = other.myHelper;
1844   myParamToColumnMap = other.myParamToColumnMap;
1845
1846   myComponents.resize( other.myComponents.size());
1847   for (int i = 0 ; i < myComponents.size(); ++i )
1848     myComponents[ i ] = new TSideFace( *other.myComponents[ i ]);
1849 }
1850
1851 //================================================================================
1852 /*!
1853  * \brief Deletes myComponents
1854  */
1855 //================================================================================
1856
1857 StdMeshers_PrismAsBlock::TSideFace::~TSideFace()
1858 {
1859   for (int i = 0 ; i < myComponents.size(); ++i )
1860     if ( myComponents[ i ] )
1861       delete myComponents[ i ];
1862 }
1863
1864 //================================================================================
1865 /*!
1866  * \brief Return geometry of the vertical curve
1867   * \param isMax - true means curve located closer to (1,1,1) block point
1868   * \retval Adaptor3d_Curve* - curve adaptor
1869  */
1870 //================================================================================
1871
1872 Adaptor3d_Curve* StdMeshers_PrismAsBlock::TSideFace::VertiCurve(const bool isMax) const
1873 {
1874   if ( !myComponents.empty() ) {
1875     if ( isMax )
1876       return myComponents.back()->VertiCurve(isMax);
1877     else
1878       return myComponents.front()->VertiCurve(isMax);
1879   }
1880   double f = myParams[0].first, l = myParams[0].second;
1881   if ( !myIsForward ) std::swap( f, l );
1882   return new TVerticalEdgeAdaptor( myParamToColumnMap, isMax ? l : f );
1883 }
1884
1885 //================================================================================
1886 /*!
1887  * \brief Return geometry of the top or bottom curve
1888   * \param isTop - 
1889   * \retval Adaptor3d_Curve* - 
1890  */
1891 //================================================================================
1892
1893 Adaptor3d_Curve* StdMeshers_PrismAsBlock::TSideFace::HorizCurve(const bool isTop) const
1894 {
1895   return new THorizontalEdgeAdaptor( this, isTop );
1896 }
1897
1898 //================================================================================
1899 /*!
1900  * \brief Return pcurves
1901   * \param pcurv - array of 4 pcurves
1902   * \retval bool - is a success
1903  */
1904 //================================================================================
1905
1906 bool StdMeshers_PrismAsBlock::TSideFace::GetPCurves(Adaptor2d_Curve2d* pcurv[4]) const
1907 {
1908   int iEdge[ 4 ] = { BOTTOM_EDGE, TOP_EDGE, V0_EDGE, V1_EDGE };
1909
1910   for ( int i = 0 ; i < 4 ; ++i ) {
1911     Handle(Geom2d_Line) line;
1912     switch ( iEdge[ i ] ) {
1913     case TOP_EDGE:
1914       line = new Geom2d_Line( gp_Pnt2d( 0, 1 ), gp::DX2d() ); break;
1915     case BOTTOM_EDGE:
1916       line = new Geom2d_Line( gp::Origin2d(), gp::DX2d() ); break;
1917     case V0_EDGE:
1918       line = new Geom2d_Line( gp::Origin2d(), gp::DY2d() ); break;
1919     case V1_EDGE:
1920       line = new Geom2d_Line( gp_Pnt2d( 1, 0 ), gp::DY2d() ); break;
1921     }
1922     pcurv[ i ] = new Geom2dAdaptor_Curve( line, 0, 1 );
1923   }
1924   return true;
1925 }
1926
1927 //================================================================================
1928 /*!
1929  * \brief Returns geometry of pcurve on a horizontal face
1930   * \param isTop - is top or bottom face
1931   * \param horFace - a horizontal face
1932   * \retval Adaptor2d_Curve2d* - curve adaptor
1933  */
1934 //================================================================================
1935
1936 Adaptor2d_Curve2d*
1937 StdMeshers_PrismAsBlock::TSideFace::HorizPCurve(const bool         isTop,
1938                                                 const TopoDS_Face& horFace) const
1939 {
1940   return new TPCurveOnHorFaceAdaptor( this, isTop, horFace );
1941 }
1942
1943 //================================================================================
1944 /*!
1945  * \brief Return a component corresponding to parameter
1946   * \param U - parameter along a horizontal size
1947   * \param localU - parameter along a horizontal size of a component
1948   * \retval TSideFace* - found component
1949  */
1950 //================================================================================
1951
1952 StdMeshers_PrismAsBlock::TSideFace*
1953 StdMeshers_PrismAsBlock::TSideFace::GetComponent(const double U,double & localU) const
1954 {
1955   localU = U;
1956   if ( myComponents.empty() )
1957     return const_cast<TSideFace*>( this );
1958
1959   int i;
1960   for ( i = 0; i < myComponents.size(); ++i )
1961     if ( U < myParams[ i ].second )
1962       break;
1963   if ( i >= myComponents.size() )
1964     i = myComponents.size() - 1;
1965
1966   double f = myParams[ i ].first, l = myParams[ i ].second;
1967   localU = ( U - f ) / ( l - f );
1968   return myComponents[ i ];
1969 }
1970
1971 //================================================================================
1972 /*!
1973  * \brief Find node columns for a parameter
1974   * \param U - parameter along a horizontal edge
1975   * \param col1 - the 1st found column
1976   * \param col2 - the 2nd found column
1977   * \retval r - normalized position of U between the found columns
1978  */
1979 //================================================================================
1980
1981 double StdMeshers_PrismAsBlock::TSideFace::GetColumns(const double      U,
1982                                                       TParam2ColumnIt & col1,
1983                                                       TParam2ColumnIt & col2) const
1984 {
1985   double u = U, r = 0;
1986   if ( !myComponents.empty() ) {
1987     TSideFace * comp = GetComponent(U,u);
1988     return comp->GetColumns( u, col1, col2 );
1989   }
1990
1991   if ( !myIsForward )
1992     u = 1 - u;
1993   double f = myParams[0].first, l = myParams[0].second;
1994   u = f + u * ( l - f );
1995
1996   col1 = col2 = getColumn( myParamToColumnMap, u );
1997   if ( ++col2 == myParamToColumnMap->end() ) {
1998     --col2;
1999     r = 0.5;
2000   }
2001   else {
2002     double uf = col1->first;
2003     double ul = col2->first;
2004     r = ( u - uf ) / ( ul - uf );
2005   }
2006   return r;
2007 }
2008
2009 //================================================================================
2010 /*!
2011  * \brief Return coordinates by normalized params
2012   * \param U - horizontal param
2013   * \param V - vertical param
2014   * \retval gp_Pnt - result point
2015  */
2016 //================================================================================
2017
2018 gp_Pnt StdMeshers_PrismAsBlock::TSideFace::Value(const Standard_Real U,
2019                                                  const Standard_Real V) const
2020 {
2021   if ( !myComponents.empty() ) {
2022     double u;
2023     TSideFace * comp = GetComponent(U,u);
2024     return comp->Value( u, V );
2025   }
2026
2027   TParam2ColumnIt u_col1, u_col2;
2028   double vR, hR = GetColumns( U, u_col1, u_col2 );
2029
2030   const SMDS_MeshNode* n1 = 0;
2031   const SMDS_MeshNode* n2 = 0;
2032   const SMDS_MeshNode* n3 = 0;
2033   const SMDS_MeshNode* n4 = 0;
2034
2035   // BEGIN issue 0020680: EDF 1252 SMESH: Bad cell created by Radial prism in center of torus
2036   // Workaround for a wrongly located point returned by mySurface.Value() for
2037   // UV located near boundary of BSpline surface.
2038   // To bypass the problem, we take point from 3D curve of edge.
2039   // It solves pb of the bloc_fiss_new.py
2040   const double tol = 1e-3;
2041   if ( V < tol || V+tol >= 1. )
2042   {
2043     n1 = V < tol ? u_col1->second.front() : u_col1->second.back();
2044     n3 = V < tol ? u_col2->second.front() : u_col2->second.back();
2045     TopoDS_Edge edge;
2046     if ( V < tol )
2047     {
2048       edge = myBaseEdge;
2049     }
2050     else
2051     {
2052       TopoDS_Shape s = myHelper->GetSubShapeByNode( n1, myHelper->GetMeshDS() );
2053       if ( s.ShapeType() != TopAbs_EDGE )
2054         s = myHelper->GetSubShapeByNode( n3, myHelper->GetMeshDS() );
2055       if ( s.ShapeType() == TopAbs_EDGE )
2056         edge = TopoDS::Edge( s );
2057     }
2058     if ( !edge.IsNull() )
2059     {
2060       double u1 = myHelper->GetNodeU( edge, n1 );
2061       double u3 = myHelper->GetNodeU( edge, n3 );
2062       double u = u1 * ( 1 - hR ) + u3 * hR;
2063       TopLoc_Location loc; double f,l;
2064       Handle(Geom_Curve) curve = BRep_Tool::Curve( edge,loc,f,l );
2065       return curve->Value( u ).Transformed( loc );
2066     }
2067   }
2068   // END issue 0020680: EDF 1252 SMESH: Bad cell created by Radial prism in center of torus
2069
2070   vR = getRAndNodes( & u_col1->second, V, n1, n2 );
2071   vR = getRAndNodes( & u_col2->second, V, n3, n4 );
2072   
2073   gp_XY uv1 = myHelper->GetNodeUV( mySurface.Face(), n1, n4);
2074   gp_XY uv2 = myHelper->GetNodeUV( mySurface.Face(), n2, n3);
2075   gp_XY uv12 = uv1 * ( 1 - vR ) + uv2 * vR;
2076
2077   gp_XY uv3 = myHelper->GetNodeUV( mySurface.Face(), n3, n2);
2078   gp_XY uv4 = myHelper->GetNodeUV( mySurface.Face(), n4, n1);
2079   gp_XY uv34 = uv3 * ( 1 - vR ) + uv4 * vR;
2080
2081   gp_XY uv = uv12 * ( 1 - hR ) + uv34 * hR;
2082
2083   gp_Pnt p = mySurface.Value( uv.X(), uv.Y() );
2084   return p;
2085 }
2086
2087
2088 //================================================================================
2089 /*!
2090  * \brief Return boundary edge
2091   * \param edge - edge index
2092   * \retval TopoDS_Edge - found edge
2093  */
2094 //================================================================================
2095
2096 TopoDS_Edge StdMeshers_PrismAsBlock::TSideFace::GetEdge(const int iEdge) const
2097 {
2098   if ( !myComponents.empty() ) {
2099     switch ( iEdge ) {
2100     case V0_EDGE : return myComponents.front()->GetEdge( iEdge );
2101     case V1_EDGE : return myComponents.back() ->GetEdge( iEdge );
2102     default: return TopoDS_Edge();
2103     }
2104   }
2105   TopoDS_Shape edge;
2106   const SMDS_MeshNode* node = 0;
2107   SMESHDS_Mesh * meshDS = myHelper->GetMesh()->GetMeshDS();
2108   TNodeColumn* column;
2109
2110   switch ( iEdge ) {
2111   case TOP_EDGE:
2112   case BOTTOM_EDGE:
2113     column = & (( ++myParamToColumnMap->begin())->second );
2114     node = ( iEdge == TOP_EDGE ) ? column->back() : column->front();
2115     edge = myHelper->GetSubShapeByNode ( node, meshDS );
2116     if ( edge.ShapeType() == TopAbs_VERTEX ) {
2117       column = & ( myParamToColumnMap->begin()->second );
2118       node = ( iEdge == TOP_EDGE ) ? column->back() : column->front();
2119     }
2120     break;
2121   case V0_EDGE:
2122   case V1_EDGE: {
2123     bool back = ( iEdge == V1_EDGE );
2124     if ( !myIsForward ) back = !back;
2125     if ( back )
2126       column = & ( myParamToColumnMap->rbegin()->second );
2127     else
2128       column = & ( myParamToColumnMap->begin()->second );
2129     if ( column->size() > 0 )
2130       edge = myHelper->GetSubShapeByNode( (*column)[ 1 ], meshDS );
2131     if ( edge.IsNull() || edge.ShapeType() == TopAbs_VERTEX )
2132       node = column->front();
2133     break;
2134   }
2135   default:;
2136   }
2137   if ( !edge.IsNull() && edge.ShapeType() == TopAbs_EDGE )
2138     return TopoDS::Edge( edge );
2139
2140   // find edge by 2 vertices
2141   TopoDS_Shape V1 = edge;
2142   TopoDS_Shape V2 = myHelper->GetSubShapeByNode( node, meshDS );
2143   if ( V2.ShapeType() == TopAbs_VERTEX && !V2.IsSame( V1 ))
2144   {
2145     TopTools_ListIteratorOfListOfShape ancestIt =
2146       myHelper->GetMesh()->GetAncestors( V1 );
2147     for ( ; ancestIt.More(); ancestIt.Next() )
2148     {
2149       const TopoDS_Shape & ancestor = ancestIt.Value();
2150       if ( ancestor.ShapeType() == TopAbs_EDGE )
2151         for ( TopExp_Explorer e( ancestor, TopAbs_VERTEX ); e.More(); e.Next() )
2152           if ( V2.IsSame( e.Current() ))
2153             return TopoDS::Edge( ancestor );
2154     }
2155   }
2156   return TopoDS_Edge();
2157 }
2158
2159 //================================================================================
2160 /*!
2161  * \brief Fill block subshapes
2162   * \param shapeMap - map to fill in
2163   * \retval int - nb inserted subshapes
2164  */
2165 //================================================================================
2166
2167 int StdMeshers_PrismAsBlock::TSideFace::InsertSubShapes(TBlockShapes& shapeMap) const
2168 {
2169   int nbInserted = 0;
2170
2171   // Insert edges
2172   vector< int > edgeIdVec;
2173   SMESH_Block::GetFaceEdgesIDs( myID, edgeIdVec );
2174
2175   for ( int i = BOTTOM_EDGE; i <=V1_EDGE ; ++i ) {
2176     TopoDS_Edge e = GetEdge( i );
2177     if ( !e.IsNull() ) {
2178       nbInserted += SMESH_Block::Insert( e, edgeIdVec[ i ], shapeMap);
2179     }
2180   }
2181
2182   // Insert corner vertices
2183
2184   TParam2ColumnIt col1, col2 ;
2185   vector< int > vertIdVec;
2186
2187   // from V0 column
2188   SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ V0_EDGE ], vertIdVec);
2189   GetColumns(0, col1, col2 );
2190   const SMDS_MeshNode* node0 = col1->second.front();
2191   const SMDS_MeshNode* node1 = col1->second.back();
2192   TopoDS_Shape v0 = myHelper->GetSubShapeByNode( node0, myHelper->GetMeshDS());
2193   TopoDS_Shape v1 = myHelper->GetSubShapeByNode( node1, myHelper->GetMeshDS());
2194   if ( v0.ShapeType() == TopAbs_VERTEX ) {
2195     nbInserted += SMESH_Block::Insert( v0, vertIdVec[ 0 ], shapeMap);
2196   }
2197   if ( v1.ShapeType() == TopAbs_VERTEX ) {
2198     nbInserted += SMESH_Block::Insert( v1, vertIdVec[ 1 ], shapeMap);
2199   }
2200   
2201   // from V1 column
2202   SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ V1_EDGE ], vertIdVec);
2203   GetColumns(1, col1, col2 );
2204   node0 = col2->second.front();
2205   node1 = col2->second.back();
2206   v0 = myHelper->GetSubShapeByNode( node0, myHelper->GetMeshDS());
2207   v1 = myHelper->GetSubShapeByNode( node1, myHelper->GetMeshDS());
2208   if ( v0.ShapeType() == TopAbs_VERTEX ) {
2209     nbInserted += SMESH_Block::Insert( v0, vertIdVec[ 0 ], shapeMap);
2210   }
2211   if ( v1.ShapeType() == TopAbs_VERTEX ) {
2212     nbInserted += SMESH_Block::Insert( v1, vertIdVec[ 1 ], shapeMap);
2213   }
2214
2215 //   TopoDS_Vertex V0, V1, Vcom;
2216 //   TopExp::Vertices( myBaseEdge, V0, V1, true );
2217 //   if ( !myIsForward ) std::swap( V0, V1 );
2218
2219 //   // bottom vertex IDs
2220 //   SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ _u0 ], vertIdVec);
2221 //   SMESH_Block::Insert( V0, vertIdVec[ 0 ], shapeMap);
2222 //   SMESH_Block::Insert( V1, vertIdVec[ 1 ], shapeMap);
2223
2224 //   TopoDS_Edge sideEdge = GetEdge( V0_EDGE );
2225 //   if ( sideEdge.IsNull() || !TopExp::CommonVertex( botEdge, sideEdge, Vcom ))
2226 //     return false;
2227
2228 //   // insert one side edge
2229 //   int edgeID;
2230 //   if ( Vcom.IsSame( V0 )) edgeID = edgeIdVec[ _v0 ];
2231 //   else                    edgeID = edgeIdVec[ _v1 ];
2232 //   SMESH_Block::Insert( sideEdge, edgeID, shapeMap);
2233
2234 //   // top vertex of the side edge
2235 //   SMESH_Block::GetEdgeVertexIDs( edgeID, vertIdVec);
2236 //   TopoDS_Vertex Vtop = TopExp::FirstVertex( sideEdge );
2237 //   if ( Vcom.IsSame( Vtop ))
2238 //     Vtop = TopExp::LastVertex( sideEdge );
2239 //   SMESH_Block::Insert( Vtop, vertIdVec[ 1 ], shapeMap);
2240
2241 //   // other side edge
2242 //   sideEdge = GetEdge( V1_EDGE );
2243 //   if ( sideEdge.IsNull() )
2244 //     return false;
2245 //   if ( edgeID = edgeIdVec[ _v1 ]) edgeID = edgeIdVec[ _v0 ];
2246 //   else                            edgeID = edgeIdVec[ _v1 ];
2247 //   SMESH_Block::Insert( sideEdge, edgeID, shapeMap);
2248   
2249 //   // top edge
2250 //   TopoDS_Edge topEdge = GetEdge( TOP_EDGE );
2251 //   SMESH_Block::Insert( topEdge, edgeIdVec[ _u1 ], shapeMap);
2252
2253 //   // top vertex of the other side edge
2254 //   if ( !TopExp::CommonVertex( topEdge, sideEdge, Vcom ))
2255 //     return false;
2256 //   SMESH_Block::GetEdgeVertexIDs( edgeID, vertIdVec );
2257 //   SMESH_Block::Insert( Vcom, vertIdVec[ 1 ], shapeMap);
2258
2259   return nbInserted;
2260 }
2261
2262 //================================================================================
2263 /*!
2264  * \brief Dump ids of nodes of sides
2265  */
2266 //================================================================================
2267
2268 void StdMeshers_PrismAsBlock::TSideFace::dumpNodes(int nbNodes) const
2269 {
2270 #ifdef _DEBUG_
2271   cout << endl << "NODES OF FACE "; SMESH_Block::DumpShapeID( myID, cout ) << endl;
2272   THorizontalEdgeAdaptor* hSize0 = (THorizontalEdgeAdaptor*) HorizCurve(0);
2273   cout << "Horiz side 0: "; hSize0->dumpNodes(nbNodes); cout << endl;
2274   THorizontalEdgeAdaptor* hSize1 = (THorizontalEdgeAdaptor*) HorizCurve(1);
2275   cout << "Horiz side 1: "; hSize1->dumpNodes(nbNodes); cout << endl;
2276   TVerticalEdgeAdaptor* vSide0 = (TVerticalEdgeAdaptor*) VertiCurve(0);
2277   cout << "Verti side 0: "; vSide0->dumpNodes(nbNodes); cout << endl;
2278   TVerticalEdgeAdaptor* vSide1 = (TVerticalEdgeAdaptor*) VertiCurve(1);
2279   cout << "Verti side 1: "; vSide1->dumpNodes(nbNodes); cout << endl;
2280   delete hSize0; delete hSize1; delete vSide0; delete vSide1;
2281 #endif
2282 }
2283
2284 //================================================================================
2285 /*!
2286  * \brief Creates TVerticalEdgeAdaptor 
2287   * \param columnsMap - node column map
2288   * \param parameter - normalized parameter
2289  */
2290 //================================================================================
2291
2292 StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::
2293 TVerticalEdgeAdaptor( const TParam2ColumnMap* columnsMap, const double parameter)
2294 {
2295   myNodeColumn = & getColumn( columnsMap, parameter )->second;
2296 }
2297
2298 //================================================================================
2299 /*!
2300  * \brief Return coordinates for the given normalized parameter
2301   * \param U - normalized parameter
2302   * \retval gp_Pnt - coordinates
2303  */
2304 //================================================================================
2305
2306 gp_Pnt StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::Value(const Standard_Real U) const
2307 {
2308   const SMDS_MeshNode* n1;
2309   const SMDS_MeshNode* n2;
2310   double r = getRAndNodes( myNodeColumn, U, n1, n2 );
2311   return gpXYZ(n1) * ( 1 - r ) + gpXYZ(n2) * r;
2312 }
2313
2314 //================================================================================
2315 /*!
2316  * \brief Dump ids of nodes
2317  */
2318 //================================================================================
2319
2320 void StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::dumpNodes(int nbNodes) const
2321 {
2322 #ifdef _DEBUG_
2323   for ( int i = 0; i < nbNodes && i < myNodeColumn->size(); ++i )
2324     cout << (*myNodeColumn)[i]->GetID() << " ";
2325   if ( nbNodes < myNodeColumn->size() )
2326     cout << myNodeColumn->back()->GetID();
2327 #endif
2328 }
2329
2330 //================================================================================
2331 /*!
2332  * \brief Return coordinates for the given normalized parameter
2333   * \param U - normalized parameter
2334   * \retval gp_Pnt - coordinates
2335  */
2336 //================================================================================
2337
2338 gp_Pnt StdMeshers_PrismAsBlock::THorizontalEdgeAdaptor::Value(const Standard_Real U) const
2339 {
2340   return mySide->TSideFace::Value( U, myV );
2341 }
2342
2343 //================================================================================
2344 /*!
2345  * \brief Dump ids of <nbNodes> first nodes and the last one
2346  */
2347 //================================================================================
2348
2349 void StdMeshers_PrismAsBlock::THorizontalEdgeAdaptor::dumpNodes(int nbNodes) const
2350 {
2351 #ifdef _DEBUG_
2352   // Not bedugged code. Last node is sometimes incorrect
2353   const TSideFace* side = mySide;
2354   double u = 0;
2355   if ( mySide->IsComplex() )
2356     side = mySide->GetComponent(0,u);
2357
2358   TParam2ColumnIt col, col2;
2359   TParam2ColumnMap* u2cols = side->GetColumns();
2360   side->GetColumns( u , col, col2 );
2361   
2362   int j, i = myV ? mySide->ColumnHeight()-1 : 0;
2363
2364   const SMDS_MeshNode* n = 0;
2365   const SMDS_MeshNode* lastN
2366     = side->IsForward() ? u2cols->rbegin()->second[ i ] : u2cols->begin()->second[ i ];
2367   for ( j = 0; j < nbNodes && n != lastN; ++j )
2368   {
2369     n = col->second[ i ];
2370     cout << n->GetID() << " ";
2371     if ( side->IsForward() )
2372       ++col;
2373     else
2374       --col;
2375   }
2376
2377   // last node
2378   u = 1;
2379   if ( mySide->IsComplex() )
2380     side = mySide->GetComponent(1,u);
2381
2382   side->GetColumns( u , col, col2 );
2383   if ( n != col->second[ i ] )
2384     cout << col->second[ i ]->GetID();
2385 #endif
2386 }
2387 //================================================================================
2388 /*!
2389  * \brief Return UV on pcurve for the given normalized parameter
2390   * \param U - normalized parameter
2391   * \retval gp_Pnt - coordinates
2392  */
2393 //================================================================================
2394
2395 gp_Pnt2d StdMeshers_PrismAsBlock::TPCurveOnHorFaceAdaptor::Value(const Standard_Real U) const
2396 {
2397   TParam2ColumnIt u_col1, u_col2;
2398   double r = mySide->GetColumns( U, u_col1, u_col2 );
2399   gp_XY uv1 = mySide->GetNodeUV( myFace, u_col1->second[ myZ ]);
2400   gp_XY uv2 = mySide->GetNodeUV( myFace, u_col2->second[ myZ ]);
2401   return uv1 * ( 1 - r ) + uv2 * r;
2402 }