1 // Copyright (C) 2007-2011 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 // SMESH SMESH : implementaion of SMESH idl descriptions
24 // File : StdMeshers_Prism_3D.cxx
26 // Created : Fri Oct 20 11:37:07 2006
27 // Author : Edward AGAPOV (eap)
29 #include "StdMeshers_Prism_3D.hxx"
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"
38 #include "utilities.h"
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>
46 #include <TopExp_Explorer.hxx>
47 #include <TopTools_ListIteratorOfListOfShape.hxx>
48 #include <TopTools_MapOfShape.hxx>
49 #include <TopTools_SequenceOfShape.hxx>
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) // {\
60 // cout << msg << " ("<< p.X() << "; " <<p.Y() << "; " <<p.Z() << ") " <<endl;\
63 typedef StdMeshers_ProjectionUtils TAssocTool;
64 typedef SMESH_Comment TCom;
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 }; //
73 //================================================================================
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
80 * it returns closest left column
82 //================================================================================
84 TParam2ColumnIt getColumn( const TParam2ColumnMap* columnsMap,
85 const double parameter )
87 TParam2ColumnIt u_col = columnsMap->upper_bound( parameter );
88 if ( u_col != columnsMap->begin() )
90 return u_col; // return left column
93 //================================================================================
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
102 //================================================================================
104 double getRAndNodes( const TNodeColumn* column,
106 const SMDS_MeshNode* & node1,
107 const SMDS_MeshNode* & node2)
109 if ( param >= 1.0 || column->size() == 1) {
110 node1 = node2 = column->back();
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 );
118 node1 = (*column)[ i ];
119 node2 = (*column)[ i + 1];
123 //================================================================================
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
130 //================================================================================
132 void splitParams( const int nbParts,
133 const TParam2ColumnMap* columnsMap,
134 vector< double > & params)
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 )
144 double partSize = ( parLast - par ) / double ( nbParts - i );
145 TParam2ColumnIt par_col = getColumn( columnsMap, par + partSize );
146 if ( par_col->first == par ) {
148 if ( par_col == last_par_col ) {
149 while ( i < nbParts - 1 )
150 params.push_back( par + partSize * i++ );
154 par = par_col->first;
155 params.push_back( par );
157 params.push_back( parLast ); // 1.
160 //================================================================================
162 * \brief Return coordinate system for z-th layer of nodes
164 //================================================================================
166 gp_Ax2 getLayerCoordSys(const int z,
167 const vector< const TNodeColumn* >& columns,
170 // gravity center of a layer
173 for ( int i = 0; i < columns.size(); ++i )
175 O += gpXYZ( (*columns[ i ])[ z ]);
176 if ( vertexCol < 0 &&
177 columns[ i ]->front()->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX )
184 int iPrev = columns.size()-1;
185 for ( int i = 0; i < columns.size(); ++i )
187 gp_Vec v1( O, gpXYZ( (*columns[ iPrev ])[ z ]));
188 gp_Vec v2( O, gpXYZ( (*columns[ i ] )[ z ]));
193 if ( vertexCol >= 0 )
195 O = gpXYZ( (*columns[ vertexCol ])[ z ]);
197 if ( xColumn < 0 || xColumn >= columns.size() )
199 // select a column for X dir
201 for ( int i = 0; i < columns.size(); ++i )
203 double dist = ( O - gpXYZ((*columns[ i ])[ z ])).SquareModulus();
204 if ( dist > maxDist )
213 gp_Vec X( O, gpXYZ( (*columns[ xColumn ])[ z ]));
215 return gp_Ax2( O, Z, X);
218 //================================================================================
220 * \brief Removes submeshes meshed with regular grid from given list
221 * \retval int - nb of removed submeshes
223 //================================================================================
225 int removeQuasiQuads(list< SMESH_subMesh* >& notQuadSubMesh)
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() )
233 SMESH_subMesh* faceSm = *smIt;
234 SMESHDS_SubMesh* faceSmDS = faceSm->GetSubMeshDS();
235 int nbQuads = faceSmDS->NbElements();
236 if ( nbQuads == 0 ) __NEXT_SM;
238 // get oredered edges
239 list< TopoDS_Edge > orderedEdges;
240 list< int > nbEdgesInWires;
242 int nbWires = SMESH_Block::GetOrderedEdges( TopoDS::Face( faceSm->GetSubShape() ),
243 V000, orderedEdges, nbEdgesInWires );
244 if ( nbWires != 1 || nbEdgesInWires.front() <= 4 )
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 )
252 if ( SMESHDS_SubMesh* edgeSmDS = mesh->MeshElements( *edge ))
253 nbSegOnEdge.push_back( edgeSmDS->NbElements() );
255 nbSegOnEdge.push_back(0);
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(); )
263 const TopoDS_Edge& e1 = *edge++;
264 const TopoDS_Edge& e2 = ( edge == orderedEdges.end() ? orderedEdges.front() : *edge );
265 if ( SMESH_Algo::IsContinuous( e1, e2 ))
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 );
273 SMDS_ElemIteratorPtr fIt = vNode->GetInverseElementIterator(SMDSAbs_Face);
274 while ( fIt->more() )
275 nbF += faceSmDS->Contains( fIt->next() );
277 list<int>::iterator nbSegIt1 = nbSegIt++;
278 if ( !vNode || nbF == 2 ) // !vNode - two edges can be meshed as one
281 if ( nbSegIt == nbSegOnEdge.end() ) nbSegIt = nbSegOnEdge.begin();
282 *nbSegIt += *nbSegIt1;
283 nbSegOnEdge.erase( nbSegIt1 );
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
298 smIt = notQuadSubMesh.erase( smIt );
303 return oldNbSM - notQuadSubMesh.size();
307 //=======================================================================
308 //function : StdMeshers_Prism_3D
310 //=======================================================================
312 StdMeshers_Prism_3D::StdMeshers_Prism_3D(int hypId, int studyId, SMESH_Gen* gen)
313 :SMESH_3D_Algo(hypId, studyId, gen)
316 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID); // 1 bit per shape type
317 myProjectTriangles = false;
320 //================================================================================
324 //================================================================================
326 StdMeshers_Prism_3D::~StdMeshers_Prism_3D()
329 //=======================================================================
330 //function : CheckHypothesis
332 //=======================================================================
334 bool StdMeshers_Prism_3D::CheckHypothesis(SMESH_Mesh& aMesh,
335 const TopoDS_Shape& aShape,
336 SMESH_Hypothesis::Hypothesis_Status& aStatus)
338 // Check shape geometry
340 aStatus = SMESH_Hypothesis::HYP_BAD_GEOMETRY;
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() ) {
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");
357 notQuadFaces.push_back( face );
360 if ( !notQuadFaces.empty() )
362 if ( notQuadFaces.size() != 2 )
363 RETURN_BAD_RESULT("Bad nb not quad faces: " << notQuadFaces.size());
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);
372 aStatus = SMESH_Hypothesis::HYP_OK;
376 //=======================================================================
379 //=======================================================================
381 bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theShape)
383 SMESH_MesherHelper helper( theMesh );
386 myHelper->IsQuadraticSubMesh( theShape );
388 // Analyse mesh and geomerty to find block subshapes and submeshes
389 if ( !myBlock.Init( myHelper, theShape ))
390 return error( myBlock.GetError());
392 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
394 int volumeID = meshDS->ShapeToIndex( theShape );
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
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 ]);
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
415 // Create nodes inside the block
417 // try to use transformation (issue 0020680)
418 vector<gp_Trsf> trsf;
419 if ( myBlock.GetLayersTransformation(trsf))
421 // loop on nodes inside the bottom face
422 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.begin();
423 for ( ; bot_column != myBotToColumnMap.end(); ++bot_column )
425 const TNode& tBotNode = bot_column->first; // bottom TNode
426 if ( tBotNode.GetPositionType() != SMDS_TOP_FACE )
427 continue; // node is not inside face
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)
434 const SMDS_MeshNode* & node = *columnNodes;
435 if ( node ) continue; // skip bottom or top node
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 );
442 } // loop on bottom nodes
444 else // use block approach
446 // loop on nodes inside the bottom face
448 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.begin();
449 for ( ; bot_column != myBotToColumnMap.end(); ++bot_column )
451 const TNode& tBotNode = bot_column->first; // bottom TNode
452 if ( tBotNode.GetPositionType() != SMDS_TOP_FACE )
453 continue; // node is not inside face
455 // column nodes; middle part of the column are zero pointers
456 TNodeColumn& column = bot_column->second;
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;
469 myShapeXYZ[ ID_BOT_FACE ] = tBotNode.GetCoords();
470 gp_XYZ botParams = tBotNode.GetParams();
472 // compute top node parameters
473 myShapeXYZ[ ID_TOP_FACE ] = gpXYZ( column.back() );
474 gp_XYZ topParams = botParams;
476 if ( column.size() > 2 ) {
477 gp_Pnt topCoords = myShapeXYZ[ ID_TOP_FACE ];
478 if ( !myBlock.ComputeParameters( topCoords, topParams, ID_TOP_FACE, topParams ))
479 return error(TCom("Can't compute normalized parameters ")
480 << "for node " << column.back()->GetID()
481 << " on the face #"<< column.back()->getshapeId() );
485 TNodeColumn::iterator columnNodes = column.begin();
486 for ( int z = 0; columnNodes != column.end(); ++columnNodes, ++z)
488 const SMDS_MeshNode* & node = *columnNodes;
489 if ( node ) continue; // skip bottom or top node
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;
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 ))
505 // compute coords for a new node
507 if ( !SMESH_Block::ShellPoint( params, myShapeXYZ, coords ))
508 return error("Can't compute coordinates by normalized parameters");
510 SHOWYXZ("TOPFacePoint ",myShapeXYZ[ ID_TOP_FACE]);
511 SHOWYXZ("BOT Node "<< tBotNode.myNode->GetID(),gpXYZ(tBotNode.myNode));
512 SHOWYXZ("ShellPoint ",coords);
515 node = meshDS->AddNode( coords.X(), coords.Y(), coords.Z() );
516 meshDS->SetNodeInVolume( node, volumeID );
518 } // loop on bottom nodes
523 SMESHDS_SubMesh* smDS = myBlock.SubMeshDS( ID_BOT_FACE );
524 if ( !smDS ) return error(COMPERR_BAD_INPUT_MESH, "Null submesh");
526 // loop on bottom mesh faces
527 SMDS_ElemIteratorPtr faceIt = smDS->GetElements();
528 while ( faceIt->more() )
530 const SMDS_MeshElement* face = faceIt->next();
531 if ( !face || face->GetType() != SMDSAbs_Face )
533 int nbNodes = face->NbNodes();
534 if ( face->IsQuadratic() )
537 // find node columns for each node
538 vector< const TNodeColumn* > columns( nbNodes );
539 for ( int i = 0; i < nbNodes; ++i )
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;
549 columns[ i ] = myBlock.GetNodeColumn( n );
551 return error(TCom("No side nodes found above node ") << n->GetID() );
555 AddPrisms( columns, myHelper );
557 } // loop on bottom mesh faces
560 myBotToColumnMap.clear();
567 //=======================================================================
568 //function : Evaluate
570 //=======================================================================
572 bool StdMeshers_Prism_3D::Evaluate(SMESH_Mesh& theMesh,
573 const TopoDS_Shape& theShape,
574 MapShapeNbElems& aResMap)
576 // find face contains only triangles
577 vector < SMESH_subMesh * >meshFaces;
578 TopTools_SequenceOfShape aFaces;
579 int NumBase = 0, i = 0, NbQFs = 0;
580 for (TopExp_Explorer exp(theShape, TopAbs_FACE); exp.More(); exp.Next()) {
582 aFaces.Append(exp.Current());
583 SMESH_subMesh *aSubMesh = theMesh.GetSubMesh(exp.Current());
584 meshFaces.push_back(aSubMesh);
585 MapShapeNbElemsItr anIt = aResMap.find(meshFaces[i-1]);
586 if( anIt==aResMap.end() ) {
587 SMESH_ComputeErrorPtr& smError = aSubMesh->GetComputeError();
588 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
591 std::vector<int> aVec = (*anIt).second;
592 int nbtri = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
593 int nbqua = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
594 if( nbtri==0 && nbqua>0 ) {
603 std::vector<int> aResVec(SMDSEntity_Last);
604 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
605 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
606 aResMap.insert(std::make_pair(sm,aResVec));
607 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
608 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
612 if(NumBase==0) NumBase = 1; // only quads => set 1 faces as base
614 // find number of 1d elems for base face
616 TopTools_MapOfShape Edges1;
617 for (TopExp_Explorer exp(aFaces.Value(NumBase), TopAbs_EDGE); exp.More(); exp.Next()) {
618 Edges1.Add(exp.Current());
619 SMESH_subMesh *sm = theMesh.GetSubMesh(exp.Current());
621 MapShapeNbElemsItr anIt = aResMap.find(sm);
622 if( anIt == aResMap.end() ) continue;
623 std::vector<int> aVec = (*anIt).second;
624 nb1d += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
627 // find face opposite to base face
629 for(i=1; i<=6; i++) {
630 if(i==NumBase) continue;
631 bool IsOpposite = true;
632 for(TopExp_Explorer exp(aFaces.Value(i), TopAbs_EDGE); exp.More(); exp.Next()) {
633 if( Edges1.Contains(exp.Current()) ) {
643 // find number of 2d elems on side faces
645 for(i=1; i<=6; i++) {
646 if( i==OppNum || i==NumBase ) continue;
647 MapShapeNbElemsItr anIt = aResMap.find( meshFaces[i-1] );
648 if( anIt == aResMap.end() ) continue;
649 std::vector<int> aVec = (*anIt).second;
650 nb2d += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
653 MapShapeNbElemsItr anIt = aResMap.find( meshFaces[NumBase-1] );
654 std::vector<int> aVec = (*anIt).second;
655 bool IsQuadratic = (aVec[SMDSEntity_Quad_Triangle]>aVec[SMDSEntity_Triangle]) ||
656 (aVec[SMDSEntity_Quad_Quadrangle]>aVec[SMDSEntity_Quadrangle]);
657 int nb2d_face0_3 = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
658 int nb2d_face0_4 = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
659 int nb0d_face0 = aVec[SMDSEntity_Node];
660 int nb1d_face0_int = ( nb2d_face0_3*3 + nb2d_face0_4*4 - nb1d ) / 2;
662 std::vector<int> aResVec(SMDSEntity_Last);
663 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
665 aResVec[SMDSEntity_Quad_Penta] = nb2d_face0_3 * ( nb2d/nb1d );
666 aResVec[SMDSEntity_Quad_Hexa] = nb2d_face0_4 * ( nb2d/nb1d );
667 aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
670 aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
671 aResVec[SMDSEntity_Penta] = nb2d_face0_3 * ( nb2d/nb1d );
672 aResVec[SMDSEntity_Hexa] = nb2d_face0_4 * ( nb2d/nb1d );
674 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
675 aResMap.insert(std::make_pair(sm,aResVec));
681 //================================================================================
683 * \brief Create prisms
684 * \param columns - columns of nodes generated from nodes of a mesh face
685 * \param helper - helper initialized by mesh and shape to add prisms to
687 //================================================================================
689 void StdMeshers_Prism_3D::AddPrisms( vector<const TNodeColumn*> & columns,
690 SMESH_MesherHelper* helper)
692 SMESHDS_Mesh * meshDS = helper->GetMeshDS();
693 int shapeID = helper->GetSubShapeID();
695 int nbNodes = columns.size();
696 int nbZ = columns[0]->size();
697 if ( nbZ < 2 ) return;
699 // find out orientation
700 bool isForward = true;
701 SMDS_VolumeTool vTool;
705 const SMDS_MeshNode* botNodes[3] = { (*columns[0])[z-1],
707 (*columns[2])[z-1] };
708 const SMDS_MeshNode* topNodes[3] = { (*columns[0])[z],
711 SMDS_VolumeOfNodes tmpVol ( botNodes[0], botNodes[1], botNodes[2],
712 topNodes[0], topNodes[1], topNodes[2]);
713 vTool.Set( &tmpVol );
714 isForward = vTool.IsForward();
718 const SMDS_MeshNode* botNodes[4] = { (*columns[0])[z-1], (*columns[1])[z-1],
719 (*columns[2])[z-1], (*columns[3])[z-1] };
720 const SMDS_MeshNode* topNodes[4] = { (*columns[0])[z], (*columns[1])[z],
721 (*columns[2])[z], (*columns[3])[z] };
722 SMDS_VolumeOfNodes tmpVol ( botNodes[0], botNodes[1], botNodes[2], botNodes[3],
723 topNodes[0], topNodes[1], topNodes[2], topNodes[3]);
724 vTool.Set( &tmpVol );
725 isForward = vTool.IsForward();
730 // vertical loop on columns
731 for ( z = 1; z < nbZ; ++z )
733 SMDS_MeshElement* vol = 0;
737 const SMDS_MeshNode* botNodes[3] = { (*columns[0])[z-1],
739 (*columns[2])[z-1] };
740 const SMDS_MeshNode* topNodes[3] = { (*columns[0])[z],
744 vol = helper->AddVolume( botNodes[0], botNodes[1], botNodes[2],
745 topNodes[0], topNodes[1], topNodes[2]);
747 vol = helper->AddVolume( topNodes[0], topNodes[1], topNodes[2],
748 botNodes[0], botNodes[1], botNodes[2]);
752 const SMDS_MeshNode* botNodes[4] = { (*columns[0])[z-1], (*columns[1])[z-1],
753 (*columns[2])[z-1], (*columns[3])[z-1] };
754 const SMDS_MeshNode* topNodes[4] = { (*columns[0])[z], (*columns[1])[z],
755 (*columns[2])[z], (*columns[3])[z] };
757 vol = helper->AddVolume( botNodes[0], botNodes[1], botNodes[2], botNodes[3],
758 topNodes[0], topNodes[1], topNodes[2], topNodes[3]);
760 vol = helper->AddVolume( topNodes[0], topNodes[1], topNodes[2], topNodes[3],
761 botNodes[0], botNodes[1], botNodes[2], botNodes[3]);
766 vector<const SMDS_MeshNode*> nodes( 2*nbNodes + 4*nbNodes);
767 vector<int> quantities( 2 + nbNodes, 4 );
768 quantities[0] = quantities[1] = nbNodes;
769 columns.resize( nbNodes + 1 );
770 columns[ nbNodes ] = columns[ 0 ];
771 for ( int i = 0; i < nbNodes; ++i ) {
772 nodes[ i ] = (*columns[ i ])[z-1]; // bottom
773 nodes[ i+nbNodes ] = (*columns[ i ])[z ]; // top
775 int di = 2*nbNodes + 4*i - 1;
776 nodes[ di ] = (*columns[i ])[z-1];
777 nodes[ di+1 ] = (*columns[i+1])[z-1];
778 nodes[ di+2 ] = (*columns[i+1])[z ];
779 nodes[ di+3 ] = (*columns[i ])[z ];
781 vol = meshDS->AddPolyhedralVolume( nodes, quantities );
783 if ( vol && shapeID > 0 )
784 meshDS->SetMeshElementOnShape( vol, shapeID );
788 //================================================================================
790 * \brief Find correspondence between bottom and top nodes
791 * If elements on the bottom and top faces are topologically different,
792 * and projection is possible and allowed, perform the projection
793 * \retval bool - is a success or not
795 //================================================================================
797 bool StdMeshers_Prism_3D::assocOrProjBottom2Top()
799 SMESH_subMesh * botSM = myBlock.SubMesh( ID_BOT_FACE );
800 SMESH_subMesh * topSM = myBlock.SubMesh( ID_TOP_FACE );
802 SMESHDS_SubMesh * botSMDS = botSM->GetSubMeshDS();
803 SMESHDS_SubMesh * topSMDS = topSM->GetSubMeshDS();
805 if ( !botSMDS || botSMDS->NbElements() == 0 )
806 return error(TCom("No elememts on face #") << botSM->GetId());
808 bool needProject = false;
810 botSMDS->NbElements() != topSMDS->NbElements() ||
811 botSMDS->NbNodes() != topSMDS->NbNodes())
813 MESSAGE("nb elem bot " << botSMDS->NbElements() << " top " << topSMDS->NbElements());
814 MESSAGE("nb node bot " << botSMDS->NbNodes() << " top " << topSMDS->NbNodes());
815 if ( myBlock.HasNotQuadElemOnTop() )
816 return error(TCom("Mesh on faces #") << botSM->GetId()
817 <<" and #"<< topSM->GetId() << " seems different" );
821 if ( 0/*needProject && !myProjectTriangles*/ )
822 return error(TCom("Mesh on faces #") << botSM->GetId()
823 <<" and #"<< topSM->GetId() << " seems different" );
824 ///RETURN_BAD_RESULT("Need to project but not allowed");
828 return projectBottomToTop();
831 TopoDS_Face botFace = TopoDS::Face( myBlock.Shape( ID_BOT_FACE ));
832 TopoDS_Face topFace = TopoDS::Face( myBlock.Shape( ID_TOP_FACE ));
833 // associate top and bottom faces
834 TAssocTool::TShapeShapeMap shape2ShapeMap;
835 if ( !TAssocTool::FindSubShapeAssociation( botFace, myBlock.Mesh(),
836 topFace, myBlock.Mesh(),
838 return error(TCom("Topology of faces #") << botSM->GetId()
839 <<" and #"<< topSM->GetId() << " seems different" );
841 // Find matching nodes of top and bottom faces
843 if ( ! TAssocTool::FindMatchingNodesOnFaces( botFace, myBlock.Mesh(),
844 topFace, myBlock.Mesh(),
845 shape2ShapeMap, n2nMap ))
846 return error(TCom("Mesh on faces #") << botSM->GetId()
847 <<" and #"<< topSM->GetId() << " seems different" );
849 // Fill myBotToColumnMap
851 int zSize = myBlock.VerticalSize();
853 TNodeNodeMap::iterator bN_tN = n2nMap.begin();
854 for ( ; bN_tN != n2nMap.end(); ++bN_tN )
856 const SMDS_MeshNode* botNode = bN_tN->first;
857 const SMDS_MeshNode* topNode = bN_tN->second;
858 if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
859 continue; // wall columns are contained in myBlock
860 // create node column
862 TNode2ColumnMap::iterator bN_col =
863 myBotToColumnMap.insert( make_pair ( bN, TNodeColumn() )).first;
864 TNodeColumn & column = bN_col->second;
865 column.resize( zSize );
866 column.front() = botNode;
867 column.back() = topNode;
872 //================================================================================
874 * \brief Remove quadrangles from the top face and
875 * create triangles there by projection from the bottom
876 * \retval bool - a success or not
878 //================================================================================
880 bool StdMeshers_Prism_3D::projectBottomToTop()
882 SMESH_subMesh * botSM = myBlock.SubMesh( ID_BOT_FACE );
883 SMESH_subMesh * topSM = myBlock.SubMesh( ID_TOP_FACE );
885 SMESHDS_SubMesh * botSMDS = botSM->GetSubMeshDS();
886 SMESHDS_SubMesh * topSMDS = topSM->GetSubMeshDS();
889 topSM->ComputeStateEngine( SMESH_subMesh::CLEAN );
891 SMESHDS_Mesh* meshDS = myBlock.MeshDS();
892 int shapeID = myHelper->GetSubShapeID();
893 int topFaceID = meshDS->ShapeToIndex( topSM->GetSubShape() );
895 // Fill myBotToColumnMap
897 int zSize = myBlock.VerticalSize();
899 SMDS_NodeIteratorPtr nIt = botSMDS->GetNodes();
900 while ( nIt->more() )
902 const SMDS_MeshNode* botNode = nIt->next();
903 if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
905 // compute bottom node params
907 gp_XYZ paramHint(-1,-1,-1);
908 if ( prevTNode.IsNeighbor( bN ))
909 paramHint = prevTNode.GetParams();
910 if ( !myBlock.ComputeParameters( bN.GetCoords(), bN.ChangeParams(),
911 ID_BOT_FACE, paramHint ))
912 return error(TCom("Can't compute normalized parameters for node ")
913 << botNode->GetID() << " on the face #"<< botSM->GetId() );
915 // compute top node coords
916 gp_XYZ topXYZ; gp_XY topUV;
917 if ( !myBlock.FacePoint( ID_TOP_FACE, bN.GetParams(), topXYZ ) ||
918 !myBlock.FaceUV ( ID_TOP_FACE, bN.GetParams(), topUV ))
919 return error(TCom("Can't compute coordinates "
920 "by normalized parameters on the face #")<< topSM->GetId() );
921 SMDS_MeshNode * topNode = meshDS->AddNode( topXYZ.X(),topXYZ.Y(),topXYZ.Z() );
922 meshDS->SetNodeOnFace( topNode, topFaceID, topUV.X(), topUV.Y() );
923 // create node column
924 TNode2ColumnMap::iterator bN_col =
925 myBotToColumnMap.insert( make_pair ( bN, TNodeColumn() )).first;
926 TNodeColumn & column = bN_col->second;
927 column.resize( zSize );
928 column.front() = botNode;
929 column.back() = topNode;
934 // loop on bottom mesh faces
935 SMDS_ElemIteratorPtr faceIt = botSMDS->GetElements();
936 while ( faceIt->more() )
938 const SMDS_MeshElement* face = faceIt->next();
939 if ( !face || face->GetType() != SMDSAbs_Face )
941 int nbNodes = face->NbNodes();
942 if ( face->IsQuadratic() )
945 // find top node in columns for each bottom node
946 vector< const SMDS_MeshNode* > nodes( nbNodes );
947 for ( int i = 0; i < nbNodes; ++i )
949 const SMDS_MeshNode* n = face->GetNode( nbNodes - i - 1 );
950 if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ) {
951 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.find( n );
952 if ( bot_column == myBotToColumnMap.end() )
953 return error(TCom("No nodes found above node ") << n->GetID() );
954 nodes[ i ] = bot_column->second.back();
957 const TNodeColumn* column = myBlock.GetNodeColumn( n );
959 return error(TCom("No side nodes found above node ") << n->GetID() );
960 nodes[ i ] = column->back();
963 // create a face, with reversed orientation
964 SMDS_MeshElement* newFace = 0;
968 newFace = myHelper->AddFace(nodes[0], nodes[1], nodes[2]);
972 newFace = myHelper->AddFace( nodes[0], nodes[1], nodes[2], nodes[3] );
976 newFace = meshDS->AddPolygonalFace( nodes );
978 if ( newFace && shapeID > 0 )
979 meshDS->SetMeshElementOnShape( newFace, shapeID );
985 //================================================================================
987 * \brief Set projection coordinates of a node to a face and it's subshapes
988 * \param faceID - the face given by in-block ID
989 * \param params - node normalized parameters
990 * \retval bool - is a success
992 //================================================================================
994 bool StdMeshers_Prism_3D::setFaceAndEdgesXYZ( const int faceID, const gp_XYZ& params, int z )
996 // find base and top edges of the face
997 enum { BASE = 0, TOP, LEFT, RIGHT };
998 vector< int > edgeVec; // 0-base, 1-top
999 SMESH_Block::GetFaceEdgesIDs( faceID, edgeVec );
1001 myBlock.EdgePoint( edgeVec[ BASE ], params, myShapeXYZ[ edgeVec[ BASE ]]);
1002 myBlock.EdgePoint( edgeVec[ TOP ], params, myShapeXYZ[ edgeVec[ TOP ]]);
1004 SHOWYXZ("\nparams ", params);
1005 SHOWYXZ("TOP is " <<edgeVec[ TOP ], myShapeXYZ[ edgeVec[ TOP]]);
1006 SHOWYXZ("BASE is "<<edgeVec[ BASE], myShapeXYZ[ edgeVec[ BASE]]);
1008 if ( faceID == SMESH_Block::ID_Fx0z || faceID == SMESH_Block::ID_Fx1z )
1010 myBlock.EdgePoint( edgeVec[ LEFT ], params, myShapeXYZ[ edgeVec[ LEFT ]]);
1011 myBlock.EdgePoint( edgeVec[ RIGHT ], params, myShapeXYZ[ edgeVec[ RIGHT ]]);
1013 SHOWYXZ("VER "<<edgeVec[ LEFT], myShapeXYZ[ edgeVec[ LEFT]]);
1014 SHOWYXZ("VER "<<edgeVec[ RIGHT], myShapeXYZ[ edgeVec[ RIGHT]]);
1016 myBlock.FacePoint( faceID, params, myShapeXYZ[ faceID ]);
1017 SHOWYXZ("FacePoint "<<faceID, myShapeXYZ[ faceID]);
1022 //================================================================================
1024 * \brief Return true if this node and other one belong to one face
1026 //================================================================================
1028 bool TNode::IsNeighbor( const TNode& other ) const
1030 if ( !other.myNode || !myNode ) return false;
1032 SMDS_ElemIteratorPtr fIt = other.myNode->GetInverseElementIterator(SMDSAbs_Face);
1033 while ( fIt->more() )
1034 if ( fIt->next()->GetNodeIndex( myNode ) >= 0 )
1039 //================================================================================
1041 * \brief Constructor. Initialization is needed
1043 //================================================================================
1045 StdMeshers_PrismAsBlock::StdMeshers_PrismAsBlock()
1050 StdMeshers_PrismAsBlock::~StdMeshers_PrismAsBlock()
1054 void StdMeshers_PrismAsBlock::Clear()
1057 myShapeIDMap.Clear();
1061 delete mySide; mySide = 0;
1063 myParam2ColumnMaps.clear();
1064 myShapeIndex2ColumnMap.clear();
1067 //================================================================================
1069 * \brief Initialization.
1070 * \param helper - helper loaded with mesh and 3D shape
1071 * \param shape3D - a closed shell or solid
1072 * \retval bool - false if a mesh or a shape are KO
1074 //================================================================================
1076 bool StdMeshers_PrismAsBlock::Init(SMESH_MesherHelper* helper,
1077 const TopoDS_Shape& shape3D)
1080 delete mySide; mySide = 0;
1082 vector< TSideFace* > sideFaces( NB_WALL_FACES, 0 );
1083 vector< pair< double, double> > params ( NB_WALL_FACES );
1084 mySide = new TSideFace( sideFaces, params );
1087 SMESHDS_Mesh* meshDS = myHelper->GetMeshDS();
1089 SMESH_Block::init();
1090 myShapeIDMap.Clear();
1091 myShapeIndex2ColumnMap.clear();
1093 int wallFaceIds[ NB_WALL_FACES ] = { // to walk around a block
1094 SMESH_Block::ID_Fx0z, SMESH_Block::ID_F1yz,
1095 SMESH_Block::ID_Fx1z, SMESH_Block::ID_F0yz
1098 myError = SMESH_ComputeError::New();
1100 // -------------------------------------------------------------
1101 // Look for top and bottom faces: not quadrangle ones or meshed
1102 // with not quadrangle elements
1103 // -------------------------------------------------------------
1105 list< SMESH_subMesh* > notQuadGeomSubMesh;
1106 list< SMESH_subMesh* > notQuadElemSubMesh;
1109 SMESH_subMesh* mainSubMesh = myHelper->GetMesh()->GetSubMeshContaining( shape3D );
1110 if ( !mainSubMesh ) return error(COMPERR_BAD_INPUT_MESH,"Null submesh of shape3D");
1112 // analyse face submeshes
1113 SMESH_subMeshIteratorPtr smIt = mainSubMesh->getDependsOnIterator(false,false);
1114 while ( smIt->more() )
1116 SMESH_subMesh* sm = smIt->next();
1117 const TopoDS_Shape& face = sm->GetSubShape();
1118 if ( face.ShapeType() != TopAbs_FACE )
1122 // is quadrangle face?
1123 list< TopoDS_Edge > orderedEdges;
1124 list< int > nbEdgesInWires;
1126 int nbWires = GetOrderedEdges( TopoDS::Face( face ),
1127 V000, orderedEdges, nbEdgesInWires );
1128 if ( nbWires != 1 || nbEdgesInWires.front() != 4 )
1129 notQuadGeomSubMesh.push_back( sm );
1131 // look for not quadrangle mesh elements
1132 if ( SMESHDS_SubMesh* smDS = sm->GetSubMeshDS() ) {
1133 bool hasNotQuad = false;
1134 SMDS_ElemIteratorPtr eIt = smDS->GetElements();
1135 while ( eIt->more() && !hasNotQuad ) {
1136 const SMDS_MeshElement* elem = eIt->next();
1137 if ( elem->GetType() == SMDSAbs_Face ) {
1138 int nbNodes = elem->NbNodes();
1139 if ( elem->IsQuadratic() )
1141 hasNotQuad = ( nbNodes != 4 );
1145 notQuadElemSubMesh.push_back( sm );
1148 return error(COMPERR_BAD_INPUT_MESH,TCom("Not meshed face #")<<sm->GetId());
1150 // check if a quadrangle face is meshed with a quadranglar grid
1151 if ( notQuadGeomSubMesh.back() != sm &&
1152 notQuadElemSubMesh.back() != sm )
1154 // count nb edges on face sides
1155 vector< int > nbEdges;
1156 nbEdges.reserve( nbEdgesInWires.front() );
1157 for ( list< TopoDS_Edge >::iterator edge = orderedEdges.begin();
1158 edge != orderedEdges.end(); ++edge )
1160 if ( SMESHDS_SubMesh* smDS = meshDS->MeshElements( *edge ))
1161 nbEdges.push_back ( smDS->NbElements() );
1163 nbEdges.push_back ( 0 );
1165 int nbQuads = sm->GetSubMeshDS()->NbElements();
1166 if ( nbEdges[0] * nbEdges[1] != nbQuads ||
1167 nbEdges[0] != nbEdges[2] ||
1168 nbEdges[1] != nbEdges[3] )
1169 notQuadElemSubMesh.push_back( sm );
1173 // ----------------------------------------------------------------------
1174 // Analyse mesh and topology of faces: choose the bottom submesh.
1175 // If there are not quadrangle geom faces, they are top and bottom ones.
1176 // Not quadrangle geom faces must be only on top and bottom.
1177 // ----------------------------------------------------------------------
1179 SMESH_subMesh * botSM = 0;
1180 SMESH_subMesh * topSM = 0;
1182 int nbNotQuad = notQuadGeomSubMesh.size();
1183 int nbNotQuadMeshed = notQuadElemSubMesh.size();
1184 bool hasNotQuad = ( nbNotQuad || nbNotQuadMeshed );
1187 if ( nbNotQuadMeshed > 2 )
1189 return error(COMPERR_BAD_INPUT_MESH,
1190 TCom("More than 2 faces with not quadrangle elements: ")
1193 int nbQuasiQuads = 0;
1194 if ( nbNotQuad > 0 && nbNotQuad != 2 )
1196 // Issue 0020843 - one of side faces is quasi-quadrilateral.
1197 // Remove from notQuadGeomSubMesh faces meshed with regular grid
1198 nbQuasiQuads = removeQuasiQuads( notQuadGeomSubMesh );
1199 nbNotQuad -= nbQuasiQuads;
1200 if ( nbNotQuad > 0 && nbNotQuad != 2 )
1201 return error(COMPERR_BAD_SHAPE,
1202 TCom("More than 2 not quadrilateral faces: ")
1206 // get found submeshes
1209 if ( nbNotQuadMeshed > 0 ) botSM = notQuadElemSubMesh.front();
1210 else botSM = notQuadGeomSubMesh.front();
1211 if ( nbNotQuadMeshed > 1 ) topSM = notQuadElemSubMesh.back();
1212 else if ( nbNotQuad > 1 ) topSM = notQuadGeomSubMesh.back();
1214 // detect other bad cases
1215 if ( nbNotQuad == 2 && nbNotQuadMeshed > 0 ) {
1217 if ( nbNotQuadMeshed == 1 )
1218 ok = ( find( notQuadGeomSubMesh.begin(),
1219 notQuadGeomSubMesh.end(), botSM ) != notQuadGeomSubMesh.end() );
1221 ok = ( notQuadGeomSubMesh == notQuadElemSubMesh );
1223 return error(COMPERR_BAD_INPUT_MESH, "Side face meshed with not quadrangle elements");
1226 myNotQuadOnTop = ( nbNotQuadMeshed > 1 );
1227 MESSAGE("myNotQuadOnTop " << myNotQuadOnTop << " nbNotQuadMeshed " << nbNotQuadMeshed);
1229 // ----------------------------------------------------------
1231 if ( nbNotQuad == 0 ) // Standard block of 6 quadrangle faces ?
1233 // SMESH_Block will perform geometry analysis, we need just to find 2
1234 // connected vertices on top and bottom
1236 TopoDS_Vertex Vbot, Vtop;
1237 if ( nbNotQuadMeshed > 0 ) // Look for vertices
1239 TopTools_IndexedMapOfShape edgeMap;
1240 TopExp::MapShapes( botSM->GetSubShape(), TopAbs_EDGE, edgeMap );
1241 // vertex 1 is any vertex of the bottom face
1242 Vbot = TopExp::FirstVertex( TopoDS::Edge( edgeMap( 1 )));
1243 // vertex 2 is end vertex of edge sharing Vbot and not belonging to the bottom face
1244 TopTools_ListIteratorOfListOfShape ancestIt = Mesh()->GetAncestors( Vbot );
1245 for ( ; Vtop.IsNull() && ancestIt.More(); ancestIt.Next() )
1247 const TopoDS_Shape & ancestor = ancestIt.Value();
1248 if ( ancestor.ShapeType() == TopAbs_EDGE && !edgeMap.FindIndex( ancestor ))
1250 TopoDS_Vertex V1, V2;
1251 TopExp::Vertices( TopoDS::Edge( ancestor ), V1, V2);
1252 if ( Vbot.IsSame ( V1 )) Vtop = V2;
1253 else if ( Vbot.IsSame ( V2 )) Vtop = V1;
1254 // check that Vtop belongs to shape3D
1255 TopExp_Explorer exp( shape3D, TopAbs_VERTEX );
1256 for ( ; exp.More(); exp.Next() )
1257 if ( Vtop.IsSame( exp.Current() ))
1264 // get shell from shape3D
1266 TopExp_Explorer exp( shape3D, TopAbs_SHELL );
1268 for ( ; exp.More(); exp.Next(), ++nbShell )
1269 shell = TopoDS::Shell( exp.Current() );
1270 // if ( nbShell != 1 )
1271 // RETURN_BAD_RESULT("There must be 1 shell in the block");
1273 // Load geometry in SMESH_Block
1274 if ( !SMESH_Block::FindBlockShapes( shell, Vbot, Vtop, myShapeIDMap )) {
1276 return error(COMPERR_BAD_SHAPE, "Can't detect top and bottom of a prism");
1279 if ( !botSM ) botSM = Mesh()->GetSubMeshContaining( myShapeIDMap( ID_BOT_FACE ));
1280 if ( !topSM ) topSM = Mesh()->GetSubMeshContaining( myShapeIDMap( ID_TOP_FACE ));
1283 } // end Standard block of 6 quadrangle faces
1284 // --------------------------------------------------------
1286 // Here the top and bottom faces are found
1287 if ( nbNotQuadMeshed == 2 ) // roughly check correspondence of horiz meshes
1289 // SMESHDS_SubMesh* topSMDS = topSM->GetSubMeshDS();
1290 // SMESHDS_SubMesh* botSMDS = botSM->GetSubMeshDS();
1291 // if ( topSMDS->NbNodes() != botSMDS->NbNodes() ||
1292 // topSMDS->NbElements() != botSMDS->NbElements() )
1293 // RETURN_BAD_RESULT("Top mesh doesn't correspond to bottom one");
1296 // ---------------------------------------------------------
1297 // If there are not quadrangle geom faces, we emulate
1298 // a block of 6 quadrangle faces.
1299 // Load SMESH_Block with faces and edges geometry
1300 // ---------------------------------------------------------
1303 // find vertex 000 - the one with smallest coordinates (for easy DEBUG :-)
1305 double minVal = DBL_MAX, minX, val;
1306 for ( TopExp_Explorer exp( botSM->GetSubShape(), TopAbs_VERTEX );
1307 exp.More(); exp.Next() )
1309 const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
1310 gp_Pnt P = BRep_Tool::Pnt( v );
1311 val = P.X() + P.Y() + P.Z();
1312 if ( val < minVal || ( val == minVal && P.X() < minX )) {
1319 // Get ordered bottom edges
1320 list< TopoDS_Edge > orderedEdges;
1322 SMESH_Block::GetOrderedEdges( TopoDS::Face( botSM->GetSubShape().Reversed() ),
1323 V000, orderedEdges, nbEInW );
1324 // if ( nbEInW.size() != 1 )
1325 // RETURN_BAD_RESULT("Wrong prism geometry");
1327 // Get Wall faces corresponding to the ordered bottom edges
1328 list< TopoDS_Face > wallFaces;
1329 if ( !GetWallFaces( Mesh(), shape3D, botSM->GetSubShape(), orderedEdges, nbEInW, wallFaces))
1330 return error(COMPERR_BAD_SHAPE, "Can't find side faces");
1332 // check that the found top and bottom faces are opposite
1334 for (TopExp_Explorer edge(botSM->GetSubShape(), TopAbs_EDGE); edge.More(); edge.Next())
1335 if ( helper->IsSubShape( edge.Current(), topSM->GetSubShape() ))
1336 return error(notQuadGeomSubMesh.empty() ? COMPERR_BAD_INPUT_MESH : COMPERR_BAD_SHAPE,
1337 "Non-quadrilateral faces are not opposite");
1340 // Protect from a distorted block (test 3D_mesh_HEXA3D/B7 on 32bit platform)
1341 // check that all wall faces have an edge common with the top face
1343 list< TopoDS_Face >::iterator faceIt = wallFaces.begin();
1344 for ( ; faceIt != wallFaces.end(); ++faceIt )
1346 bool hasCommon = false;
1347 for (TopExp_Explorer edge(*faceIt, TopAbs_EDGE); !hasCommon && edge.More(); edge.Next())
1348 if ( helper->IsSubShape( edge.Current(), topSM->GetSubShape() ))
1351 return error(COMPERR_BAD_SHAPE);
1355 // Find columns of wall nodes and calculate edges' lengths
1356 // --------------------------------------------------------
1358 myParam2ColumnMaps.clear();
1359 myParam2ColumnMaps.resize( orderedEdges.size() ); // total nb edges
1361 int iE, nbEdges = nbEInW.front(); // nb outer edges
1362 vector< double > edgeLength( nbEdges );
1363 map< double, int > len2edgeMap;
1365 list< TopoDS_Edge >::iterator edgeIt = orderedEdges.begin();
1366 list< TopoDS_Face >::iterator faceIt = wallFaces.begin();
1367 for ( iE = 0; iE < nbEdges; ++edgeIt, ++faceIt )
1369 TParam2ColumnMap & faceColumns = myParam2ColumnMaps[ iE ];
1370 if ( !myHelper->LoadNodeColumns( faceColumns, *faceIt, *edgeIt, meshDS ))
1371 return error(COMPERR_BAD_INPUT_MESH, TCom("Can't find regular quadrangle mesh ")
1372 << "on a side face #" << MeshDS()->ShapeToIndex( *faceIt ));
1374 SHOWYXZ("\np1 F "<<iE, gpXYZ(faceColumns.begin()->second.front() ));
1375 SHOWYXZ("p2 F "<<iE, gpXYZ(faceColumns.rbegin()->second.front() ));
1376 SHOWYXZ("V First "<<iE, BRep_Tool::Pnt( TopExp::FirstVertex(*edgeIt,true )));
1378 edgeLength[ iE ] = SMESH_Algo::EdgeLength( *edgeIt );
1380 if ( nbEdges < NB_WALL_FACES ) // fill map used to split faces
1382 SMESHDS_SubMesh* smDS = meshDS->MeshElements( *edgeIt);
1384 return error(COMPERR_BAD_INPUT_MESH, TCom("Null submesh on the edge #")
1385 << MeshDS()->ShapeToIndex( *edgeIt ));
1386 // assure length uniqueness
1387 edgeLength[ iE ] *= smDS->NbNodes() + edgeLength[ iE ] / ( 1000 + iE );
1388 len2edgeMap[ edgeLength[ iE ]] = iE;
1392 // Load columns of internal edges (forming holes)
1393 // and fill map ShapeIndex to TParam2ColumnMap for them
1394 for ( ; edgeIt != orderedEdges.end() ; ++edgeIt, ++faceIt )
1396 TParam2ColumnMap & faceColumns = myParam2ColumnMaps[ iE ];
1397 if ( !myHelper->LoadNodeColumns( faceColumns, *faceIt, *edgeIt, meshDS ))
1398 return error(COMPERR_BAD_INPUT_MESH, TCom("Can't find regular quadrangle mesh ")
1399 << "on a side face #" << MeshDS()->ShapeToIndex( *faceIt ));
1401 int id = MeshDS()->ShapeToIndex( *edgeIt );
1402 bool isForward = true; // meaningless for intenal wires
1403 myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1404 // columns for vertices
1406 const SMDS_MeshNode* n0 = faceColumns.begin()->second.front();
1407 id = n0->getshapeId();
1408 myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1410 const SMDS_MeshNode* n1 = faceColumns.rbegin()->second.front();
1411 id = n1->getshapeId();
1412 myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1413 // SHOWYXZ("\np1 F "<<iE, gpXYZ(faceColumns.begin()->second.front() ));
1414 // SHOWYXZ("p2 F "<<iE, gpXYZ(faceColumns.rbegin()->second.front() ));
1415 // SHOWYXZ("V First "<<iE, BRep_Tool::Pnt( TopExp::FirstVertex(*edgeIt,true )));
1419 // Create 4 wall faces of a block
1420 // -------------------------------
1422 if ( nbEdges <= NB_WALL_FACES ) // ************* Split faces if necessary
1424 map< int, int > iE2nbSplit;
1425 if ( nbEdges != NB_WALL_FACES ) // define how to split
1427 if ( len2edgeMap.size() != nbEdges )
1428 RETURN_BAD_RESULT("Uniqueness of edge lengths not assured");
1429 map< double, int >::reverse_iterator maxLen_i = len2edgeMap.rbegin();
1430 map< double, int >::reverse_iterator midLen_i = ++len2edgeMap.rbegin();
1431 double maxLen = maxLen_i->first;
1432 double midLen = ( len2edgeMap.size() == 1 ) ? 0 : midLen_i->first;
1433 switch ( nbEdges ) {
1434 case 1: // 0-th edge is split into 4 parts
1435 iE2nbSplit.insert( make_pair( 0, 4 )); break;
1436 case 2: // either the longest edge is split into 3 parts, or both edges into halves
1437 if ( maxLen / 3 > midLen / 2 ) {
1438 iE2nbSplit.insert( make_pair( maxLen_i->second, 3 ));
1441 iE2nbSplit.insert( make_pair( maxLen_i->second, 2 ));
1442 iE2nbSplit.insert( make_pair( midLen_i->second, 2 ));
1446 // split longest into halves
1447 iE2nbSplit.insert( make_pair( maxLen_i->second, 2 ));
1450 // Create TSideFace's
1451 faceIt = wallFaces.begin();
1452 edgeIt = orderedEdges.begin();
1454 for ( iE = 0; iE < nbEdges; ++edgeIt, ++faceIt )
1457 map< int, int >::iterator i_nb = iE2nbSplit.find( iE );
1458 if ( i_nb != iE2nbSplit.end() ) {
1460 int nbSplit = i_nb->second;
1461 vector< double > params;
1462 splitParams( nbSplit, &myParam2ColumnMaps[ iE ], params );
1463 bool isForward = ( edgeIt->Orientation() == TopAbs_FORWARD );
1464 for ( int i = 0; i < nbSplit; ++i ) {
1465 double f = ( isForward ? params[ i ] : params[ nbSplit - i-1 ]);
1466 double l = ( isForward ? params[ i+1 ] : params[ nbSplit - i ]);
1467 TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1469 &myParam2ColumnMaps[ iE ], f, l );
1470 mySide->SetComponent( iSide++, comp );
1474 TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1476 &myParam2ColumnMaps[ iE ]);
1477 mySide->SetComponent( iSide++, comp );
1482 else { // **************************** Unite faces
1484 // unite first faces
1485 int nbExraFaces = nbEdges - 3;
1487 double u0 = 0, sumLen = 0;
1488 for ( iE = 0; iE < nbExraFaces; ++iE )
1489 sumLen += edgeLength[ iE ];
1491 vector< TSideFace* > components( nbExraFaces );
1492 vector< pair< double, double> > params( nbExraFaces );
1493 faceIt = wallFaces.begin();
1494 edgeIt = orderedEdges.begin();
1495 for ( iE = 0; iE < nbExraFaces; ++edgeIt, ++faceIt )
1497 components[ iE ] = new TSideFace( myHelper, wallFaceIds[ iSide ],
1499 &myParam2ColumnMaps[ iE ]);
1500 double u1 = u0 + edgeLength[ iE ] / sumLen;
1501 params[ iE ] = make_pair( u0 , u1 );
1505 mySide->SetComponent( iSide++, new TSideFace( components, params ));
1507 // fill the rest faces
1508 for ( ; iE < nbEdges; ++faceIt, ++edgeIt )
1510 TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1512 &myParam2ColumnMaps[ iE ]);
1513 mySide->SetComponent( iSide++, comp );
1519 // Fill geometry fields of SMESH_Block
1520 // ------------------------------------
1522 TopoDS_Face botF = TopoDS::Face( botSM->GetSubShape() );
1523 TopoDS_Face topF = TopoDS::Face( topSM->GetSubShape() );
1525 vector< int > botEdgeIdVec;
1526 SMESH_Block::GetFaceEdgesIDs( ID_BOT_FACE, botEdgeIdVec );
1528 bool isForward[NB_WALL_FACES] = { true, true, true, true };
1529 Adaptor2d_Curve2d* botPcurves[NB_WALL_FACES];
1530 Adaptor2d_Curve2d* topPcurves[NB_WALL_FACES];
1532 for ( int iF = 0; iF < NB_WALL_FACES; ++iF )
1534 TSideFace * sideFace = mySide->GetComponent( iF );
1536 RETURN_BAD_RESULT("NULL TSideFace");
1537 int fID = sideFace->FaceID();
1539 // fill myShapeIDMap
1540 if ( sideFace->InsertSubShapes( myShapeIDMap ) != 8 &&
1541 !sideFace->IsComplex())
1542 MESSAGE( ": Warning : InsertSubShapes() < 8 on side " << iF );
1544 // side faces geometry
1545 Adaptor2d_Curve2d* pcurves[NB_WALL_FACES];
1546 if ( !sideFace->GetPCurves( pcurves ))
1547 RETURN_BAD_RESULT("TSideFace::GetPCurves() failed");
1549 SMESH_Block::TFace& tFace = myFace[ fID - ID_FirstF ];
1550 tFace.Set( fID, sideFace->Surface(), pcurves, isForward );
1552 SHOWYXZ( endl<<"F "<< iF << " id " << fID << " FRW " << sideFace->IsForward(), sideFace->Value(0,0));
1553 // edges 3D geometry
1554 vector< int > edgeIdVec;
1555 SMESH_Block::GetFaceEdgesIDs( fID, edgeIdVec );
1556 for ( int isMax = 0; isMax < 2; ++isMax ) {
1558 int eID = edgeIdVec[ isMax ];
1559 SMESH_Block::TEdge& tEdge = myEdge[ eID - ID_FirstE ];
1560 tEdge.Set( eID, sideFace->HorizCurve(isMax), true);
1561 SHOWYXZ(eID<<" HOR"<<isMax<<"(0)", sideFace->HorizCurve(isMax)->Value(0));
1562 SHOWYXZ(eID<<" HOR"<<isMax<<"(1)", sideFace->HorizCurve(isMax)->Value(1));
1565 int eID = edgeIdVec[ isMax+2 ];
1566 SMESH_Block::TEdge& tEdge = myEdge[ eID - ID_FirstE ];
1567 tEdge.Set( eID, sideFace->VertiCurve(isMax), true);
1568 SHOWYXZ(eID<<" VER"<<isMax<<"(0)", sideFace->VertiCurve(isMax)->Value(0));
1569 SHOWYXZ(eID<<" VER"<<isMax<<"(1)", sideFace->VertiCurve(isMax)->Value(1));
1572 vector< int > vertexIdVec;
1573 SMESH_Block::GetEdgeVertexIDs( eID, vertexIdVec );
1574 myPnt[ vertexIdVec[0] - ID_FirstV ] = tEdge.GetCurve()->Value(0).XYZ();
1575 myPnt[ vertexIdVec[1] - ID_FirstV ] = tEdge.GetCurve()->Value(1).XYZ();
1578 // pcurves on horizontal faces
1579 for ( iE = 0; iE < NB_WALL_FACES; ++iE ) {
1580 if ( edgeIdVec[ BOTTOM_EDGE ] == botEdgeIdVec[ iE ] ) {
1581 botPcurves[ iE ] = sideFace->HorizPCurve( false, botF );
1582 topPcurves[ iE ] = sideFace->HorizPCurve( true, topF );
1586 //sideFace->dumpNodes( 4 ); // debug
1588 // horizontal faces geometry
1590 SMESH_Block::TFace& tFace = myFace[ ID_BOT_FACE - ID_FirstF ];
1591 tFace.Set( ID_BOT_FACE, new BRepAdaptor_Surface( botF ), botPcurves, isForward );
1592 SMESH_Block::Insert( botF, ID_BOT_FACE, myShapeIDMap );
1595 SMESH_Block::TFace& tFace = myFace[ ID_TOP_FACE - ID_FirstF ];
1596 tFace.Set( ID_TOP_FACE, new BRepAdaptor_Surface( topF ), topPcurves, isForward );
1597 SMESH_Block::Insert( topF, ID_TOP_FACE, myShapeIDMap );
1600 // Fill map ShapeIndex to TParam2ColumnMap
1601 // ----------------------------------------
1603 list< TSideFace* > fList;
1604 list< TSideFace* >::iterator fListIt;
1605 fList.push_back( mySide );
1606 for ( fListIt = fList.begin(); fListIt != fList.end(); ++fListIt)
1608 int nb = (*fListIt)->NbComponents();
1609 for ( int i = 0; i < nb; ++i ) {
1610 if ( TSideFace* comp = (*fListIt)->GetComponent( i ))
1611 fList.push_back( comp );
1613 if ( TParam2ColumnMap* cols = (*fListIt)->GetColumns()) {
1614 // columns for a base edge
1615 int id = MeshDS()->ShapeToIndex( (*fListIt)->BaseEdge() );
1616 bool isForward = (*fListIt)->IsForward();
1617 myShapeIndex2ColumnMap[ id ] = make_pair( cols, isForward );
1619 // columns for vertices
1620 const SMDS_MeshNode* n0 = cols->begin()->second.front();
1621 id = n0->getshapeId();
1622 myShapeIndex2ColumnMap[ id ] = make_pair( cols, isForward );
1624 const SMDS_MeshNode* n1 = cols->rbegin()->second.front();
1625 id = n1->getshapeId();
1626 myShapeIndex2ColumnMap[ id ] = make_pair( cols, !isForward );
1630 // gp_XYZ testPar(0.25, 0.25, 0), testCoord;
1631 // if ( !FacePoint( ID_BOT_FACE, testPar, testCoord ))
1632 // RETURN_BAD_RESULT("TEST FacePoint() FAILED");
1633 // SHOWYXZ("IN TEST PARAM" , testPar);
1634 // SHOWYXZ("OUT TEST CORD" , testCoord);
1635 // if ( !ComputeParameters( testCoord, testPar , ID_BOT_FACE))
1636 // RETURN_BAD_RESULT("TEST ComputeParameters() FAILED");
1637 // SHOWYXZ("OUT TEST PARAM" , testPar);
1642 //================================================================================
1644 * \brief Return pointer to column of nodes
1645 * \param node - bottom node from which the returned column goes up
1646 * \retval const TNodeColumn* - the found column
1648 //================================================================================
1650 const TNodeColumn* StdMeshers_PrismAsBlock::GetNodeColumn(const SMDS_MeshNode* node) const
1652 int sID = node->getshapeId();
1654 map<int, pair< TParam2ColumnMap*, bool > >::const_iterator col_frw =
1655 myShapeIndex2ColumnMap.find( sID );
1656 if ( col_frw != myShapeIndex2ColumnMap.end() ) {
1657 const TParam2ColumnMap* cols = col_frw->second.first;
1658 TParam2ColumnIt u_col = cols->begin();
1659 for ( ; u_col != cols->end(); ++u_col )
1660 if ( u_col->second[ 0 ] == node )
1661 return & u_col->second;
1666 //=======================================================================
1667 //function : GetLayersTransformation
1668 //purpose : Return transformations to get coordinates of nodes of each layer
1669 // by nodes of the bottom. Layer is a set of nodes at a certain step
1670 // from bottom to top.
1671 //=======================================================================
1673 bool StdMeshers_PrismAsBlock::GetLayersTransformation(vector<gp_Trsf> & trsf) const
1675 const int zSize = VerticalSize();
1676 if ( zSize < 3 ) return true;
1677 trsf.resize( zSize - 2 );
1679 // Select some node columns by which we will define coordinate system of layers
1681 vector< const TNodeColumn* > columns;
1683 const TopoDS_Shape& baseFace = Shape(ID_BOT_FACE);
1684 list< TopoDS_Edge > orderedEdges;
1685 list< int > nbEdgesInWires;
1686 GetOrderedEdges( TopoDS::Face( baseFace ), TopoDS_Vertex(), orderedEdges, nbEdgesInWires );
1688 list< TopoDS_Edge >::iterator edgeIt = orderedEdges.begin();
1689 for ( int iE = 0; iE < nbEdgesInWires.front(); ++iE, ++edgeIt )
1691 if ( BRep_Tool::Degenerated( *edgeIt )) continue;
1692 const TParam2ColumnMap& u2colMap =
1693 GetParam2ColumnMap( myHelper->GetMeshDS()->ShapeToIndex( *edgeIt ), isReverse );
1694 isReverse = ( edgeIt->Orientation() == TopAbs_REVERSED );
1695 double f = u2colMap.begin()->first, l = u2colMap.rbegin()->first;
1696 if ( isReverse ) swap ( f, l );
1697 const int nbCol = 5;
1698 for ( int i = 0; i < nbCol; ++i )
1700 double u = f + i/double(nbCol) * ( l - f );
1701 const TNodeColumn* col = & getColumn( & u2colMap, u )->second;
1702 if ( columns.empty() || col != columns.back() )
1703 columns.push_back( col );
1708 // Find tolerance to check transformations
1713 for ( int i = 0; i < columns.size(); ++i )
1714 bndBox.Add( gpXYZ( columns[i]->front() ));
1715 tol2 = bndBox.SquareExtent() * 1e-5;
1718 // Compute transformations
1721 gp_Trsf fromCsZ, toCs0;
1722 gp_Ax3 cs0 = getLayerCoordSys(0, columns, xCol );
1723 //double dist0 = cs0.Location().Distance( gpXYZ( (*columns[0])[0]));
1724 toCs0.SetTransformation( cs0 );
1725 for ( int z = 1; z < zSize-1; ++z )
1727 gp_Ax3 csZ = getLayerCoordSys(z, columns, xCol );
1728 //double distZ = csZ.Location().Distance( gpXYZ( (*columns[0])[z]));
1729 fromCsZ.SetTransformation( csZ );
1731 gp_Trsf& t = trsf[ z-1 ];
1732 t = fromCsZ * toCs0;
1733 //t.SetScaleFactor( distZ/dist0 ); - it does not work properly, wrong base point
1735 // check a transformation
1736 for ( int i = 0; i < columns.size(); ++i )
1738 gp_Pnt p0 = gpXYZ( (*columns[i])[0] );
1739 gp_Pnt pz = gpXYZ( (*columns[i])[z] );
1740 t.Transforms( p0.ChangeCoord() );
1741 if ( p0.SquareDistance( pz ) > tol2 )
1748 //================================================================================
1750 * \brief Check curve orientation of a bootom edge
1751 * \param meshDS - mesh DS
1752 * \param columnsMap - node columns map of side face
1753 * \param bottomEdge - the bootom edge
1754 * \param sideFaceID - side face in-block ID
1755 * \retval bool - true if orientation coinside with in-block froward orientation
1757 //================================================================================
1759 bool StdMeshers_PrismAsBlock::IsForwardEdge(SMESHDS_Mesh* meshDS,
1760 const TParam2ColumnMap& columnsMap,
1761 const TopoDS_Edge & bottomEdge,
1762 const int sideFaceID)
1764 bool isForward = false;
1765 if ( SMESH_MesherHelper::IsClosedEdge( bottomEdge ))
1767 isForward = ( bottomEdge.Orientation() == TopAbs_FORWARD );
1771 const TNodeColumn& firstCol = columnsMap.begin()->second;
1772 const SMDS_MeshNode* bottomNode = firstCol[0];
1773 TopoDS_Shape firstVertex = SMESH_MesherHelper::GetSubShapeByNode( bottomNode, meshDS );
1774 isForward = ( firstVertex.IsSame( TopExp::FirstVertex( bottomEdge, true )));
1776 // on 2 of 4 sides first vertex is end
1777 if ( sideFaceID == ID_Fx1z || sideFaceID == ID_F0yz )
1778 isForward = !isForward;
1782 //================================================================================
1784 * \brief Find wall faces by bottom edges
1785 * \param mesh - the mesh
1786 * \param mainShape - the prism
1787 * \param bottomFace - the bottom face
1788 * \param bottomEdges - edges bounding the bottom face
1789 * \param wallFaces - faces list to fill in
1791 //================================================================================
1793 bool StdMeshers_PrismAsBlock::GetWallFaces( SMESH_Mesh* mesh,
1794 const TopoDS_Shape & mainShape,
1795 const TopoDS_Shape & bottomFace,
1796 std::list< TopoDS_Edge >& bottomEdges,
1797 std::list< int > & nbEInW,
1798 std::list< TopoDS_Face >& wallFaces)
1802 TopTools_IndexedMapOfShape faceMap;
1803 TopExp::MapShapes( mainShape, TopAbs_FACE, faceMap );
1805 list< TopoDS_Edge >::iterator edge = bottomEdges.begin();
1806 std::list< int >::iterator nbE = nbEInW.begin();
1808 while ( edge != bottomEdges.end() )
1811 if ( BRep_Tool::Degenerated( *edge ))
1813 edge = bottomEdges.erase( edge );
1819 PShapeIteratorPtr fIt = myHelper->GetAncestors( *edge, *mesh, TopAbs_FACE );
1820 while ( fIt->more() )
1822 const TopoDS_Shape* face = fIt->next();
1823 if ( !bottomFace.IsSame( *face ) && // not bottom
1824 faceMap.FindIndex( *face )) // belongs to the prism
1826 wallFaces.push_back( TopoDS::Face( *face ));
1838 return ( wallFaces.size() == bottomEdges.size() );
1841 //================================================================================
1843 * \brief Constructor
1844 * \param faceID - in-block ID
1845 * \param face - geom face
1846 * \param columnsMap - map of node columns
1847 * \param first - first normalized param
1848 * \param last - last normalized param
1850 //================================================================================
1852 StdMeshers_PrismAsBlock::TSideFace::TSideFace(SMESH_MesherHelper* helper,
1854 const TopoDS_Face& face,
1855 const TopoDS_Edge& baseEdge,
1856 TParam2ColumnMap* columnsMap,
1860 myParamToColumnMap( columnsMap ),
1861 myBaseEdge( baseEdge ),
1864 mySurface.Initialize( face );
1865 myParams.resize( 1 );
1866 myParams[ 0 ] = make_pair( first, last );
1867 myIsForward = StdMeshers_PrismAsBlock::IsForwardEdge( myHelper->GetMeshDS(),
1868 *myParamToColumnMap,
1872 //================================================================================
1874 * \brief Constructor of complex side face
1876 //================================================================================
1878 StdMeshers_PrismAsBlock::TSideFace::
1879 TSideFace(const vector< TSideFace* >& components,
1880 const vector< pair< double, double> > & params)
1881 :myID( components[0] ? components[0]->myID : 0 ),
1882 myParamToColumnMap( 0 ),
1884 myIsForward( true ),
1885 myComponents( components ),
1886 myHelper( components[0] ? components[0]->myHelper : 0 )
1888 //================================================================================
1890 * \brief Copy constructor
1891 * \param other - other side
1893 //================================================================================
1895 StdMeshers_PrismAsBlock::TSideFace::TSideFace( const TSideFace& other )
1898 mySurface = other.mySurface;
1899 myBaseEdge = other.myBaseEdge;
1900 myParams = other.myParams;
1901 myIsForward = other.myIsForward;
1902 myHelper = other.myHelper;
1903 myParamToColumnMap = other.myParamToColumnMap;
1905 myComponents.resize( other.myComponents.size());
1906 for (int i = 0 ; i < myComponents.size(); ++i )
1907 myComponents[ i ] = new TSideFace( *other.myComponents[ i ]);
1910 //================================================================================
1912 * \brief Deletes myComponents
1914 //================================================================================
1916 StdMeshers_PrismAsBlock::TSideFace::~TSideFace()
1918 for (int i = 0 ; i < myComponents.size(); ++i )
1919 if ( myComponents[ i ] )
1920 delete myComponents[ i ];
1923 //================================================================================
1925 * \brief Return geometry of the vertical curve
1926 * \param isMax - true means curve located closer to (1,1,1) block point
1927 * \retval Adaptor3d_Curve* - curve adaptor
1929 //================================================================================
1931 Adaptor3d_Curve* StdMeshers_PrismAsBlock::TSideFace::VertiCurve(const bool isMax) const
1933 if ( !myComponents.empty() ) {
1935 return myComponents.back()->VertiCurve(isMax);
1937 return myComponents.front()->VertiCurve(isMax);
1939 double f = myParams[0].first, l = myParams[0].second;
1940 if ( !myIsForward ) std::swap( f, l );
1941 return new TVerticalEdgeAdaptor( myParamToColumnMap, isMax ? l : f );
1944 //================================================================================
1946 * \brief Return geometry of the top or bottom curve
1948 * \retval Adaptor3d_Curve* -
1950 //================================================================================
1952 Adaptor3d_Curve* StdMeshers_PrismAsBlock::TSideFace::HorizCurve(const bool isTop) const
1954 return new THorizontalEdgeAdaptor( this, isTop );
1957 //================================================================================
1959 * \brief Return pcurves
1960 * \param pcurv - array of 4 pcurves
1961 * \retval bool - is a success
1963 //================================================================================
1965 bool StdMeshers_PrismAsBlock::TSideFace::GetPCurves(Adaptor2d_Curve2d* pcurv[4]) const
1967 int iEdge[ 4 ] = { BOTTOM_EDGE, TOP_EDGE, V0_EDGE, V1_EDGE };
1969 for ( int i = 0 ; i < 4 ; ++i ) {
1970 Handle(Geom2d_Line) line;
1971 switch ( iEdge[ i ] ) {
1973 line = new Geom2d_Line( gp_Pnt2d( 0, 1 ), gp::DX2d() ); break;
1975 line = new Geom2d_Line( gp::Origin2d(), gp::DX2d() ); break;
1977 line = new Geom2d_Line( gp::Origin2d(), gp::DY2d() ); break;
1979 line = new Geom2d_Line( gp_Pnt2d( 1, 0 ), gp::DY2d() ); break;
1981 pcurv[ i ] = new Geom2dAdaptor_Curve( line, 0, 1 );
1986 //================================================================================
1988 * \brief Returns geometry of pcurve on a horizontal face
1989 * \param isTop - is top or bottom face
1990 * \param horFace - a horizontal face
1991 * \retval Adaptor2d_Curve2d* - curve adaptor
1993 //================================================================================
1996 StdMeshers_PrismAsBlock::TSideFace::HorizPCurve(const bool isTop,
1997 const TopoDS_Face& horFace) const
1999 return new TPCurveOnHorFaceAdaptor( this, isTop, horFace );
2002 //================================================================================
2004 * \brief Return a component corresponding to parameter
2005 * \param U - parameter along a horizontal size
2006 * \param localU - parameter along a horizontal size of a component
2007 * \retval TSideFace* - found component
2009 //================================================================================
2011 StdMeshers_PrismAsBlock::TSideFace*
2012 StdMeshers_PrismAsBlock::TSideFace::GetComponent(const double U,double & localU) const
2015 if ( myComponents.empty() )
2016 return const_cast<TSideFace*>( this );
2019 for ( i = 0; i < myComponents.size(); ++i )
2020 if ( U < myParams[ i ].second )
2022 if ( i >= myComponents.size() )
2023 i = myComponents.size() - 1;
2025 double f = myParams[ i ].first, l = myParams[ i ].second;
2026 localU = ( U - f ) / ( l - f );
2027 return myComponents[ i ];
2030 //================================================================================
2032 * \brief Find node columns for a parameter
2033 * \param U - parameter along a horizontal edge
2034 * \param col1 - the 1st found column
2035 * \param col2 - the 2nd found column
2036 * \retval r - normalized position of U between the found columns
2038 //================================================================================
2040 double StdMeshers_PrismAsBlock::TSideFace::GetColumns(const double U,
2041 TParam2ColumnIt & col1,
2042 TParam2ColumnIt & col2) const
2044 double u = U, r = 0;
2045 if ( !myComponents.empty() ) {
2046 TSideFace * comp = GetComponent(U,u);
2047 return comp->GetColumns( u, col1, col2 );
2052 double f = myParams[0].first, l = myParams[0].second;
2053 u = f + u * ( l - f );
2055 col1 = col2 = getColumn( myParamToColumnMap, u );
2056 if ( ++col2 == myParamToColumnMap->end() ) {
2061 double uf = col1->first;
2062 double ul = col2->first;
2063 r = ( u - uf ) / ( ul - uf );
2068 //================================================================================
2070 * \brief Return coordinates by normalized params
2071 * \param U - horizontal param
2072 * \param V - vertical param
2073 * \retval gp_Pnt - result point
2075 //================================================================================
2077 gp_Pnt StdMeshers_PrismAsBlock::TSideFace::Value(const Standard_Real U,
2078 const Standard_Real V) const
2080 if ( !myComponents.empty() ) {
2082 TSideFace * comp = GetComponent(U,u);
2083 return comp->Value( u, V );
2086 TParam2ColumnIt u_col1, u_col2;
2087 double vR, hR = GetColumns( U, u_col1, u_col2 );
2089 const SMDS_MeshNode* n1 = 0;
2090 const SMDS_MeshNode* n2 = 0;
2091 const SMDS_MeshNode* n3 = 0;
2092 const SMDS_MeshNode* n4 = 0;
2094 // BEGIN issue 0020680: EDF 1252 SMESH: Bad cell created by Radial prism in center of torus
2095 // Workaround for a wrongly located point returned by mySurface.Value() for
2096 // UV located near boundary of BSpline surface.
2097 // To bypass the problem, we take point from 3D curve of edge.
2098 // It solves pb of the bloc_fiss_new.py
2099 const double tol = 1e-3;
2100 if ( V < tol || V+tol >= 1. )
2102 n1 = V < tol ? u_col1->second.front() : u_col1->second.back();
2103 n3 = V < tol ? u_col2->second.front() : u_col2->second.back();
2111 TopoDS_Shape s = myHelper->GetSubShapeByNode( n1, myHelper->GetMeshDS() );
2112 if ( s.ShapeType() != TopAbs_EDGE )
2113 s = myHelper->GetSubShapeByNode( n3, myHelper->GetMeshDS() );
2114 if ( s.ShapeType() == TopAbs_EDGE )
2115 edge = TopoDS::Edge( s );
2117 if ( !edge.IsNull() )
2119 double u1 = myHelper->GetNodeU( edge, n1 );
2120 double u3 = myHelper->GetNodeU( edge, n3 );
2121 double u = u1 * ( 1 - hR ) + u3 * hR;
2122 TopLoc_Location loc; double f,l;
2123 Handle(Geom_Curve) curve = BRep_Tool::Curve( edge,loc,f,l );
2124 return curve->Value( u ).Transformed( loc );
2127 // END issue 0020680: EDF 1252 SMESH: Bad cell created by Radial prism in center of torus
2129 vR = getRAndNodes( & u_col1->second, V, n1, n2 );
2130 vR = getRAndNodes( & u_col2->second, V, n3, n4 );
2132 gp_XY uv1 = myHelper->GetNodeUV( mySurface.Face(), n1, n4);
2133 gp_XY uv2 = myHelper->GetNodeUV( mySurface.Face(), n2, n3);
2134 gp_XY uv12 = uv1 * ( 1 - vR ) + uv2 * vR;
2136 gp_XY uv3 = myHelper->GetNodeUV( mySurface.Face(), n3, n2);
2137 gp_XY uv4 = myHelper->GetNodeUV( mySurface.Face(), n4, n1);
2138 gp_XY uv34 = uv3 * ( 1 - vR ) + uv4 * vR;
2140 gp_XY uv = uv12 * ( 1 - hR ) + uv34 * hR;
2142 gp_Pnt p = mySurface.Value( uv.X(), uv.Y() );
2147 //================================================================================
2149 * \brief Return boundary edge
2150 * \param edge - edge index
2151 * \retval TopoDS_Edge - found edge
2153 //================================================================================
2155 TopoDS_Edge StdMeshers_PrismAsBlock::TSideFace::GetEdge(const int iEdge) const
2157 if ( !myComponents.empty() ) {
2159 case V0_EDGE : return myComponents.front()->GetEdge( iEdge );
2160 case V1_EDGE : return myComponents.back() ->GetEdge( iEdge );
2161 default: return TopoDS_Edge();
2165 const SMDS_MeshNode* node = 0;
2166 SMESHDS_Mesh * meshDS = myHelper->GetMesh()->GetMeshDS();
2167 TNodeColumn* column;
2172 column = & (( ++myParamToColumnMap->begin())->second );
2173 node = ( iEdge == TOP_EDGE ) ? column->back() : column->front();
2174 edge = myHelper->GetSubShapeByNode ( node, meshDS );
2175 if ( edge.ShapeType() == TopAbs_VERTEX ) {
2176 column = & ( myParamToColumnMap->begin()->second );
2177 node = ( iEdge == TOP_EDGE ) ? column->back() : column->front();
2182 bool back = ( iEdge == V1_EDGE );
2183 if ( !myIsForward ) back = !back;
2185 column = & ( myParamToColumnMap->rbegin()->second );
2187 column = & ( myParamToColumnMap->begin()->second );
2188 if ( column->size() > 0 )
2189 edge = myHelper->GetSubShapeByNode( (*column)[ 1 ], meshDS );
2190 if ( edge.IsNull() || edge.ShapeType() == TopAbs_VERTEX )
2191 node = column->front();
2196 if ( !edge.IsNull() && edge.ShapeType() == TopAbs_EDGE )
2197 return TopoDS::Edge( edge );
2199 // find edge by 2 vertices
2200 TopoDS_Shape V1 = edge;
2201 TopoDS_Shape V2 = myHelper->GetSubShapeByNode( node, meshDS );
2202 if ( V2.ShapeType() == TopAbs_VERTEX && !V2.IsSame( V1 ))
2204 TopoDS_Shape ancestor = myHelper->GetCommonAncestor( V1, V2, *myHelper->GetMesh(), TopAbs_EDGE);
2205 if ( !ancestor.IsNull() )
2206 return TopoDS::Edge( ancestor );
2208 return TopoDS_Edge();
2211 //================================================================================
2213 * \brief Fill block subshapes
2214 * \param shapeMap - map to fill in
2215 * \retval int - nb inserted subshapes
2217 //================================================================================
2219 int StdMeshers_PrismAsBlock::TSideFace::InsertSubShapes(TBlockShapes& shapeMap) const
2224 vector< int > edgeIdVec;
2225 SMESH_Block::GetFaceEdgesIDs( myID, edgeIdVec );
2227 for ( int i = BOTTOM_EDGE; i <=V1_EDGE ; ++i ) {
2228 TopoDS_Edge e = GetEdge( i );
2229 if ( !e.IsNull() ) {
2230 nbInserted += SMESH_Block::Insert( e, edgeIdVec[ i ], shapeMap);
2234 // Insert corner vertices
2236 TParam2ColumnIt col1, col2 ;
2237 vector< int > vertIdVec;
2240 SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ V0_EDGE ], vertIdVec);
2241 GetColumns(0, col1, col2 );
2242 const SMDS_MeshNode* node0 = col1->second.front();
2243 const SMDS_MeshNode* node1 = col1->second.back();
2244 TopoDS_Shape v0 = myHelper->GetSubShapeByNode( node0, myHelper->GetMeshDS());
2245 TopoDS_Shape v1 = myHelper->GetSubShapeByNode( node1, myHelper->GetMeshDS());
2246 if ( v0.ShapeType() == TopAbs_VERTEX ) {
2247 nbInserted += SMESH_Block::Insert( v0, vertIdVec[ 0 ], shapeMap);
2249 if ( v1.ShapeType() == TopAbs_VERTEX ) {
2250 nbInserted += SMESH_Block::Insert( v1, vertIdVec[ 1 ], shapeMap);
2254 SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ V1_EDGE ], vertIdVec);
2255 GetColumns(1, col1, col2 );
2256 node0 = col2->second.front();
2257 node1 = col2->second.back();
2258 v0 = myHelper->GetSubShapeByNode( node0, myHelper->GetMeshDS());
2259 v1 = myHelper->GetSubShapeByNode( node1, myHelper->GetMeshDS());
2260 if ( v0.ShapeType() == TopAbs_VERTEX ) {
2261 nbInserted += SMESH_Block::Insert( v0, vertIdVec[ 0 ], shapeMap);
2263 if ( v1.ShapeType() == TopAbs_VERTEX ) {
2264 nbInserted += SMESH_Block::Insert( v1, vertIdVec[ 1 ], shapeMap);
2267 // TopoDS_Vertex V0, V1, Vcom;
2268 // TopExp::Vertices( myBaseEdge, V0, V1, true );
2269 // if ( !myIsForward ) std::swap( V0, V1 );
2271 // // bottom vertex IDs
2272 // SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ _u0 ], vertIdVec);
2273 // SMESH_Block::Insert( V0, vertIdVec[ 0 ], shapeMap);
2274 // SMESH_Block::Insert( V1, vertIdVec[ 1 ], shapeMap);
2276 // TopoDS_Edge sideEdge = GetEdge( V0_EDGE );
2277 // if ( sideEdge.IsNull() || !TopExp::CommonVertex( botEdge, sideEdge, Vcom ))
2280 // // insert one side edge
2282 // if ( Vcom.IsSame( V0 )) edgeID = edgeIdVec[ _v0 ];
2283 // else edgeID = edgeIdVec[ _v1 ];
2284 // SMESH_Block::Insert( sideEdge, edgeID, shapeMap);
2286 // // top vertex of the side edge
2287 // SMESH_Block::GetEdgeVertexIDs( edgeID, vertIdVec);
2288 // TopoDS_Vertex Vtop = TopExp::FirstVertex( sideEdge );
2289 // if ( Vcom.IsSame( Vtop ))
2290 // Vtop = TopExp::LastVertex( sideEdge );
2291 // SMESH_Block::Insert( Vtop, vertIdVec[ 1 ], shapeMap);
2293 // // other side edge
2294 // sideEdge = GetEdge( V1_EDGE );
2295 // if ( sideEdge.IsNull() )
2297 // if ( edgeID = edgeIdVec[ _v1 ]) edgeID = edgeIdVec[ _v0 ];
2298 // else edgeID = edgeIdVec[ _v1 ];
2299 // SMESH_Block::Insert( sideEdge, edgeID, shapeMap);
2302 // TopoDS_Edge topEdge = GetEdge( TOP_EDGE );
2303 // SMESH_Block::Insert( topEdge, edgeIdVec[ _u1 ], shapeMap);
2305 // // top vertex of the other side edge
2306 // if ( !TopExp::CommonVertex( topEdge, sideEdge, Vcom ))
2308 // SMESH_Block::GetEdgeVertexIDs( edgeID, vertIdVec );
2309 // SMESH_Block::Insert( Vcom, vertIdVec[ 1 ], shapeMap);
2314 //================================================================================
2316 * \brief Dump ids of nodes of sides
2318 //================================================================================
2320 void StdMeshers_PrismAsBlock::TSideFace::dumpNodes(int nbNodes) const
2323 cout << endl << "NODES OF FACE "; SMESH_Block::DumpShapeID( myID, cout ) << endl;
2324 THorizontalEdgeAdaptor* hSize0 = (THorizontalEdgeAdaptor*) HorizCurve(0);
2325 cout << "Horiz side 0: "; hSize0->dumpNodes(nbNodes); cout << endl;
2326 THorizontalEdgeAdaptor* hSize1 = (THorizontalEdgeAdaptor*) HorizCurve(1);
2327 cout << "Horiz side 1: "; hSize1->dumpNodes(nbNodes); cout << endl;
2328 TVerticalEdgeAdaptor* vSide0 = (TVerticalEdgeAdaptor*) VertiCurve(0);
2329 cout << "Verti side 0: "; vSide0->dumpNodes(nbNodes); cout << endl;
2330 TVerticalEdgeAdaptor* vSide1 = (TVerticalEdgeAdaptor*) VertiCurve(1);
2331 cout << "Verti side 1: "; vSide1->dumpNodes(nbNodes); cout << endl;
2332 delete hSize0; delete hSize1; delete vSide0; delete vSide1;
2336 //================================================================================
2338 * \brief Creates TVerticalEdgeAdaptor
2339 * \param columnsMap - node column map
2340 * \param parameter - normalized parameter
2342 //================================================================================
2344 StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::
2345 TVerticalEdgeAdaptor( const TParam2ColumnMap* columnsMap, const double parameter)
2347 myNodeColumn = & getColumn( columnsMap, parameter )->second;
2350 //================================================================================
2352 * \brief Return coordinates for the given normalized parameter
2353 * \param U - normalized parameter
2354 * \retval gp_Pnt - coordinates
2356 //================================================================================
2358 gp_Pnt StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::Value(const Standard_Real U) const
2360 const SMDS_MeshNode* n1;
2361 const SMDS_MeshNode* n2;
2362 double r = getRAndNodes( myNodeColumn, U, n1, n2 );
2363 return gpXYZ(n1) * ( 1 - r ) + gpXYZ(n2) * r;
2366 //================================================================================
2368 * \brief Dump ids of nodes
2370 //================================================================================
2372 void StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::dumpNodes(int nbNodes) const
2375 for ( int i = 0; i < nbNodes && i < myNodeColumn->size(); ++i )
2376 cout << (*myNodeColumn)[i]->GetID() << " ";
2377 if ( nbNodes < myNodeColumn->size() )
2378 cout << myNodeColumn->back()->GetID();
2382 //================================================================================
2384 * \brief Return coordinates for the given normalized parameter
2385 * \param U - normalized parameter
2386 * \retval gp_Pnt - coordinates
2388 //================================================================================
2390 gp_Pnt StdMeshers_PrismAsBlock::THorizontalEdgeAdaptor::Value(const Standard_Real U) const
2392 return mySide->TSideFace::Value( U, myV );
2395 //================================================================================
2397 * \brief Dump ids of <nbNodes> first nodes and the last one
2399 //================================================================================
2401 void StdMeshers_PrismAsBlock::THorizontalEdgeAdaptor::dumpNodes(int nbNodes) const
2404 // Not bedugged code. Last node is sometimes incorrect
2405 const TSideFace* side = mySide;
2407 if ( mySide->IsComplex() )
2408 side = mySide->GetComponent(0,u);
2410 TParam2ColumnIt col, col2;
2411 TParam2ColumnMap* u2cols = side->GetColumns();
2412 side->GetColumns( u , col, col2 );
2414 int j, i = myV ? mySide->ColumnHeight()-1 : 0;
2416 const SMDS_MeshNode* n = 0;
2417 const SMDS_MeshNode* lastN
2418 = side->IsForward() ? u2cols->rbegin()->second[ i ] : u2cols->begin()->second[ i ];
2419 for ( j = 0; j < nbNodes && n != lastN; ++j )
2421 n = col->second[ i ];
2422 cout << n->GetID() << " ";
2423 if ( side->IsForward() )
2431 if ( mySide->IsComplex() )
2432 side = mySide->GetComponent(1,u);
2434 side->GetColumns( u , col, col2 );
2435 if ( n != col->second[ i ] )
2436 cout << col->second[ i ]->GetID();
2439 //================================================================================
2441 * \brief Return UV on pcurve for the given normalized parameter
2442 * \param U - normalized parameter
2443 * \retval gp_Pnt - coordinates
2445 //================================================================================
2447 gp_Pnt2d StdMeshers_PrismAsBlock::TPCurveOnHorFaceAdaptor::Value(const Standard_Real U) const
2449 TParam2ColumnIt u_col1, u_col2;
2450 double r = mySide->GetColumns( U, u_col1, u_col2 );
2451 gp_XY uv1 = mySide->GetNodeUV( myFace, u_col1->second[ myZ ]);
2452 gp_XY uv2 = mySide->GetNodeUV( myFace, u_col2->second[ myZ ]);
2453 return uv1 * ( 1 - r ) + uv2 * r;