1 // Copyright (C) 2007-2012 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
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.
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.
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
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 // File : StdMeshers_Prism_3D.cxx
25 // Created : Fri Oct 20 11:37:07 2006
26 // Author : Edward AGAPOV (eap)
28 #include "StdMeshers_Prism_3D.hxx"
30 #include "StdMeshers_ProjectionUtils.hxx"
31 #include "SMESH_MesherHelper.hxx"
32 #include "SMDS_VolumeTool.hxx"
33 #include "SMDS_VolumeOfNodes.hxx"
34 #include "SMDS_EdgePosition.hxx"
35 #include "SMESH_Comment.hxx"
37 #include "utilities.h"
39 #include <BRep_Tool.hxx>
40 #include <Bnd_B3d.hxx>
41 #include <Geom2dAdaptor_Curve.hxx>
42 #include <Geom2d_Line.hxx>
43 #include <Geom_Curve.hxx>
45 #include <TopExp_Explorer.hxx>
46 #include <TopTools_ListIteratorOfListOfShape.hxx>
47 #include <TopTools_MapOfShape.hxx>
48 #include <TopTools_SequenceOfShape.hxx>
55 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
56 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
57 #define SHOWYXZ(msg, xyz) // {\
59 // cout << msg << " ("<< p.X() << "; " <<p.Y() << "; " <<p.Z() << ") " <<endl;\
62 typedef StdMeshers_ProjectionUtils TAssocTool;
63 typedef SMESH_Comment TCom;
65 enum { ID_BOT_FACE = SMESH_Block::ID_Fxy0,
66 ID_TOP_FACE = SMESH_Block::ID_Fxy1,
67 BOTTOM_EDGE = 0, TOP_EDGE, V0_EDGE, V1_EDGE, // edge IDs in face
68 NB_WALL_FACES = 4 }; //
72 //================================================================================
74 * \brief Return iterator pointing to node column for the given parameter
75 * \param columnsMap - node column map
76 * \param parameter - parameter
77 * \retval TParam2ColumnMap::iterator - result
79 * it returns closest left column
81 //================================================================================
83 TParam2ColumnIt getColumn( const TParam2ColumnMap* columnsMap,
84 const double parameter )
86 TParam2ColumnIt u_col = columnsMap->upper_bound( parameter );
87 if ( u_col != columnsMap->begin() )
89 return u_col; // return left column
92 //================================================================================
94 * \brief Return nodes around given parameter and a ratio
95 * \param column - node column
96 * \param param - parameter
97 * \param node1 - lower node
98 * \param node2 - upper node
99 * \retval double - ratio
101 //================================================================================
103 double getRAndNodes( const TNodeColumn* column,
105 const SMDS_MeshNode* & node1,
106 const SMDS_MeshNode* & node2)
108 if ( param >= 1.0 || column->size() == 1) {
109 node1 = node2 = column->back();
113 int i = int( param * ( column->size() - 1 ));
114 double u0 = double( i )/ double( column->size() - 1 );
115 double r = ( param - u0 ) * ( column->size() - 1 );
117 node1 = (*column)[ i ];
118 node2 = (*column)[ i + 1];
122 //================================================================================
124 * \brief Compute boundary parameters of face parts
125 * \param nbParts - nb of parts to split columns into
126 * \param columnsMap - node columns of the face to split
127 * \param params - computed parameters
129 //================================================================================
131 void splitParams( const int nbParts,
132 const TParam2ColumnMap* columnsMap,
133 vector< double > & params)
136 params.reserve( nbParts + 1 );
137 TParam2ColumnIt last_par_col = --columnsMap->end();
138 double par = columnsMap->begin()->first; // 0.
139 double parLast = last_par_col->first;
140 params.push_back( par );
141 for ( int i = 0; i < nbParts - 1; ++ i )
143 double partSize = ( parLast - par ) / double ( nbParts - i );
144 TParam2ColumnIt par_col = getColumn( columnsMap, par + partSize );
145 if ( par_col->first == par ) {
147 if ( par_col == last_par_col ) {
148 while ( i < nbParts - 1 )
149 params.push_back( par + partSize * i++ );
153 par = par_col->first;
154 params.push_back( par );
156 params.push_back( parLast ); // 1.
159 //================================================================================
161 * \brief Return coordinate system for z-th layer of nodes
163 //================================================================================
165 gp_Ax2 getLayerCoordSys(const int z,
166 const vector< const TNodeColumn* >& columns,
169 // gravity center of a layer
172 for ( int i = 0; i < columns.size(); ++i )
174 O += gpXYZ( (*columns[ i ])[ z ]);
175 if ( vertexCol < 0 &&
176 columns[ i ]->front()->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX )
183 int iPrev = columns.size()-1;
184 for ( int i = 0; i < columns.size(); ++i )
186 gp_Vec v1( O, gpXYZ( (*columns[ iPrev ])[ z ]));
187 gp_Vec v2( O, gpXYZ( (*columns[ i ] )[ z ]));
192 if ( vertexCol >= 0 )
194 O = gpXYZ( (*columns[ vertexCol ])[ z ]);
196 if ( xColumn < 0 || xColumn >= columns.size() )
198 // select a column for X dir
200 for ( int i = 0; i < columns.size(); ++i )
202 double dist = ( O - gpXYZ((*columns[ i ])[ z ])).SquareModulus();
203 if ( dist > maxDist )
212 gp_Vec X( O, gpXYZ( (*columns[ xColumn ])[ z ]));
214 return gp_Ax2( O, Z, X);
217 //================================================================================
219 * \brief Removes submeshes meshed with regular grid from given list
220 * \retval int - nb of removed submeshes
222 //================================================================================
224 int removeQuasiQuads(list< SMESH_subMesh* >& notQuadSubMesh)
226 int oldNbSM = notQuadSubMesh.size();
227 SMESHDS_Mesh* mesh = notQuadSubMesh.front()->GetFather()->GetMeshDS();
228 list< SMESH_subMesh* >::iterator smIt = notQuadSubMesh.begin();
229 #define __NEXT_SM { ++smIt; continue; }
230 while ( smIt != notQuadSubMesh.end() )
232 SMESH_subMesh* faceSm = *smIt;
233 SMESHDS_SubMesh* faceSmDS = faceSm->GetSubMeshDS();
234 int nbQuads = faceSmDS->NbElements();
235 if ( nbQuads == 0 ) __NEXT_SM;
237 // get oredered edges
238 list< TopoDS_Edge > orderedEdges;
239 list< int > nbEdgesInWires;
241 int nbWires = SMESH_Block::GetOrderedEdges( TopoDS::Face( faceSm->GetSubShape() ),
242 V000, orderedEdges, nbEdgesInWires );
243 if ( nbWires != 1 || nbEdgesInWires.front() <= 4 )
246 // get nb of segements on edges
247 list<int> nbSegOnEdge;
248 list< TopoDS_Edge >::iterator edge = orderedEdges.begin();
249 for ( ; edge != orderedEdges.end(); ++edge )
251 if ( SMESHDS_SubMesh* edgeSmDS = mesh->MeshElements( *edge ))
252 nbSegOnEdge.push_back( edgeSmDS->NbElements() );
254 nbSegOnEdge.push_back(0);
257 // unite nbSegOnEdge of continues edges
258 int nbEdges = nbEdgesInWires.front();
259 list<int>::iterator nbSegIt = nbSegOnEdge.begin();
260 for ( edge = orderedEdges.begin(); edge != orderedEdges.end(); )
262 const TopoDS_Edge& e1 = *edge++;
263 const TopoDS_Edge& e2 = ( edge == orderedEdges.end() ? orderedEdges.front() : *edge );
264 if ( SMESH_Algo::IsContinuous( e1, e2 ))
266 // common vertex of continues edges must be shared by two 2D mesh elems of geom face
267 TopoDS_Vertex vCommon = TopExp::LastVertex( e1, true );
268 const SMDS_MeshNode* vNode = SMESH_Algo::VertexNode( vCommon, mesh );
272 SMDS_ElemIteratorPtr fIt = vNode->GetInverseElementIterator(SMDSAbs_Face);
273 while ( fIt->more() )
274 nbF += faceSmDS->Contains( fIt->next() );
276 list<int>::iterator nbSegIt1 = nbSegIt++;
277 if ( !vNode || nbF == 2 ) // !vNode - two edges can be meshed as one
280 if ( nbSegIt == nbSegOnEdge.end() ) nbSegIt = nbSegOnEdge.begin();
281 *nbSegIt += *nbSegIt1;
282 nbSegOnEdge.erase( nbSegIt1 );
291 vector<int> nbSegVec( nbSegOnEdge.begin(), nbSegOnEdge.end());
292 if ( nbSegVec.size() == 4 &&
293 nbSegVec[0] == nbSegVec[2] &&
294 nbSegVec[1] == nbSegVec[3] &&
295 nbSegVec[0] * nbSegVec[1] == nbQuads
297 smIt = notQuadSubMesh.erase( smIt );
302 return oldNbSM - notQuadSubMesh.size();
306 //=======================================================================
307 //function : StdMeshers_Prism_3D
309 //=======================================================================
311 StdMeshers_Prism_3D::StdMeshers_Prism_3D(int hypId, int studyId, SMESH_Gen* gen)
312 :SMESH_3D_Algo(hypId, studyId, gen)
315 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID); // 1 bit per shape type
316 myProjectTriangles = false;
319 //================================================================================
323 //================================================================================
325 StdMeshers_Prism_3D::~StdMeshers_Prism_3D()
328 //=======================================================================
329 //function : CheckHypothesis
331 //=======================================================================
333 bool StdMeshers_Prism_3D::CheckHypothesis(SMESH_Mesh& aMesh,
334 const TopoDS_Shape& aShape,
335 SMESH_Hypothesis::Hypothesis_Status& aStatus)
337 // Check shape geometry
339 aStatus = SMESH_Hypothesis::HYP_BAD_GEOMETRY;
341 // find not quadrangle faces
342 list< TopoDS_Shape > notQuadFaces;
343 int nbEdge, nbWire, nbFace = 0;
344 TopExp_Explorer exp( aShape, TopAbs_FACE );
345 for ( ; exp.More(); exp.Next() ) {
347 const TopoDS_Shape& face = exp.Current();
348 nbEdge = TAssocTool::Count( face, TopAbs_EDGE, 0 );
349 nbWire = TAssocTool::Count( face, TopAbs_WIRE, 0 );
350 if ( nbEdge!= 4 || nbWire!= 1 ) {
351 if ( !notQuadFaces.empty() ) {
352 if ( TAssocTool::Count( notQuadFaces.back(), TopAbs_EDGE, 0 ) != nbEdge ||
353 TAssocTool::Count( notQuadFaces.back(), TopAbs_WIRE, 0 ) != nbWire )
354 RETURN_BAD_RESULT("Different not quad faces");
356 notQuadFaces.push_back( face );
359 if ( !notQuadFaces.empty() )
361 if ( notQuadFaces.size() != 2 )
362 RETURN_BAD_RESULT("Bad nb not quad faces: " << notQuadFaces.size());
364 // check total nb faces
365 nbEdge = TAssocTool::Count( notQuadFaces.back(), TopAbs_EDGE, 0 );
366 if ( nbFace != nbEdge + 2 )
367 RETURN_BAD_RESULT("Bad nb of faces: " << nbFace << " but must be " << nbEdge + 2);
371 aStatus = SMESH_Hypothesis::HYP_OK;
375 //=======================================================================
378 //=======================================================================
380 bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theShape)
382 SMESH_MesherHelper helper( theMesh );
385 myHelper->IsQuadraticSubMesh( theShape );
387 // Analyse mesh and geomerty to find block sub-shapes and submeshes
388 if ( !myBlock.Init( myHelper, theShape ))
389 return error( myBlock.GetError());
391 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
393 int volumeID = meshDS->ShapeToIndex( theShape );
396 // To compute coordinates of a node inside a block, it is necessary to know
397 // 1. normalized parameters of the node by which
398 // 2. coordinates of node projections on all block sub-shapes are computed
400 // So we fill projections on vertices at once as they are same for all nodes
401 myShapeXYZ.resize( myBlock.NbSubShapes() );
402 for ( int iV = SMESH_Block::ID_FirstV; iV < SMESH_Block::ID_FirstE; ++iV ) {
403 myBlock.VertexPoint( iV, myShapeXYZ[ iV ]);
404 SHOWYXZ("V point " <<iV << " ", myShapeXYZ[ iV ]);
407 // Projections on the top and bottom faces are taken from nodes existing
408 // on these faces; find correspondence between bottom and top nodes
409 myBotToColumnMap.clear();
410 if ( !assocOrProjBottom2Top() ) // it also fills myBotToColumnMap
414 // Create nodes inside the block
416 // try to use transformation (issue 0020680)
417 vector<gp_Trsf> trsf;
418 if ( myBlock.GetLayersTransformation(trsf))
420 // loop on nodes inside the bottom face
421 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.begin();
422 for ( ; bot_column != myBotToColumnMap.end(); ++bot_column )
424 const TNode& tBotNode = bot_column->first; // bottom TNode
425 if ( tBotNode.GetPositionType() != SMDS_TOP_FACE )
426 continue; // node is not inside face
428 // column nodes; middle part of the column are zero pointers
429 TNodeColumn& column = bot_column->second;
430 TNodeColumn::iterator columnNodes = column.begin();
431 for ( int z = 0; columnNodes != column.end(); ++columnNodes, ++z)
433 const SMDS_MeshNode* & node = *columnNodes;
434 if ( node ) continue; // skip bottom or top node
436 gp_XYZ coords = tBotNode.GetCoords();
437 trsf[z-1].Transforms( coords );
438 node = meshDS->AddNode( coords.X(), coords.Y(), coords.Z() );
439 meshDS->SetNodeInVolume( node, volumeID );
441 } // loop on bottom nodes
443 else // use block approach
445 // loop on nodes inside the bottom face
447 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.begin();
448 for ( ; bot_column != myBotToColumnMap.end(); ++bot_column )
450 const TNode& tBotNode = bot_column->first; // bottom TNode
451 if ( tBotNode.GetPositionType() != SMDS_TOP_FACE )
452 continue; // node is not inside face
454 // column nodes; middle part of the column are zero pointers
455 TNodeColumn& column = bot_column->second;
457 // compute bottom node parameters
458 gp_XYZ paramHint(-1,-1,-1);
459 if ( prevBNode.IsNeighbor( tBotNode ))
460 paramHint = prevBNode.GetParams();
461 if ( !myBlock.ComputeParameters( tBotNode.GetCoords(), tBotNode.ChangeParams(),
462 ID_BOT_FACE, paramHint ))
463 return error(TCom("Can't compute normalized parameters for node ")
464 << tBotNode.myNode->GetID() << " on the face #"
465 << myBlock.SubMesh( ID_BOT_FACE )->GetId() );
466 prevBNode = tBotNode;
468 myShapeXYZ[ ID_BOT_FACE ] = tBotNode.GetCoords();
469 gp_XYZ botParams = tBotNode.GetParams();
471 // compute top node parameters
472 myShapeXYZ[ ID_TOP_FACE ] = gpXYZ( column.back() );
473 gp_XYZ topParams = botParams;
475 if ( column.size() > 2 ) {
476 gp_Pnt topCoords = myShapeXYZ[ ID_TOP_FACE ];
477 if ( !myBlock.ComputeParameters( topCoords, topParams, ID_TOP_FACE, topParams ))
478 return error(TCom("Can't compute normalized parameters ")
479 << "for node " << column.back()->GetID()
480 << " on the face #"<< column.back()->getshapeId() );
484 TNodeColumn::iterator columnNodes = column.begin();
485 for ( int z = 0; columnNodes != column.end(); ++columnNodes, ++z)
487 const SMDS_MeshNode* & node = *columnNodes;
488 if ( node ) continue; // skip bottom or top node
490 // params of a node to create
491 double rz = (double) z / (double) ( column.size() - 1 );
492 gp_XYZ params = botParams * ( 1 - rz ) + topParams * rz;
494 // set coords on all faces and nodes
495 const int nbSideFaces = 4;
496 int sideFaceIDs[nbSideFaces] = { SMESH_Block::ID_Fx0z,
497 SMESH_Block::ID_Fx1z,
498 SMESH_Block::ID_F0yz,
499 SMESH_Block::ID_F1yz };
500 for ( int iF = 0; iF < nbSideFaces; ++iF )
501 if ( !setFaceAndEdgesXYZ( sideFaceIDs[ iF ], params, z ))
504 // compute coords for a new node
506 if ( !SMESH_Block::ShellPoint( params, myShapeXYZ, coords ))
507 return error("Can't compute coordinates by normalized parameters");
509 SHOWYXZ("TOPFacePoint ",myShapeXYZ[ ID_TOP_FACE]);
510 SHOWYXZ("BOT Node "<< tBotNode.myNode->GetID(),gpXYZ(tBotNode.myNode));
511 SHOWYXZ("ShellPoint ",coords);
514 node = meshDS->AddNode( coords.X(), coords.Y(), coords.Z() );
515 meshDS->SetNodeInVolume( node, volumeID );
517 } // loop on bottom nodes
522 SMESHDS_SubMesh* smDS = myBlock.SubMeshDS( ID_BOT_FACE );
523 if ( !smDS ) return error(COMPERR_BAD_INPUT_MESH, "Null submesh");
525 // loop on bottom mesh faces
526 SMDS_ElemIteratorPtr faceIt = smDS->GetElements();
527 while ( faceIt->more() )
529 const SMDS_MeshElement* face = faceIt->next();
530 if ( !face || face->GetType() != SMDSAbs_Face )
532 int nbNodes = face->NbNodes();
533 if ( face->IsQuadratic() )
536 // find node columns for each node
537 vector< const TNodeColumn* > columns( nbNodes );
538 for ( int i = 0; i < nbNodes; ++i )
540 const SMDS_MeshNode* n = face->GetNode( i );
541 if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ) {
542 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.find( n );
543 if ( bot_column == myBotToColumnMap.end() )
544 return error(TCom("No nodes found above node ") << n->GetID() );
545 columns[ i ] = & bot_column->second;
548 columns[ i ] = myBlock.GetNodeColumn( n );
550 return error(TCom("No side nodes found above node ") << n->GetID() );
554 AddPrisms( columns, myHelper );
556 } // loop on bottom mesh faces
559 myBotToColumnMap.clear();
566 //=======================================================================
567 //function : Evaluate
569 //=======================================================================
571 bool StdMeshers_Prism_3D::Evaluate(SMESH_Mesh& theMesh,
572 const TopoDS_Shape& theShape,
573 MapShapeNbElems& aResMap)
575 // find face contains only triangles
576 vector < SMESH_subMesh * >meshFaces;
577 TopTools_SequenceOfShape aFaces;
578 int NumBase = 0, i = 0, NbQFs = 0;
579 for (TopExp_Explorer exp(theShape, TopAbs_FACE); exp.More(); exp.Next()) {
581 aFaces.Append(exp.Current());
582 SMESH_subMesh *aSubMesh = theMesh.GetSubMesh(exp.Current());
583 meshFaces.push_back(aSubMesh);
584 MapShapeNbElemsItr anIt = aResMap.find(meshFaces[i-1]);
585 if( anIt==aResMap.end() ) {
586 SMESH_ComputeErrorPtr& smError = aSubMesh->GetComputeError();
587 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
590 std::vector<int> aVec = (*anIt).second;
591 int nbtri = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
592 int nbqua = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
593 if( nbtri==0 && nbqua>0 ) {
602 std::vector<int> aResVec(SMDSEntity_Last);
603 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
604 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
605 aResMap.insert(std::make_pair(sm,aResVec));
606 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
607 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
611 if(NumBase==0) NumBase = 1; // only quads => set 1 faces as base
613 // find number of 1d elems for base face
615 TopTools_MapOfShape Edges1;
616 for (TopExp_Explorer exp(aFaces.Value(NumBase), TopAbs_EDGE); exp.More(); exp.Next()) {
617 Edges1.Add(exp.Current());
618 SMESH_subMesh *sm = theMesh.GetSubMesh(exp.Current());
620 MapShapeNbElemsItr anIt = aResMap.find(sm);
621 if( anIt == aResMap.end() ) continue;
622 std::vector<int> aVec = (*anIt).second;
623 nb1d += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
626 // find face opposite to base face
628 for(i=1; i<=6; i++) {
629 if(i==NumBase) continue;
630 bool IsOpposite = true;
631 for(TopExp_Explorer exp(aFaces.Value(i), TopAbs_EDGE); exp.More(); exp.Next()) {
632 if( Edges1.Contains(exp.Current()) ) {
642 // find number of 2d elems on side faces
644 for(i=1; i<=6; i++) {
645 if( i==OppNum || i==NumBase ) continue;
646 MapShapeNbElemsItr anIt = aResMap.find( meshFaces[i-1] );
647 if( anIt == aResMap.end() ) continue;
648 std::vector<int> aVec = (*anIt).second;
649 nb2d += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
652 MapShapeNbElemsItr anIt = aResMap.find( meshFaces[NumBase-1] );
653 std::vector<int> aVec = (*anIt).second;
654 bool IsQuadratic = (aVec[SMDSEntity_Quad_Triangle]>aVec[SMDSEntity_Triangle]) ||
655 (aVec[SMDSEntity_Quad_Quadrangle]>aVec[SMDSEntity_Quadrangle]);
656 int nb2d_face0_3 = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
657 int nb2d_face0_4 = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
658 int nb0d_face0 = aVec[SMDSEntity_Node];
659 int nb1d_face0_int = ( nb2d_face0_3*3 + nb2d_face0_4*4 - nb1d ) / 2;
661 std::vector<int> aResVec(SMDSEntity_Last);
662 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
664 aResVec[SMDSEntity_Quad_Penta] = nb2d_face0_3 * ( nb2d/nb1d );
665 aResVec[SMDSEntity_Quad_Hexa] = nb2d_face0_4 * ( nb2d/nb1d );
666 aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
669 aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
670 aResVec[SMDSEntity_Penta] = nb2d_face0_3 * ( nb2d/nb1d );
671 aResVec[SMDSEntity_Hexa] = nb2d_face0_4 * ( nb2d/nb1d );
673 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
674 aResMap.insert(std::make_pair(sm,aResVec));
680 //================================================================================
682 * \brief Create prisms
683 * \param columns - columns of nodes generated from nodes of a mesh face
684 * \param helper - helper initialized by mesh and shape to add prisms to
686 //================================================================================
688 void StdMeshers_Prism_3D::AddPrisms( vector<const TNodeColumn*> & columns,
689 SMESH_MesherHelper* helper)
691 int nbNodes = columns.size();
692 int nbZ = columns[0]->size();
693 if ( nbZ < 2 ) return;
695 // find out orientation
696 bool isForward = true;
697 SMDS_VolumeTool vTool;
701 SMDS_VolumeOfNodes tmpPenta ( (*columns[0])[z-1], // bottom
704 (*columns[0])[z], // top
707 vTool.Set( &tmpPenta );
708 isForward = vTool.IsForward();
712 SMDS_VolumeOfNodes tmpHex( (*columns[0])[z-1], (*columns[1])[z-1], // bottom
713 (*columns[2])[z-1], (*columns[3])[z-1],
714 (*columns[0])[z], (*columns[1])[z], // top
715 (*columns[2])[z], (*columns[3])[z] );
716 vTool.Set( &tmpHex );
717 isForward = vTool.IsForward();
721 const int di = (nbNodes+1) / 3;
722 SMDS_VolumeOfNodes tmpVol ( (*columns[0] )[z-1],
723 (*columns[di] )[z-1],
724 (*columns[2*di])[z-1],
727 (*columns[2*di])[z] );
728 vTool.Set( &tmpVol );
729 isForward = vTool.IsForward();
732 // vertical loop on columns
734 helper->SetElementsOnShape( true );
738 case 3: { // ---------- pentahedra
739 const int i1 = isForward ? 1 : 2;
740 const int i2 = isForward ? 2 : 1;
741 for ( z = 1; z < nbZ; ++z )
742 helper->AddVolume( (*columns[0 ])[z-1], // bottom
745 (*columns[0 ])[z], // top
750 case 4: { // ---------- hexahedra
751 const int i1 = isForward ? 1 : 3;
752 const int i3 = isForward ? 3 : 1;
753 for ( z = 1; z < nbZ; ++z )
754 helper->AddVolume( (*columns[0])[z-1], (*columns[i1])[z-1], // bottom
755 (*columns[2])[z-1], (*columns[i3])[z-1],
756 (*columns[0])[z], (*columns[i1])[z], // top
757 (*columns[2])[z], (*columns[i3])[z] );
760 case 6: { // ---------- octahedra
761 const int iBase1 = isForward ? -1 : 0;
762 const int iBase2 = isForward ? 0 :-1;
763 for ( z = 1; z < nbZ; ++z )
764 helper->AddVolume( (*columns[0])[z+iBase1], (*columns[1])[z+iBase1], // bottom or top
765 (*columns[2])[z+iBase1], (*columns[3])[z+iBase1],
766 (*columns[4])[z+iBase1], (*columns[5])[z+iBase1],
767 (*columns[0])[z+iBase2], (*columns[1])[z+iBase2], // top or bottom
768 (*columns[2])[z+iBase2], (*columns[3])[z+iBase2],
769 (*columns[4])[z+iBase2], (*columns[5])[z+iBase2] );
772 default: // ---------- polyhedra
773 vector<int> quantities( 2 + nbNodes, 4 );
774 quantities[0] = quantities[1] = nbNodes;
775 columns.resize( nbNodes + 1 );
776 columns[ nbNodes ] = columns[ 0 ];
777 const int i1 = isForward ? 1 : 3;
778 const int i3 = isForward ? 3 : 1;
779 const int iBase1 = isForward ? -1 : 0;
780 const int iBase2 = isForward ? 0 :-1;
781 vector<const SMDS_MeshNode*> nodes( 2*nbNodes + 4*nbNodes);
782 for ( z = 1; z < nbZ; ++z )
784 for ( int i = 0; i < nbNodes; ++i ) {
785 nodes[ i ] = (*columns[ i ])[z+iBase1]; // bottom or top
786 nodes[ 2*nbNodes-i-1 ] = (*columns[ i ])[z+iBase2]; // top or bottom
788 int di = 2*nbNodes + 4*i;
789 nodes[ di+0 ] = (*columns[i ])[z ];
790 nodes[ di+i1] = (*columns[i+1])[z ];
791 nodes[ di+2 ] = (*columns[i+1])[z-1];
792 nodes[ di+i3] = (*columns[i ])[z-1];
794 helper->AddPolyhedralVolume( nodes, quantities );
797 } // switch ( nbNodes )
800 //================================================================================
802 * \brief Find correspondence between bottom and top nodes
803 * If elements on the bottom and top faces are topologically different,
804 * and projection is possible and allowed, perform the projection
805 * \retval bool - is a success or not
807 //================================================================================
809 bool StdMeshers_Prism_3D::assocOrProjBottom2Top()
811 SMESH_subMesh * botSM = myBlock.SubMesh( ID_BOT_FACE );
812 SMESH_subMesh * topSM = myBlock.SubMesh( ID_TOP_FACE );
814 SMESHDS_SubMesh * botSMDS = botSM->GetSubMeshDS();
815 SMESHDS_SubMesh * topSMDS = topSM->GetSubMeshDS();
817 if ( !botSMDS || botSMDS->NbElements() == 0 )
818 return error(TCom("No elememts on face #") << botSM->GetId());
820 bool needProject = false;
822 botSMDS->NbElements() != topSMDS->NbElements() ||
823 botSMDS->NbNodes() != topSMDS->NbNodes())
825 MESSAGE("nb elem bot " << botSMDS->NbElements() << " top " << topSMDS->NbElements());
826 MESSAGE("nb node bot " << botSMDS->NbNodes() << " top " << topSMDS->NbNodes());
827 if ( myBlock.HasNotQuadElemOnTop() )
828 return error(TCom("Mesh on faces #") << botSM->GetId()
829 <<" and #"<< topSM->GetId() << " seems different" );
833 if ( 0/*needProject && !myProjectTriangles*/ )
834 return error(TCom("Mesh on faces #") << botSM->GetId()
835 <<" and #"<< topSM->GetId() << " seems different" );
836 ///RETURN_BAD_RESULT("Need to project but not allowed");
840 return projectBottomToTop();
843 TopoDS_Face botFace = TopoDS::Face( myBlock.Shape( ID_BOT_FACE ));
844 TopoDS_Face topFace = TopoDS::Face( myBlock.Shape( ID_TOP_FACE ));
845 // associate top and bottom faces
846 TAssocTool::TShapeShapeMap shape2ShapeMap;
847 if ( !TAssocTool::FindSubShapeAssociation( botFace, myBlock.Mesh(),
848 topFace, myBlock.Mesh(),
850 return error(TCom("Topology of faces #") << botSM->GetId()
851 <<" and #"<< topSM->GetId() << " seems different" );
853 // Find matching nodes of top and bottom faces
855 if ( ! TAssocTool::FindMatchingNodesOnFaces( botFace, myBlock.Mesh(),
856 topFace, myBlock.Mesh(),
857 shape2ShapeMap, n2nMap ))
858 return error(TCom("Mesh on faces #") << botSM->GetId()
859 <<" and #"<< topSM->GetId() << " seems different" );
861 // Fill myBotToColumnMap
863 int zSize = myBlock.VerticalSize();
865 TNodeNodeMap::iterator bN_tN = n2nMap.begin();
866 for ( ; bN_tN != n2nMap.end(); ++bN_tN )
868 const SMDS_MeshNode* botNode = bN_tN->first;
869 const SMDS_MeshNode* topNode = bN_tN->second;
870 if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
871 continue; // wall columns are contained in myBlock
872 // create node column
874 TNode2ColumnMap::iterator bN_col =
875 myBotToColumnMap.insert( make_pair ( bN, TNodeColumn() )).first;
876 TNodeColumn & column = bN_col->second;
877 column.resize( zSize );
878 column.front() = botNode;
879 column.back() = topNode;
884 //================================================================================
886 * \brief Remove quadrangles from the top face and
887 * create triangles there by projection from the bottom
888 * \retval bool - a success or not
890 //================================================================================
892 bool StdMeshers_Prism_3D::projectBottomToTop()
894 SMESH_subMesh * botSM = myBlock.SubMesh( ID_BOT_FACE );
895 SMESH_subMesh * topSM = myBlock.SubMesh( ID_TOP_FACE );
897 SMESHDS_SubMesh * botSMDS = botSM->GetSubMeshDS();
898 SMESHDS_SubMesh * topSMDS = topSM->GetSubMeshDS();
901 topSM->ComputeStateEngine( SMESH_subMesh::CLEAN );
903 SMESHDS_Mesh* meshDS = myBlock.MeshDS();
904 int shapeID = myHelper->GetSubShapeID();
905 int topFaceID = meshDS->ShapeToIndex( topSM->GetSubShape() );
907 // Fill myBotToColumnMap
909 int zSize = myBlock.VerticalSize();
911 SMDS_NodeIteratorPtr nIt = botSMDS->GetNodes();
912 while ( nIt->more() )
914 const SMDS_MeshNode* botNode = nIt->next();
915 if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
917 // compute bottom node params
919 gp_XYZ paramHint(-1,-1,-1);
920 if ( prevTNode.IsNeighbor( bN ))
921 paramHint = prevTNode.GetParams();
922 if ( !myBlock.ComputeParameters( bN.GetCoords(), bN.ChangeParams(),
923 ID_BOT_FACE, paramHint ))
924 return error(TCom("Can't compute normalized parameters for node ")
925 << botNode->GetID() << " on the face #"<< botSM->GetId() );
927 // compute top node coords
928 gp_XYZ topXYZ; gp_XY topUV;
929 if ( !myBlock.FacePoint( ID_TOP_FACE, bN.GetParams(), topXYZ ) ||
930 !myBlock.FaceUV ( ID_TOP_FACE, bN.GetParams(), topUV ))
931 return error(TCom("Can't compute coordinates "
932 "by normalized parameters on the face #")<< topSM->GetId() );
933 SMDS_MeshNode * topNode = meshDS->AddNode( topXYZ.X(),topXYZ.Y(),topXYZ.Z() );
934 meshDS->SetNodeOnFace( topNode, topFaceID, topUV.X(), topUV.Y() );
935 // create node column
936 TNode2ColumnMap::iterator bN_col =
937 myBotToColumnMap.insert( make_pair ( bN, TNodeColumn() )).first;
938 TNodeColumn & column = bN_col->second;
939 column.resize( zSize );
940 column.front() = botNode;
941 column.back() = topNode;
946 // loop on bottom mesh faces
947 SMDS_ElemIteratorPtr faceIt = botSMDS->GetElements();
948 while ( faceIt->more() )
950 const SMDS_MeshElement* face = faceIt->next();
951 if ( !face || face->GetType() != SMDSAbs_Face )
953 int nbNodes = face->NbNodes();
954 if ( face->IsQuadratic() )
957 // find top node in columns for each bottom node
958 vector< const SMDS_MeshNode* > nodes( nbNodes );
959 for ( int i = 0; i < nbNodes; ++i )
961 const SMDS_MeshNode* n = face->GetNode( nbNodes - i - 1 );
962 if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ) {
963 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.find( n );
964 if ( bot_column == myBotToColumnMap.end() )
965 return error(TCom("No nodes found above node ") << n->GetID() );
966 nodes[ i ] = bot_column->second.back();
969 const TNodeColumn* column = myBlock.GetNodeColumn( n );
971 return error(TCom("No side nodes found above node ") << n->GetID() );
972 nodes[ i ] = column->back();
975 // create a face, with reversed orientation
976 SMDS_MeshElement* newFace = 0;
980 newFace = myHelper->AddFace(nodes[0], nodes[1], nodes[2]);
984 newFace = myHelper->AddFace( nodes[0], nodes[1], nodes[2], nodes[3] );
988 newFace = meshDS->AddPolygonalFace( nodes );
990 if ( newFace && shapeID > 0 )
991 meshDS->SetMeshElementOnShape( newFace, shapeID );
997 //================================================================================
999 * \brief Set projection coordinates of a node to a face and it's sub-shapes
1000 * \param faceID - the face given by in-block ID
1001 * \param params - node normalized parameters
1002 * \retval bool - is a success
1004 //================================================================================
1006 bool StdMeshers_Prism_3D::setFaceAndEdgesXYZ( const int faceID, const gp_XYZ& params, int z )
1008 // find base and top edges of the face
1009 enum { BASE = 0, TOP, LEFT, RIGHT };
1010 vector< int > edgeVec; // 0-base, 1-top
1011 SMESH_Block::GetFaceEdgesIDs( faceID, edgeVec );
1013 myBlock.EdgePoint( edgeVec[ BASE ], params, myShapeXYZ[ edgeVec[ BASE ]]);
1014 myBlock.EdgePoint( edgeVec[ TOP ], params, myShapeXYZ[ edgeVec[ TOP ]]);
1016 SHOWYXZ("\nparams ", params);
1017 SHOWYXZ("TOP is " <<edgeVec[ TOP ], myShapeXYZ[ edgeVec[ TOP]]);
1018 SHOWYXZ("BASE is "<<edgeVec[ BASE], myShapeXYZ[ edgeVec[ BASE]]);
1020 if ( faceID == SMESH_Block::ID_Fx0z || faceID == SMESH_Block::ID_Fx1z )
1022 myBlock.EdgePoint( edgeVec[ LEFT ], params, myShapeXYZ[ edgeVec[ LEFT ]]);
1023 myBlock.EdgePoint( edgeVec[ RIGHT ], params, myShapeXYZ[ edgeVec[ RIGHT ]]);
1025 SHOWYXZ("VER "<<edgeVec[ LEFT], myShapeXYZ[ edgeVec[ LEFT]]);
1026 SHOWYXZ("VER "<<edgeVec[ RIGHT], myShapeXYZ[ edgeVec[ RIGHT]]);
1028 myBlock.FacePoint( faceID, params, myShapeXYZ[ faceID ]);
1029 SHOWYXZ("FacePoint "<<faceID, myShapeXYZ[ faceID]);
1034 //================================================================================
1036 * \brief Return true if this node and other one belong to one face
1038 //================================================================================
1040 bool TNode::IsNeighbor( const TNode& other ) const
1042 if ( !other.myNode || !myNode ) return false;
1044 SMDS_ElemIteratorPtr fIt = other.myNode->GetInverseElementIterator(SMDSAbs_Face);
1045 while ( fIt->more() )
1046 if ( fIt->next()->GetNodeIndex( myNode ) >= 0 )
1051 //================================================================================
1053 * \brief Constructor. Initialization is needed
1055 //================================================================================
1057 StdMeshers_PrismAsBlock::StdMeshers_PrismAsBlock()
1062 StdMeshers_PrismAsBlock::~StdMeshers_PrismAsBlock()
1066 void StdMeshers_PrismAsBlock::Clear()
1069 myShapeIDMap.Clear();
1073 delete mySide; mySide = 0;
1075 myParam2ColumnMaps.clear();
1076 myShapeIndex2ColumnMap.clear();
1079 //================================================================================
1081 * \brief Initialization.
1082 * \param helper - helper loaded with mesh and 3D shape
1083 * \param shape3D - a closed shell or solid
1084 * \retval bool - false if a mesh or a shape are KO
1086 //================================================================================
1088 bool StdMeshers_PrismAsBlock::Init(SMESH_MesherHelper* helper,
1089 const TopoDS_Shape& shape3D)
1092 delete mySide; mySide = 0;
1094 vector< TSideFace* > sideFaces( NB_WALL_FACES, 0 );
1095 vector< pair< double, double> > params ( NB_WALL_FACES );
1096 mySide = new TSideFace( sideFaces, params );
1099 SMESHDS_Mesh* meshDS = myHelper->GetMeshDS();
1101 SMESH_Block::init();
1102 myShapeIDMap.Clear();
1103 myShapeIndex2ColumnMap.clear();
1105 int wallFaceIds[ NB_WALL_FACES ] = { // to walk around a block
1106 SMESH_Block::ID_Fx0z, SMESH_Block::ID_F1yz,
1107 SMESH_Block::ID_Fx1z, SMESH_Block::ID_F0yz
1110 myError = SMESH_ComputeError::New();
1112 // -------------------------------------------------------------
1113 // Look for top and bottom faces: not quadrangle ones or meshed
1114 // with not quadrangle elements
1115 // -------------------------------------------------------------
1117 list< SMESH_subMesh* > notQuadGeomSubMesh;
1118 list< SMESH_subMesh* > notQuadElemSubMesh;
1121 SMESH_subMesh* mainSubMesh = myHelper->GetMesh()->GetSubMeshContaining( shape3D );
1122 if ( !mainSubMesh ) return error(COMPERR_BAD_INPUT_MESH,"Null submesh of shape3D");
1124 // analyse face submeshes
1125 SMESH_subMeshIteratorPtr smIt = mainSubMesh->getDependsOnIterator(false,false);
1126 while ( smIt->more() )
1128 SMESH_subMesh* sm = smIt->next();
1129 const TopoDS_Shape& face = sm->GetSubShape();
1130 if ( face.ShapeType() != TopAbs_FACE )
1134 // is quadrangle face?
1135 list< TopoDS_Edge > orderedEdges;
1136 list< int > nbEdgesInWires;
1138 int nbWires = GetOrderedEdges( TopoDS::Face( face ),
1139 V000, orderedEdges, nbEdgesInWires );
1140 if ( nbWires != 1 || nbEdgesInWires.front() != 4 )
1141 notQuadGeomSubMesh.push_back( sm );
1143 // look for not quadrangle mesh elements
1144 if ( SMESHDS_SubMesh* smDS = sm->GetSubMeshDS() ) {
1145 bool hasNotQuad = false;
1146 SMDS_ElemIteratorPtr eIt = smDS->GetElements();
1147 while ( eIt->more() && !hasNotQuad ) {
1148 const SMDS_MeshElement* elem = eIt->next();
1149 if ( elem->GetType() == SMDSAbs_Face ) {
1150 int nbNodes = elem->NbNodes();
1151 if ( elem->IsQuadratic() )
1153 hasNotQuad = ( nbNodes != 4 );
1157 notQuadElemSubMesh.push_back( sm );
1160 return error(COMPERR_BAD_INPUT_MESH,TCom("Not meshed face #")<<sm->GetId());
1162 // check if a quadrangle face is meshed with a quadranglar grid
1163 if ( notQuadGeomSubMesh.back() != sm &&
1164 notQuadElemSubMesh.back() != sm )
1166 // count nb edges on face sides
1167 vector< int > nbEdges;
1168 nbEdges.reserve( nbEdgesInWires.front() );
1169 for ( list< TopoDS_Edge >::iterator edge = orderedEdges.begin();
1170 edge != orderedEdges.end(); ++edge )
1172 if ( SMESHDS_SubMesh* smDS = meshDS->MeshElements( *edge ))
1173 nbEdges.push_back ( smDS->NbElements() );
1175 nbEdges.push_back ( 0 );
1177 int nbQuads = sm->GetSubMeshDS()->NbElements();
1178 if ( nbEdges[0] * nbEdges[1] != nbQuads ||
1179 nbEdges[0] != nbEdges[2] ||
1180 nbEdges[1] != nbEdges[3] )
1181 notQuadElemSubMesh.push_back( sm );
1185 // ----------------------------------------------------------------------
1186 // Analyse mesh and topology of faces: choose the bottom submesh.
1187 // If there are not quadrangle geom faces, they are top and bottom ones.
1188 // Not quadrangle geom faces must be only on top and bottom.
1189 // ----------------------------------------------------------------------
1191 SMESH_subMesh * botSM = 0;
1192 SMESH_subMesh * topSM = 0;
1194 int nbNotQuad = notQuadGeomSubMesh.size();
1195 int nbNotQuadMeshed = notQuadElemSubMesh.size();
1196 bool hasNotQuad = ( nbNotQuad || nbNotQuadMeshed );
1199 if ( nbNotQuadMeshed > 2 )
1201 return error(COMPERR_BAD_INPUT_MESH,
1202 TCom("More than 2 faces with not quadrangle elements: ")
1205 int nbQuasiQuads = 0;
1206 if ( nbNotQuad > 0 && nbNotQuad != 2 )
1208 // Issue 0020843 - one of side faces is quasi-quadrilateral.
1209 // Remove from notQuadGeomSubMesh faces meshed with regular grid
1210 nbQuasiQuads = removeQuasiQuads( notQuadGeomSubMesh );
1211 nbNotQuad -= nbQuasiQuads;
1212 if ( nbNotQuad > 0 && nbNotQuad != 2 )
1213 return error(COMPERR_BAD_SHAPE,
1214 TCom("More than 2 not quadrilateral faces: ")
1218 // get found submeshes
1221 if ( nbNotQuadMeshed > 0 ) botSM = notQuadElemSubMesh.front();
1222 else botSM = notQuadGeomSubMesh.front();
1223 if ( nbNotQuadMeshed > 1 ) topSM = notQuadElemSubMesh.back();
1224 else if ( nbNotQuad > 1 ) topSM = notQuadGeomSubMesh.back();
1226 // detect other bad cases
1227 if ( nbNotQuad == 2 && nbNotQuadMeshed > 0 ) {
1229 if ( nbNotQuadMeshed == 1 )
1230 ok = ( find( notQuadGeomSubMesh.begin(),
1231 notQuadGeomSubMesh.end(), botSM ) != notQuadGeomSubMesh.end() );
1233 ok = ( notQuadGeomSubMesh == notQuadElemSubMesh );
1235 return error(COMPERR_BAD_INPUT_MESH, "Side face meshed with not quadrangle elements");
1238 myNotQuadOnTop = ( nbNotQuadMeshed > 1 );
1239 MESSAGE("myNotQuadOnTop " << myNotQuadOnTop << " nbNotQuadMeshed " << nbNotQuadMeshed);
1241 // ----------------------------------------------------------
1243 if ( nbNotQuad == 0 ) // Standard block of 6 quadrangle faces ?
1245 // SMESH_Block will perform geometry analysis, we need just to find 2
1246 // connected vertices on top and bottom
1248 TopoDS_Vertex Vbot, Vtop;
1249 if ( nbNotQuadMeshed > 0 ) // Look for vertices
1251 TopTools_IndexedMapOfShape edgeMap;
1252 TopExp::MapShapes( botSM->GetSubShape(), TopAbs_EDGE, edgeMap );
1253 // vertex 1 is any vertex of the bottom face
1254 Vbot = TopExp::FirstVertex( TopoDS::Edge( edgeMap( 1 )));
1255 // vertex 2 is end vertex of edge sharing Vbot and not belonging to the bottom face
1256 TopTools_ListIteratorOfListOfShape ancestIt = Mesh()->GetAncestors( Vbot );
1257 for ( ; Vtop.IsNull() && ancestIt.More(); ancestIt.Next() )
1259 const TopoDS_Shape & ancestor = ancestIt.Value();
1260 if ( ancestor.ShapeType() == TopAbs_EDGE && !edgeMap.FindIndex( ancestor ))
1262 TopoDS_Vertex V1, V2;
1263 TopExp::Vertices( TopoDS::Edge( ancestor ), V1, V2);
1264 if ( Vbot.IsSame ( V1 )) Vtop = V2;
1265 else if ( Vbot.IsSame ( V2 )) Vtop = V1;
1266 // check that Vtop belongs to shape3D
1267 TopExp_Explorer exp( shape3D, TopAbs_VERTEX );
1268 for ( ; exp.More(); exp.Next() )
1269 if ( Vtop.IsSame( exp.Current() ))
1276 // get shell from shape3D
1278 TopExp_Explorer exp( shape3D, TopAbs_SHELL );
1280 for ( ; exp.More(); exp.Next(), ++nbShell )
1281 shell = TopoDS::Shell( exp.Current() );
1282 // if ( nbShell != 1 )
1283 // RETURN_BAD_RESULT("There must be 1 shell in the block");
1285 // Load geometry in SMESH_Block
1286 if ( !SMESH_Block::FindBlockShapes( shell, Vbot, Vtop, myShapeIDMap )) {
1288 return error(COMPERR_BAD_SHAPE, "Can't detect top and bottom of a prism");
1291 if ( !botSM ) botSM = Mesh()->GetSubMeshContaining( myShapeIDMap( ID_BOT_FACE ));
1292 if ( !topSM ) topSM = Mesh()->GetSubMeshContaining( myShapeIDMap( ID_TOP_FACE ));
1295 } // end Standard block of 6 quadrangle faces
1296 // --------------------------------------------------------
1298 // Here the top and bottom faces are found
1299 if ( nbNotQuadMeshed == 2 ) // roughly check correspondence of horiz meshes
1301 // SMESHDS_SubMesh* topSMDS = topSM->GetSubMeshDS();
1302 // SMESHDS_SubMesh* botSMDS = botSM->GetSubMeshDS();
1303 // if ( topSMDS->NbNodes() != botSMDS->NbNodes() ||
1304 // topSMDS->NbElements() != botSMDS->NbElements() )
1305 // RETURN_BAD_RESULT("Top mesh doesn't correspond to bottom one");
1308 // ---------------------------------------------------------
1309 // If there are not quadrangle geom faces, we emulate
1310 // a block of 6 quadrangle faces.
1311 // Load SMESH_Block with faces and edges geometry
1312 // ---------------------------------------------------------
1315 // find vertex 000 - the one with smallest coordinates (for easy DEBUG :-)
1317 double minVal = DBL_MAX, minX, val;
1318 for ( TopExp_Explorer exp( botSM->GetSubShape(), TopAbs_VERTEX );
1319 exp.More(); exp.Next() )
1321 const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
1322 gp_Pnt P = BRep_Tool::Pnt( v );
1323 val = P.X() + P.Y() + P.Z();
1324 if ( val < minVal || ( val == minVal && P.X() < minX )) {
1331 // Get ordered bottom edges
1332 list< TopoDS_Edge > orderedEdges;
1334 SMESH_Block::GetOrderedEdges( TopoDS::Face( botSM->GetSubShape().Reversed() ),
1335 V000, orderedEdges, nbEInW );
1336 // if ( nbEInW.size() != 1 )
1337 // RETURN_BAD_RESULT("Wrong prism geometry");
1339 // Get Wall faces corresponding to the ordered bottom edges
1340 list< TopoDS_Face > wallFaces;
1341 if ( !GetWallFaces( Mesh(), shape3D, botSM->GetSubShape(), orderedEdges, nbEInW, wallFaces))
1342 return error(COMPERR_BAD_SHAPE, "Can't find side faces");
1344 // check that the found top and bottom faces are opposite
1346 for (TopExp_Explorer edge(botSM->GetSubShape(), TopAbs_EDGE); edge.More(); edge.Next())
1347 if ( helper->IsSubShape( edge.Current(), topSM->GetSubShape() ))
1348 return error(notQuadGeomSubMesh.empty() ? COMPERR_BAD_INPUT_MESH : COMPERR_BAD_SHAPE,
1349 "Non-quadrilateral faces are not opposite");
1352 // Protect from a distorted block (test 3D_mesh_HEXA3D/B7 on 32bit platform)
1353 // check that all wall faces have an edge common with the top face
1355 list< TopoDS_Face >::iterator faceIt = wallFaces.begin();
1356 for ( ; faceIt != wallFaces.end(); ++faceIt )
1358 bool hasCommon = false;
1359 for (TopExp_Explorer edge(*faceIt, TopAbs_EDGE); !hasCommon && edge.More(); edge.Next())
1360 if ( helper->IsSubShape( edge.Current(), topSM->GetSubShape() ))
1363 return error(COMPERR_BAD_SHAPE);
1367 // Find columns of wall nodes and calculate edges' lengths
1368 // --------------------------------------------------------
1370 myParam2ColumnMaps.clear();
1371 myParam2ColumnMaps.resize( orderedEdges.size() ); // total nb edges
1373 int iE, nbEdges = nbEInW.front(); // nb outer edges
1374 vector< double > edgeLength( nbEdges );
1375 map< double, int > len2edgeMap;
1377 list< TopoDS_Edge >::iterator edgeIt = orderedEdges.begin();
1378 list< TopoDS_Face >::iterator faceIt = wallFaces.begin();
1379 for ( iE = 0; iE < nbEdges; ++edgeIt, ++faceIt )
1381 TParam2ColumnMap & faceColumns = myParam2ColumnMaps[ iE ];
1382 if ( !myHelper->LoadNodeColumns( faceColumns, *faceIt, *edgeIt, meshDS ))
1383 return error(COMPERR_BAD_INPUT_MESH, TCom("Can't find regular quadrangle mesh ")
1384 << "on a side face #" << MeshDS()->ShapeToIndex( *faceIt ));
1386 SHOWYXZ("\np1 F "<<iE, gpXYZ(faceColumns.begin()->second.front() ));
1387 SHOWYXZ("p2 F "<<iE, gpXYZ(faceColumns.rbegin()->second.front() ));
1388 SHOWYXZ("V First "<<iE, BRep_Tool::Pnt( TopExp::FirstVertex(*edgeIt,true )));
1390 edgeLength[ iE ] = SMESH_Algo::EdgeLength( *edgeIt );
1392 if ( nbEdges < NB_WALL_FACES ) // fill map used to split faces
1394 SMESHDS_SubMesh* smDS = meshDS->MeshElements( *edgeIt);
1396 return error(COMPERR_BAD_INPUT_MESH, TCom("Null submesh on the edge #")
1397 << MeshDS()->ShapeToIndex( *edgeIt ));
1398 // assure length uniqueness
1399 edgeLength[ iE ] *= smDS->NbNodes() + edgeLength[ iE ] / ( 1000 + iE );
1400 len2edgeMap[ edgeLength[ iE ]] = iE;
1404 // Load columns of internal edges (forming holes)
1405 // and fill map ShapeIndex to TParam2ColumnMap for them
1406 for ( ; edgeIt != orderedEdges.end() ; ++edgeIt, ++faceIt )
1408 TParam2ColumnMap & faceColumns = myParam2ColumnMaps[ iE ];
1409 if ( !myHelper->LoadNodeColumns( faceColumns, *faceIt, *edgeIt, meshDS ))
1410 return error(COMPERR_BAD_INPUT_MESH, TCom("Can't find regular quadrangle mesh ")
1411 << "on a side face #" << MeshDS()->ShapeToIndex( *faceIt ));
1413 int id = MeshDS()->ShapeToIndex( *edgeIt );
1414 bool isForward = true; // meaningless for intenal wires
1415 myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1416 // columns for vertices
1418 const SMDS_MeshNode* n0 = faceColumns.begin()->second.front();
1419 id = n0->getshapeId();
1420 myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1422 const SMDS_MeshNode* n1 = faceColumns.rbegin()->second.front();
1423 id = n1->getshapeId();
1424 myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1425 // SHOWYXZ("\np1 F "<<iE, gpXYZ(faceColumns.begin()->second.front() ));
1426 // SHOWYXZ("p2 F "<<iE, gpXYZ(faceColumns.rbegin()->second.front() ));
1427 // SHOWYXZ("V First "<<iE, BRep_Tool::Pnt( TopExp::FirstVertex(*edgeIt,true )));
1431 // Create 4 wall faces of a block
1432 // -------------------------------
1434 if ( nbEdges <= NB_WALL_FACES ) // ************* Split faces if necessary
1436 map< int, int > iE2nbSplit;
1437 if ( nbEdges != NB_WALL_FACES ) // define how to split
1439 if ( len2edgeMap.size() != nbEdges )
1440 RETURN_BAD_RESULT("Uniqueness of edge lengths not assured");
1441 map< double, int >::reverse_iterator maxLen_i = len2edgeMap.rbegin();
1442 map< double, int >::reverse_iterator midLen_i = ++len2edgeMap.rbegin();
1443 double maxLen = maxLen_i->first;
1444 double midLen = ( len2edgeMap.size() == 1 ) ? 0 : midLen_i->first;
1445 switch ( nbEdges ) {
1446 case 1: // 0-th edge is split into 4 parts
1447 iE2nbSplit.insert( make_pair( 0, 4 )); break;
1448 case 2: // either the longest edge is split into 3 parts, or both edges into halves
1449 if ( maxLen / 3 > midLen / 2 ) {
1450 iE2nbSplit.insert( make_pair( maxLen_i->second, 3 ));
1453 iE2nbSplit.insert( make_pair( maxLen_i->second, 2 ));
1454 iE2nbSplit.insert( make_pair( midLen_i->second, 2 ));
1458 // split longest into halves
1459 iE2nbSplit.insert( make_pair( maxLen_i->second, 2 ));
1462 // Create TSideFace's
1463 faceIt = wallFaces.begin();
1464 edgeIt = orderedEdges.begin();
1466 for ( iE = 0; iE < nbEdges; ++edgeIt, ++faceIt )
1469 map< int, int >::iterator i_nb = iE2nbSplit.find( iE );
1470 if ( i_nb != iE2nbSplit.end() ) {
1472 int nbSplit = i_nb->second;
1473 vector< double > params;
1474 splitParams( nbSplit, &myParam2ColumnMaps[ iE ], params );
1475 const bool isForward =
1476 StdMeshers_PrismAsBlock::IsForwardEdge( myHelper->GetMeshDS(),
1477 myParam2ColumnMaps[iE],
1478 *edgeIt, SMESH_Block::ID_Fx0z );
1479 for ( int i = 0; i < nbSplit; ++i ) {
1480 double f = ( isForward ? params[ i ] : params[ nbSplit - i-1 ]);
1481 double l = ( isForward ? params[ i+1 ] : params[ nbSplit - i ]);
1482 TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1484 &myParam2ColumnMaps[ iE ], f, l );
1485 mySide->SetComponent( iSide++, comp );
1489 TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1491 &myParam2ColumnMaps[ iE ]);
1492 mySide->SetComponent( iSide++, comp );
1497 else { // **************************** Unite faces
1499 // unite first faces
1500 int nbExraFaces = nbEdges - 3;
1502 double u0 = 0, sumLen = 0;
1503 for ( iE = 0; iE < nbExraFaces; ++iE )
1504 sumLen += edgeLength[ iE ];
1506 vector< TSideFace* > components( nbExraFaces );
1507 vector< pair< double, double> > params( nbExraFaces );
1508 faceIt = wallFaces.begin();
1509 edgeIt = orderedEdges.begin();
1510 for ( iE = 0; iE < nbExraFaces; ++edgeIt, ++faceIt )
1512 components[ iE ] = new TSideFace( myHelper, wallFaceIds[ iSide ],
1514 &myParam2ColumnMaps[ iE ]);
1515 double u1 = u0 + edgeLength[ iE ] / sumLen;
1516 params[ iE ] = make_pair( u0 , u1 );
1520 mySide->SetComponent( iSide++, new TSideFace( components, params ));
1522 // fill the rest faces
1523 for ( ; iE < nbEdges; ++faceIt, ++edgeIt )
1525 TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1527 &myParam2ColumnMaps[ iE ]);
1528 mySide->SetComponent( iSide++, comp );
1534 // Fill geometry fields of SMESH_Block
1535 // ------------------------------------
1537 TopoDS_Face botF = TopoDS::Face( botSM->GetSubShape() );
1538 TopoDS_Face topF = TopoDS::Face( topSM->GetSubShape() );
1540 vector< int > botEdgeIdVec;
1541 SMESH_Block::GetFaceEdgesIDs( ID_BOT_FACE, botEdgeIdVec );
1543 bool isForward[NB_WALL_FACES] = { true, true, true, true };
1544 Adaptor2d_Curve2d* botPcurves[NB_WALL_FACES];
1545 Adaptor2d_Curve2d* topPcurves[NB_WALL_FACES];
1547 for ( int iF = 0; iF < NB_WALL_FACES; ++iF )
1549 TSideFace * sideFace = mySide->GetComponent( iF );
1551 RETURN_BAD_RESULT("NULL TSideFace");
1552 int fID = sideFace->FaceID();
1554 // fill myShapeIDMap
1555 if ( sideFace->InsertSubShapes( myShapeIDMap ) != 8 &&
1556 !sideFace->IsComplex())
1557 MESSAGE( ": Warning : InsertSubShapes() < 8 on side " << iF );
1559 // side faces geometry
1560 Adaptor2d_Curve2d* pcurves[NB_WALL_FACES];
1561 if ( !sideFace->GetPCurves( pcurves ))
1562 RETURN_BAD_RESULT("TSideFace::GetPCurves() failed");
1564 SMESH_Block::TFace& tFace = myFace[ fID - ID_FirstF ];
1565 tFace.Set( fID, sideFace->Surface(), pcurves, isForward );
1567 SHOWYXZ( endl<<"F "<< iF << " id " << fID << " FRW " << sideFace->IsForward(), sideFace->Value(0,0));
1568 // edges 3D geometry
1569 vector< int > edgeIdVec;
1570 SMESH_Block::GetFaceEdgesIDs( fID, edgeIdVec );
1571 for ( int isMax = 0; isMax < 2; ++isMax ) {
1573 int eID = edgeIdVec[ isMax ];
1574 SMESH_Block::TEdge& tEdge = myEdge[ eID - ID_FirstE ];
1575 tEdge.Set( eID, sideFace->HorizCurve(isMax), true);
1576 SHOWYXZ(eID<<" HOR"<<isMax<<"(0)", sideFace->HorizCurve(isMax)->Value(0));
1577 SHOWYXZ(eID<<" HOR"<<isMax<<"(1)", sideFace->HorizCurve(isMax)->Value(1));
1580 int eID = edgeIdVec[ isMax+2 ];
1581 SMESH_Block::TEdge& tEdge = myEdge[ eID - ID_FirstE ];
1582 tEdge.Set( eID, sideFace->VertiCurve(isMax), true);
1583 SHOWYXZ(eID<<" VER"<<isMax<<"(0)", sideFace->VertiCurve(isMax)->Value(0));
1584 SHOWYXZ(eID<<" VER"<<isMax<<"(1)", sideFace->VertiCurve(isMax)->Value(1));
1587 vector< int > vertexIdVec;
1588 SMESH_Block::GetEdgeVertexIDs( eID, vertexIdVec );
1589 myPnt[ vertexIdVec[0] - ID_FirstV ] = tEdge.GetCurve()->Value(0).XYZ();
1590 myPnt[ vertexIdVec[1] - ID_FirstV ] = tEdge.GetCurve()->Value(1).XYZ();
1593 // pcurves on horizontal faces
1594 for ( iE = 0; iE < NB_WALL_FACES; ++iE ) {
1595 if ( edgeIdVec[ BOTTOM_EDGE ] == botEdgeIdVec[ iE ] ) {
1596 botPcurves[ iE ] = sideFace->HorizPCurve( false, botF );
1597 topPcurves[ iE ] = sideFace->HorizPCurve( true, topF );
1601 //sideFace->dumpNodes( 4 ); // debug
1603 // horizontal faces geometry
1605 SMESH_Block::TFace& tFace = myFace[ ID_BOT_FACE - ID_FirstF ];
1606 tFace.Set( ID_BOT_FACE, new BRepAdaptor_Surface( botF ), botPcurves, isForward );
1607 SMESH_Block::Insert( botF, ID_BOT_FACE, myShapeIDMap );
1610 SMESH_Block::TFace& tFace = myFace[ ID_TOP_FACE - ID_FirstF ];
1611 tFace.Set( ID_TOP_FACE, new BRepAdaptor_Surface( topF ), topPcurves, isForward );
1612 SMESH_Block::Insert( topF, ID_TOP_FACE, myShapeIDMap );
1615 // Fill map ShapeIndex to TParam2ColumnMap
1616 // ----------------------------------------
1618 list< TSideFace* > fList;
1619 list< TSideFace* >::iterator fListIt;
1620 fList.push_back( mySide );
1621 for ( fListIt = fList.begin(); fListIt != fList.end(); ++fListIt)
1623 int nb = (*fListIt)->NbComponents();
1624 for ( int i = 0; i < nb; ++i ) {
1625 if ( TSideFace* comp = (*fListIt)->GetComponent( i ))
1626 fList.push_back( comp );
1628 if ( TParam2ColumnMap* cols = (*fListIt)->GetColumns()) {
1629 // columns for a base edge
1630 int id = MeshDS()->ShapeToIndex( (*fListIt)->BaseEdge() );
1631 bool isForward = (*fListIt)->IsForward();
1632 myShapeIndex2ColumnMap[ id ] = make_pair( cols, isForward );
1634 // columns for vertices
1635 const SMDS_MeshNode* n0 = cols->begin()->second.front();
1636 id = n0->getshapeId();
1637 myShapeIndex2ColumnMap[ id ] = make_pair( cols, isForward );
1639 const SMDS_MeshNode* n1 = cols->rbegin()->second.front();
1640 id = n1->getshapeId();
1641 myShapeIndex2ColumnMap[ id ] = make_pair( cols, !isForward );
1645 // gp_XYZ testPar(0.25, 0.25, 0), testCoord;
1646 // if ( !FacePoint( ID_BOT_FACE, testPar, testCoord ))
1647 // RETURN_BAD_RESULT("TEST FacePoint() FAILED");
1648 // SHOWYXZ("IN TEST PARAM" , testPar);
1649 // SHOWYXZ("OUT TEST CORD" , testCoord);
1650 // if ( !ComputeParameters( testCoord, testPar , ID_BOT_FACE))
1651 // RETURN_BAD_RESULT("TEST ComputeParameters() FAILED");
1652 // SHOWYXZ("OUT TEST PARAM" , testPar);
1657 //================================================================================
1659 * \brief Return pointer to column of nodes
1660 * \param node - bottom node from which the returned column goes up
1661 * \retval const TNodeColumn* - the found column
1663 //================================================================================
1665 const TNodeColumn* StdMeshers_PrismAsBlock::GetNodeColumn(const SMDS_MeshNode* node) const
1667 int sID = node->getshapeId();
1669 map<int, pair< TParam2ColumnMap*, bool > >::const_iterator col_frw =
1670 myShapeIndex2ColumnMap.find( sID );
1671 if ( col_frw != myShapeIndex2ColumnMap.end() ) {
1672 const TParam2ColumnMap* cols = col_frw->second.first;
1673 TParam2ColumnIt u_col = cols->begin();
1674 for ( ; u_col != cols->end(); ++u_col )
1675 if ( u_col->second[ 0 ] == node )
1676 return & u_col->second;
1681 //=======================================================================
1682 //function : GetLayersTransformation
1683 //purpose : Return transformations to get coordinates of nodes of each layer
1684 // by nodes of the bottom. Layer is a set of nodes at a certain step
1685 // from bottom to top.
1686 //=======================================================================
1688 bool StdMeshers_PrismAsBlock::GetLayersTransformation(vector<gp_Trsf> & trsf) const
1690 const int zSize = VerticalSize();
1691 if ( zSize < 3 ) return true;
1692 trsf.resize( zSize - 2 );
1694 // Select some node columns by which we will define coordinate system of layers
1696 vector< const TNodeColumn* > columns;
1698 const TopoDS_Shape& baseFace = Shape(ID_BOT_FACE);
1699 list< TopoDS_Edge > orderedEdges;
1700 list< int > nbEdgesInWires;
1701 GetOrderedEdges( TopoDS::Face( baseFace ), TopoDS_Vertex(), orderedEdges, nbEdgesInWires );
1703 list< TopoDS_Edge >::iterator edgeIt = orderedEdges.begin();
1704 for ( int iE = 0; iE < nbEdgesInWires.front(); ++iE, ++edgeIt )
1706 if ( BRep_Tool::Degenerated( *edgeIt )) continue;
1707 const TParam2ColumnMap* u2colMap =
1708 GetParam2ColumnMap( myHelper->GetMeshDS()->ShapeToIndex( *edgeIt ), isReverse );
1709 if ( !u2colMap ) return false;
1710 isReverse = ( edgeIt->Orientation() == TopAbs_REVERSED );
1711 double f = u2colMap->begin()->first, l = u2colMap->rbegin()->first;
1712 if ( isReverse ) swap ( f, l );
1713 const int nbCol = 5;
1714 for ( int i = 0; i < nbCol; ++i )
1716 double u = f + i/double(nbCol) * ( l - f );
1717 const TNodeColumn* col = & getColumn( u2colMap, u )->second;
1718 if ( columns.empty() || col != columns.back() )
1719 columns.push_back( col );
1724 // Find tolerance to check transformations
1729 for ( int i = 0; i < columns.size(); ++i )
1730 bndBox.Add( gpXYZ( columns[i]->front() ));
1731 tol2 = bndBox.SquareExtent() * 1e-5;
1734 // Compute transformations
1737 gp_Trsf fromCsZ, toCs0;
1738 gp_Ax3 cs0 = getLayerCoordSys(0, columns, xCol );
1739 //double dist0 = cs0.Location().Distance( gpXYZ( (*columns[0])[0]));
1740 toCs0.SetTransformation( cs0 );
1741 for ( int z = 1; z < zSize-1; ++z )
1743 gp_Ax3 csZ = getLayerCoordSys(z, columns, xCol );
1744 //double distZ = csZ.Location().Distance( gpXYZ( (*columns[0])[z]));
1745 fromCsZ.SetTransformation( csZ );
1747 gp_Trsf& t = trsf[ z-1 ];
1748 t = fromCsZ * toCs0;
1749 //t.SetScaleFactor( distZ/dist0 ); - it does not work properly, wrong base point
1751 // check a transformation
1752 for ( int i = 0; i < columns.size(); ++i )
1754 gp_Pnt p0 = gpXYZ( (*columns[i])[0] );
1755 gp_Pnt pz = gpXYZ( (*columns[i])[z] );
1756 t.Transforms( p0.ChangeCoord() );
1757 if ( p0.SquareDistance( pz ) > tol2 )
1764 //================================================================================
1766 * \brief Check curve orientation of a bootom edge
1767 * \param meshDS - mesh DS
1768 * \param columnsMap - node columns map of side face
1769 * \param bottomEdge - the bootom edge
1770 * \param sideFaceID - side face in-block ID
1771 * \retval bool - true if orientation coinside with in-block forward orientation
1773 //================================================================================
1775 bool StdMeshers_PrismAsBlock::IsForwardEdge(SMESHDS_Mesh* meshDS,
1776 const TParam2ColumnMap& columnsMap,
1777 const TopoDS_Edge & bottomEdge,
1778 const int sideFaceID)
1780 bool isForward = false;
1781 if ( SMESH_MesherHelper::IsClosedEdge( bottomEdge ))
1783 isForward = ( bottomEdge.Orientation() == TopAbs_FORWARD );
1787 const TNodeColumn& firstCol = columnsMap.begin()->second;
1788 const SMDS_MeshNode* bottomNode = firstCol[0];
1789 TopoDS_Shape firstVertex = SMESH_MesherHelper::GetSubShapeByNode( bottomNode, meshDS );
1790 isForward = ( firstVertex.IsSame( TopExp::FirstVertex( bottomEdge, true )));
1792 // on 2 of 4 sides first vertex is end
1793 if ( sideFaceID == ID_Fx1z || sideFaceID == ID_F0yz )
1794 isForward = !isForward;
1798 //================================================================================
1800 * \brief Find wall faces by bottom edges
1801 * \param mesh - the mesh
1802 * \param mainShape - the prism
1803 * \param bottomFace - the bottom face
1804 * \param bottomEdges - edges bounding the bottom face
1805 * \param wallFaces - faces list to fill in
1807 //================================================================================
1809 bool StdMeshers_PrismAsBlock::GetWallFaces( SMESH_Mesh* mesh,
1810 const TopoDS_Shape & mainShape,
1811 const TopoDS_Shape & bottomFace,
1812 std::list< TopoDS_Edge >& bottomEdges,
1813 std::list< int > & nbEInW,
1814 std::list< TopoDS_Face >& wallFaces)
1818 TopTools_IndexedMapOfShape faceMap;
1819 TopExp::MapShapes( mainShape, TopAbs_FACE, faceMap );
1821 list< TopoDS_Edge >::iterator edge = bottomEdges.begin();
1822 std::list< int >::iterator nbE = nbEInW.begin();
1824 while ( edge != bottomEdges.end() )
1827 if ( BRep_Tool::Degenerated( *edge ))
1829 edge = bottomEdges.erase( edge );
1835 PShapeIteratorPtr fIt = myHelper->GetAncestors( *edge, *mesh, TopAbs_FACE );
1836 while ( fIt->more() )
1838 const TopoDS_Shape* face = fIt->next();
1839 if ( !bottomFace.IsSame( *face ) && // not bottom
1840 faceMap.FindIndex( *face )) // belongs to the prism
1842 wallFaces.push_back( TopoDS::Face( *face ));
1854 return ( wallFaces.size() == bottomEdges.size() );
1857 //================================================================================
1859 * \brief Constructor
1860 * \param faceID - in-block ID
1861 * \param face - geom face
1862 * \param columnsMap - map of node columns
1863 * \param first - first normalized param
1864 * \param last - last normalized param
1866 //================================================================================
1868 StdMeshers_PrismAsBlock::TSideFace::TSideFace(SMESH_MesherHelper* helper,
1870 const TopoDS_Face& face,
1871 const TopoDS_Edge& baseEdge,
1872 TParam2ColumnMap* columnsMap,
1876 myParamToColumnMap( columnsMap ),
1877 myBaseEdge( baseEdge ),
1880 mySurface.Initialize( face );
1881 myParams.resize( 1 );
1882 myParams[ 0 ] = make_pair( first, last );
1883 myIsForward = StdMeshers_PrismAsBlock::IsForwardEdge( myHelper->GetMeshDS(),
1884 *myParamToColumnMap,
1888 //================================================================================
1890 * \brief Constructor of complex side face
1892 //================================================================================
1894 StdMeshers_PrismAsBlock::TSideFace::
1895 TSideFace(const vector< TSideFace* >& components,
1896 const vector< pair< double, double> > & params)
1897 :myID( components[0] ? components[0]->myID : 0 ),
1898 myParamToColumnMap( 0 ),
1900 myIsForward( true ),
1901 myComponents( components ),
1902 myHelper( components[0] ? components[0]->myHelper : 0 )
1904 //================================================================================
1906 * \brief Copy constructor
1907 * \param other - other side
1909 //================================================================================
1911 StdMeshers_PrismAsBlock::TSideFace::TSideFace( const TSideFace& other )
1914 mySurface = other.mySurface;
1915 myBaseEdge = other.myBaseEdge;
1916 myParams = other.myParams;
1917 myIsForward = other.myIsForward;
1918 myHelper = other.myHelper;
1919 myParamToColumnMap = other.myParamToColumnMap;
1921 myComponents.resize( other.myComponents.size());
1922 for (int i = 0 ; i < myComponents.size(); ++i )
1923 myComponents[ i ] = new TSideFace( *other.myComponents[ i ]);
1926 //================================================================================
1928 * \brief Deletes myComponents
1930 //================================================================================
1932 StdMeshers_PrismAsBlock::TSideFace::~TSideFace()
1934 for (int i = 0 ; i < myComponents.size(); ++i )
1935 if ( myComponents[ i ] )
1936 delete myComponents[ i ];
1939 //================================================================================
1941 * \brief Return geometry of the vertical curve
1942 * \param isMax - true means curve located closer to (1,1,1) block point
1943 * \retval Adaptor3d_Curve* - curve adaptor
1945 //================================================================================
1947 Adaptor3d_Curve* StdMeshers_PrismAsBlock::TSideFace::VertiCurve(const bool isMax) const
1949 if ( !myComponents.empty() ) {
1951 return myComponents.back()->VertiCurve(isMax);
1953 return myComponents.front()->VertiCurve(isMax);
1955 double f = myParams[0].first, l = myParams[0].second;
1956 if ( !myIsForward ) std::swap( f, l );
1957 return new TVerticalEdgeAdaptor( myParamToColumnMap, isMax ? l : f );
1960 //================================================================================
1962 * \brief Return geometry of the top or bottom curve
1964 * \retval Adaptor3d_Curve* -
1966 //================================================================================
1968 Adaptor3d_Curve* StdMeshers_PrismAsBlock::TSideFace::HorizCurve(const bool isTop) const
1970 return new THorizontalEdgeAdaptor( this, isTop );
1973 //================================================================================
1975 * \brief Return pcurves
1976 * \param pcurv - array of 4 pcurves
1977 * \retval bool - is a success
1979 //================================================================================
1981 bool StdMeshers_PrismAsBlock::TSideFace::GetPCurves(Adaptor2d_Curve2d* pcurv[4]) const
1983 int iEdge[ 4 ] = { BOTTOM_EDGE, TOP_EDGE, V0_EDGE, V1_EDGE };
1985 for ( int i = 0 ; i < 4 ; ++i ) {
1986 Handle(Geom2d_Line) line;
1987 switch ( iEdge[ i ] ) {
1989 line = new Geom2d_Line( gp_Pnt2d( 0, 1 ), gp::DX2d() ); break;
1991 line = new Geom2d_Line( gp::Origin2d(), gp::DX2d() ); break;
1993 line = new Geom2d_Line( gp::Origin2d(), gp::DY2d() ); break;
1995 line = new Geom2d_Line( gp_Pnt2d( 1, 0 ), gp::DY2d() ); break;
1997 pcurv[ i ] = new Geom2dAdaptor_Curve( line, 0, 1 );
2002 //================================================================================
2004 * \brief Returns geometry of pcurve on a horizontal face
2005 * \param isTop - is top or bottom face
2006 * \param horFace - a horizontal face
2007 * \retval Adaptor2d_Curve2d* - curve adaptor
2009 //================================================================================
2012 StdMeshers_PrismAsBlock::TSideFace::HorizPCurve(const bool isTop,
2013 const TopoDS_Face& horFace) const
2015 return new TPCurveOnHorFaceAdaptor( this, isTop, horFace );
2018 //================================================================================
2020 * \brief Return a component corresponding to parameter
2021 * \param U - parameter along a horizontal size
2022 * \param localU - parameter along a horizontal size of a component
2023 * \retval TSideFace* - found component
2025 //================================================================================
2027 StdMeshers_PrismAsBlock::TSideFace*
2028 StdMeshers_PrismAsBlock::TSideFace::GetComponent(const double U,double & localU) const
2031 if ( myComponents.empty() )
2032 return const_cast<TSideFace*>( this );
2035 for ( i = 0; i < myComponents.size(); ++i )
2036 if ( U < myParams[ i ].second )
2038 if ( i >= myComponents.size() )
2039 i = myComponents.size() - 1;
2041 double f = myParams[ i ].first, l = myParams[ i ].second;
2042 localU = ( U - f ) / ( l - f );
2043 return myComponents[ i ];
2046 //================================================================================
2048 * \brief Find node columns for a parameter
2049 * \param U - parameter along a horizontal edge
2050 * \param col1 - the 1st found column
2051 * \param col2 - the 2nd found column
2052 * \retval r - normalized position of U between the found columns
2054 //================================================================================
2056 double StdMeshers_PrismAsBlock::TSideFace::GetColumns(const double U,
2057 TParam2ColumnIt & col1,
2058 TParam2ColumnIt & col2) const
2060 double u = U, r = 0;
2061 if ( !myComponents.empty() ) {
2062 TSideFace * comp = GetComponent(U,u);
2063 return comp->GetColumns( u, col1, col2 );
2068 double f = myParams[0].first, l = myParams[0].second;
2069 u = f + u * ( l - f );
2071 col1 = col2 = getColumn( myParamToColumnMap, u );
2072 if ( ++col2 == myParamToColumnMap->end() ) {
2077 double uf = col1->first;
2078 double ul = col2->first;
2079 r = ( u - uf ) / ( ul - uf );
2084 //================================================================================
2086 * \brief Return coordinates by normalized params
2087 * \param U - horizontal param
2088 * \param V - vertical param
2089 * \retval gp_Pnt - result point
2091 //================================================================================
2093 gp_Pnt StdMeshers_PrismAsBlock::TSideFace::Value(const Standard_Real U,
2094 const Standard_Real V) const
2096 if ( !myComponents.empty() ) {
2098 TSideFace * comp = GetComponent(U,u);
2099 return comp->Value( u, V );
2102 TParam2ColumnIt u_col1, u_col2;
2103 double vR, hR = GetColumns( U, u_col1, u_col2 );
2105 const SMDS_MeshNode* n1 = 0;
2106 const SMDS_MeshNode* n2 = 0;
2107 const SMDS_MeshNode* n3 = 0;
2108 const SMDS_MeshNode* n4 = 0;
2110 // BEGIN issue 0020680: EDF 1252 SMESH: Bad cell created by Radial prism in center of torus
2111 // Workaround for a wrongly located point returned by mySurface.Value() for
2112 // UV located near boundary of BSpline surface.
2113 // To bypass the problem, we take point from 3D curve of edge.
2114 // It solves pb of the bloc_fiss_new.py
2115 const double tol = 1e-3;
2116 if ( V < tol || V+tol >= 1. )
2118 n1 = V < tol ? u_col1->second.front() : u_col1->second.back();
2119 n3 = V < tol ? u_col2->second.front() : u_col2->second.back();
2127 TopoDS_Shape s = myHelper->GetSubShapeByNode( n1, myHelper->GetMeshDS() );
2128 if ( s.ShapeType() != TopAbs_EDGE )
2129 s = myHelper->GetSubShapeByNode( n3, myHelper->GetMeshDS() );
2130 if ( s.ShapeType() == TopAbs_EDGE )
2131 edge = TopoDS::Edge( s );
2133 if ( !edge.IsNull() )
2135 double u1 = myHelper->GetNodeU( edge, n1 );
2136 double u3 = myHelper->GetNodeU( edge, n3 );
2137 double u = u1 * ( 1 - hR ) + u3 * hR;
2138 TopLoc_Location loc; double f,l;
2139 Handle(Geom_Curve) curve = BRep_Tool::Curve( edge,loc,f,l );
2140 return curve->Value( u ).Transformed( loc );
2143 // END issue 0020680: EDF 1252 SMESH: Bad cell created by Radial prism in center of torus
2145 vR = getRAndNodes( & u_col1->second, V, n1, n2 );
2146 vR = getRAndNodes( & u_col2->second, V, n3, n4 );
2148 gp_XY uv1 = myHelper->GetNodeUV( mySurface.Face(), n1, n4);
2149 gp_XY uv2 = myHelper->GetNodeUV( mySurface.Face(), n2, n3);
2150 gp_XY uv12 = uv1 * ( 1 - vR ) + uv2 * vR;
2152 gp_XY uv3 = myHelper->GetNodeUV( mySurface.Face(), n3, n2);
2153 gp_XY uv4 = myHelper->GetNodeUV( mySurface.Face(), n4, n1);
2154 gp_XY uv34 = uv3 * ( 1 - vR ) + uv4 * vR;
2156 gp_XY uv = uv12 * ( 1 - hR ) + uv34 * hR;
2158 gp_Pnt p = mySurface.Value( uv.X(), uv.Y() );
2163 //================================================================================
2165 * \brief Return boundary edge
2166 * \param edge - edge index
2167 * \retval TopoDS_Edge - found edge
2169 //================================================================================
2171 TopoDS_Edge StdMeshers_PrismAsBlock::TSideFace::GetEdge(const int iEdge) const
2173 if ( !myComponents.empty() ) {
2175 case V0_EDGE : return myComponents.front()->GetEdge( iEdge );
2176 case V1_EDGE : return myComponents.back() ->GetEdge( iEdge );
2177 default: return TopoDS_Edge();
2181 const SMDS_MeshNode* node = 0;
2182 SMESHDS_Mesh * meshDS = myHelper->GetMesh()->GetMeshDS();
2183 TNodeColumn* column;
2188 column = & (( ++myParamToColumnMap->begin())->second );
2189 node = ( iEdge == TOP_EDGE ) ? column->back() : column->front();
2190 edge = myHelper->GetSubShapeByNode ( node, meshDS );
2191 if ( edge.ShapeType() == TopAbs_VERTEX ) {
2192 column = & ( myParamToColumnMap->begin()->second );
2193 node = ( iEdge == TOP_EDGE ) ? column->back() : column->front();
2198 bool back = ( iEdge == V1_EDGE );
2199 if ( !myIsForward ) back = !back;
2201 column = & ( myParamToColumnMap->rbegin()->second );
2203 column = & ( myParamToColumnMap->begin()->second );
2204 if ( column->size() > 0 )
2205 edge = myHelper->GetSubShapeByNode( (*column)[ 1 ], meshDS );
2206 if ( edge.IsNull() || edge.ShapeType() == TopAbs_VERTEX )
2207 node = column->front();
2212 if ( !edge.IsNull() && edge.ShapeType() == TopAbs_EDGE )
2213 return TopoDS::Edge( edge );
2215 // find edge by 2 vertices
2216 TopoDS_Shape V1 = edge;
2217 TopoDS_Shape V2 = myHelper->GetSubShapeByNode( node, meshDS );
2218 if ( V2.ShapeType() == TopAbs_VERTEX && !V2.IsSame( V1 ))
2220 TopoDS_Shape ancestor = myHelper->GetCommonAncestor( V1, V2, *myHelper->GetMesh(), TopAbs_EDGE);
2221 if ( !ancestor.IsNull() )
2222 return TopoDS::Edge( ancestor );
2224 return TopoDS_Edge();
2227 //================================================================================
2229 * \brief Fill block sub-shapes
2230 * \param shapeMap - map to fill in
2231 * \retval int - nb inserted sub-shapes
2233 //================================================================================
2235 int StdMeshers_PrismAsBlock::TSideFace::InsertSubShapes(TBlockShapes& shapeMap) const
2240 vector< int > edgeIdVec;
2241 SMESH_Block::GetFaceEdgesIDs( myID, edgeIdVec );
2243 for ( int i = BOTTOM_EDGE; i <=V1_EDGE ; ++i ) {
2244 TopoDS_Edge e = GetEdge( i );
2245 if ( !e.IsNull() ) {
2246 nbInserted += SMESH_Block::Insert( e, edgeIdVec[ i ], shapeMap);
2250 // Insert corner vertices
2252 TParam2ColumnIt col1, col2 ;
2253 vector< int > vertIdVec;
2256 SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ V0_EDGE ], vertIdVec);
2257 GetColumns(0, col1, col2 );
2258 const SMDS_MeshNode* node0 = col1->second.front();
2259 const SMDS_MeshNode* node1 = col1->second.back();
2260 TopoDS_Shape v0 = myHelper->GetSubShapeByNode( node0, myHelper->GetMeshDS());
2261 TopoDS_Shape v1 = myHelper->GetSubShapeByNode( node1, myHelper->GetMeshDS());
2262 if ( v0.ShapeType() == TopAbs_VERTEX ) {
2263 nbInserted += SMESH_Block::Insert( v0, vertIdVec[ 0 ], shapeMap);
2265 if ( v1.ShapeType() == TopAbs_VERTEX ) {
2266 nbInserted += SMESH_Block::Insert( v1, vertIdVec[ 1 ], shapeMap);
2270 SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ V1_EDGE ], vertIdVec);
2271 GetColumns(1, col1, col2 );
2272 node0 = col2->second.front();
2273 node1 = col2->second.back();
2274 v0 = myHelper->GetSubShapeByNode( node0, myHelper->GetMeshDS());
2275 v1 = myHelper->GetSubShapeByNode( node1, myHelper->GetMeshDS());
2276 if ( v0.ShapeType() == TopAbs_VERTEX ) {
2277 nbInserted += SMESH_Block::Insert( v0, vertIdVec[ 0 ], shapeMap);
2279 if ( v1.ShapeType() == TopAbs_VERTEX ) {
2280 nbInserted += SMESH_Block::Insert( v1, vertIdVec[ 1 ], shapeMap);
2283 // TopoDS_Vertex V0, V1, Vcom;
2284 // TopExp::Vertices( myBaseEdge, V0, V1, true );
2285 // if ( !myIsForward ) std::swap( V0, V1 );
2287 // // bottom vertex IDs
2288 // SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ _u0 ], vertIdVec);
2289 // SMESH_Block::Insert( V0, vertIdVec[ 0 ], shapeMap);
2290 // SMESH_Block::Insert( V1, vertIdVec[ 1 ], shapeMap);
2292 // TopoDS_Edge sideEdge = GetEdge( V0_EDGE );
2293 // if ( sideEdge.IsNull() || !TopExp::CommonVertex( botEdge, sideEdge, Vcom ))
2296 // // insert one side edge
2298 // if ( Vcom.IsSame( V0 )) edgeID = edgeIdVec[ _v0 ];
2299 // else edgeID = edgeIdVec[ _v1 ];
2300 // SMESH_Block::Insert( sideEdge, edgeID, shapeMap);
2302 // // top vertex of the side edge
2303 // SMESH_Block::GetEdgeVertexIDs( edgeID, vertIdVec);
2304 // TopoDS_Vertex Vtop = TopExp::FirstVertex( sideEdge );
2305 // if ( Vcom.IsSame( Vtop ))
2306 // Vtop = TopExp::LastVertex( sideEdge );
2307 // SMESH_Block::Insert( Vtop, vertIdVec[ 1 ], shapeMap);
2309 // // other side edge
2310 // sideEdge = GetEdge( V1_EDGE );
2311 // if ( sideEdge.IsNull() )
2313 // if ( edgeID = edgeIdVec[ _v1 ]) edgeID = edgeIdVec[ _v0 ];
2314 // else edgeID = edgeIdVec[ _v1 ];
2315 // SMESH_Block::Insert( sideEdge, edgeID, shapeMap);
2318 // TopoDS_Edge topEdge = GetEdge( TOP_EDGE );
2319 // SMESH_Block::Insert( topEdge, edgeIdVec[ _u1 ], shapeMap);
2321 // // top vertex of the other side edge
2322 // if ( !TopExp::CommonVertex( topEdge, sideEdge, Vcom ))
2324 // SMESH_Block::GetEdgeVertexIDs( edgeID, vertIdVec );
2325 // SMESH_Block::Insert( Vcom, vertIdVec[ 1 ], shapeMap);
2330 //================================================================================
2332 * \brief Dump ids of nodes of sides
2334 //================================================================================
2336 void StdMeshers_PrismAsBlock::TSideFace::dumpNodes(int nbNodes) const
2339 cout << endl << "NODES OF FACE "; SMESH_Block::DumpShapeID( myID, cout ) << endl;
2340 THorizontalEdgeAdaptor* hSize0 = (THorizontalEdgeAdaptor*) HorizCurve(0);
2341 cout << "Horiz side 0: "; hSize0->dumpNodes(nbNodes); cout << endl;
2342 THorizontalEdgeAdaptor* hSize1 = (THorizontalEdgeAdaptor*) HorizCurve(1);
2343 cout << "Horiz side 1: "; hSize1->dumpNodes(nbNodes); cout << endl;
2344 TVerticalEdgeAdaptor* vSide0 = (TVerticalEdgeAdaptor*) VertiCurve(0);
2345 cout << "Verti side 0: "; vSide0->dumpNodes(nbNodes); cout << endl;
2346 TVerticalEdgeAdaptor* vSide1 = (TVerticalEdgeAdaptor*) VertiCurve(1);
2347 cout << "Verti side 1: "; vSide1->dumpNodes(nbNodes); cout << endl;
2348 delete hSize0; delete hSize1; delete vSide0; delete vSide1;
2352 //================================================================================
2354 * \brief Creates TVerticalEdgeAdaptor
2355 * \param columnsMap - node column map
2356 * \param parameter - normalized parameter
2358 //================================================================================
2360 StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::
2361 TVerticalEdgeAdaptor( const TParam2ColumnMap* columnsMap, const double parameter)
2363 myNodeColumn = & getColumn( columnsMap, parameter )->second;
2366 //================================================================================
2368 * \brief Return coordinates for the given normalized parameter
2369 * \param U - normalized parameter
2370 * \retval gp_Pnt - coordinates
2372 //================================================================================
2374 gp_Pnt StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::Value(const Standard_Real U) const
2376 const SMDS_MeshNode* n1;
2377 const SMDS_MeshNode* n2;
2378 double r = getRAndNodes( myNodeColumn, U, n1, n2 );
2379 return gpXYZ(n1) * ( 1 - r ) + gpXYZ(n2) * r;
2382 //================================================================================
2384 * \brief Dump ids of nodes
2386 //================================================================================
2388 void StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::dumpNodes(int nbNodes) const
2391 for ( int i = 0; i < nbNodes && i < myNodeColumn->size(); ++i )
2392 cout << (*myNodeColumn)[i]->GetID() << " ";
2393 if ( nbNodes < myNodeColumn->size() )
2394 cout << myNodeColumn->back()->GetID();
2398 //================================================================================
2400 * \brief Return coordinates for the given normalized parameter
2401 * \param U - normalized parameter
2402 * \retval gp_Pnt - coordinates
2404 //================================================================================
2406 gp_Pnt StdMeshers_PrismAsBlock::THorizontalEdgeAdaptor::Value(const Standard_Real U) const
2408 return mySide->TSideFace::Value( U, myV );
2411 //================================================================================
2413 * \brief Dump ids of <nbNodes> first nodes and the last one
2415 //================================================================================
2417 void StdMeshers_PrismAsBlock::THorizontalEdgeAdaptor::dumpNodes(int nbNodes) const
2420 // Not bedugged code. Last node is sometimes incorrect
2421 const TSideFace* side = mySide;
2423 if ( mySide->IsComplex() )
2424 side = mySide->GetComponent(0,u);
2426 TParam2ColumnIt col, col2;
2427 TParam2ColumnMap* u2cols = side->GetColumns();
2428 side->GetColumns( u , col, col2 );
2430 int j, i = myV ? mySide->ColumnHeight()-1 : 0;
2432 const SMDS_MeshNode* n = 0;
2433 const SMDS_MeshNode* lastN
2434 = side->IsForward() ? u2cols->rbegin()->second[ i ] : u2cols->begin()->second[ i ];
2435 for ( j = 0; j < nbNodes && n != lastN; ++j )
2437 n = col->second[ i ];
2438 cout << n->GetID() << " ";
2439 if ( side->IsForward() )
2447 if ( mySide->IsComplex() )
2448 side = mySide->GetComponent(1,u);
2450 side->GetColumns( u , col, col2 );
2451 if ( n != col->second[ i ] )
2452 cout << col->second[ i ]->GetID();
2455 //================================================================================
2457 * \brief Return UV on pcurve for the given normalized parameter
2458 * \param U - normalized parameter
2459 * \retval gp_Pnt - coordinates
2461 //================================================================================
2463 gp_Pnt2d StdMeshers_PrismAsBlock::TPCurveOnHorFaceAdaptor::Value(const Standard_Real U) const
2465 TParam2ColumnIt u_col1, u_col2;
2466 double r = mySide->GetColumns( U, u_col1, u_col2 );
2467 gp_XY uv1 = mySide->GetNodeUV( myFace, u_col1->second[ myZ ]);
2468 gp_XY uv2 = mySide->GetNodeUV( myFace, u_col2->second[ myZ ]);
2469 return uv1 * ( 1 - r ) + uv2 * r;