1 // Copyright (C) 2007-2008 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
22 // SMESH SMESH : implementaion of SMESH idl descriptions
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);
218 //=======================================================================
219 //function : StdMeshers_Prism_3D
221 //=======================================================================
223 StdMeshers_Prism_3D::StdMeshers_Prism_3D(int hypId, int studyId, SMESH_Gen* gen)
224 :SMESH_3D_Algo(hypId, studyId, gen)
227 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID); // 1 bit per shape type
228 myProjectTriangles = false;
231 //================================================================================
235 //================================================================================
237 StdMeshers_Prism_3D::~StdMeshers_Prism_3D()
240 //=======================================================================
241 //function : CheckHypothesis
243 //=======================================================================
245 bool StdMeshers_Prism_3D::CheckHypothesis(SMESH_Mesh& aMesh,
246 const TopoDS_Shape& aShape,
247 SMESH_Hypothesis::Hypothesis_Status& aStatus)
249 // Check shape geometry
251 aStatus = SMESH_Hypothesis::HYP_BAD_GEOMETRY;
253 // find not quadrangle faces
254 list< TopoDS_Shape > notQuadFaces;
255 int nbEdge, nbWire, nbFace = 0;
256 TopExp_Explorer exp( aShape, TopAbs_FACE );
257 for ( ; exp.More(); exp.Next() ) {
259 const TopoDS_Shape& face = exp.Current();
260 nbEdge = TAssocTool::Count( face, TopAbs_EDGE, 0 );
261 nbWire = TAssocTool::Count( face, TopAbs_WIRE, 0 );
262 if ( nbEdge!= 4 || nbWire!= 1 ) {
263 if ( !notQuadFaces.empty() ) {
264 if ( TAssocTool::Count( notQuadFaces.back(), TopAbs_EDGE, 0 ) != nbEdge ||
265 TAssocTool::Count( notQuadFaces.back(), TopAbs_WIRE, 0 ) != nbWire )
266 RETURN_BAD_RESULT("Different not quad faces");
268 notQuadFaces.push_back( face );
271 if ( !notQuadFaces.empty() )
273 if ( notQuadFaces.size() != 2 )
274 RETURN_BAD_RESULT("Bad nb not quad faces: " << notQuadFaces.size());
276 // check total nb faces
277 nbEdge = TAssocTool::Count( notQuadFaces.back(), TopAbs_EDGE, 0 );
278 if ( nbFace != nbEdge + 2 )
279 RETURN_BAD_RESULT("Bad nb of faces: " << nbFace << " but must be " << nbEdge + 2);
283 aStatus = SMESH_Hypothesis::HYP_OK;
287 //=======================================================================
290 //=======================================================================
292 bool StdMeshers_Prism_3D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theShape)
294 SMESH_MesherHelper helper( theMesh );
297 myHelper->IsQuadraticSubMesh( theShape );
299 // Analyse mesh and geomerty to find block subshapes and submeshes
300 if ( !myBlock.Init( myHelper, theShape ))
301 return error( myBlock.GetError());
303 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
305 int volumeID = meshDS->ShapeToIndex( theShape );
308 // To compute coordinates of a node inside a block, it is necessary to know
309 // 1. normalized parameters of the node by which
310 // 2. coordinates of node projections on all block sub-shapes are computed
312 // So we fill projections on vertices at once as they are same for all nodes
313 myShapeXYZ.resize( myBlock.NbSubShapes() );
314 for ( int iV = SMESH_Block::ID_FirstV; iV < SMESH_Block::ID_FirstE; ++iV ) {
315 myBlock.VertexPoint( iV, myShapeXYZ[ iV ]);
316 SHOWYXZ("V point " <<iV << " ", myShapeXYZ[ iV ]);
319 // Projections on the top and bottom faces are taken from nodes existing
320 // on these faces; find correspondence between bottom and top nodes
321 myBotToColumnMap.clear();
322 if ( !assocOrProjBottom2Top() ) // it also fill myBotToColumnMap
326 // Create nodes inside the block
328 // try to use transformation (issue 0020680)
329 vector<gp_Trsf> trsf;
330 if ( myBlock.GetLayersTransformation(trsf))
332 // loop on nodes inside the bottom face
333 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.begin();
334 for ( ; bot_column != myBotToColumnMap.end(); ++bot_column )
336 const TNode& tBotNode = bot_column->first; // bottom TNode
337 if ( tBotNode.GetPositionType() != SMDS_TOP_FACE )
338 continue; // node is not inside face
340 // column nodes; middle part of the column are zero pointers
341 TNodeColumn& column = bot_column->second;
342 TNodeColumn::iterator columnNodes = column.begin();
343 for ( int z = 0; columnNodes != column.end(); ++columnNodes, ++z)
345 const SMDS_MeshNode* & node = *columnNodes;
346 if ( node ) continue; // skip bottom or top node
348 gp_XYZ coords = tBotNode.GetCoords();
349 trsf[z-1].Transforms( coords );
350 node = meshDS->AddNode( coords.X(), coords.Y(), coords.Z() );
351 meshDS->SetNodeInVolume( node, volumeID );
353 } // loop on bottom nodes
355 else // use block approach
357 // loop on nodes inside the bottom face
359 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.begin();
360 for ( ; bot_column != myBotToColumnMap.end(); ++bot_column )
362 const TNode& tBotNode = bot_column->first; // bottom TNode
363 if ( tBotNode.GetPositionType() != SMDS_TOP_FACE )
364 continue; // node is not inside face
366 // column nodes; middle part of the column are zero pointers
367 TNodeColumn& column = bot_column->second;
369 // compute bottom node parameters
370 gp_XYZ paramHint(-1,-1,-1);
371 if ( prevBNode.IsNeighbor( tBotNode ))
372 paramHint = prevBNode.GetParams();
373 if ( !myBlock.ComputeParameters( tBotNode.GetCoords(), tBotNode.ChangeParams(),
374 ID_BOT_FACE, paramHint ))
375 return error(TCom("Can't compute normalized parameters for node ")
376 << tBotNode.myNode->GetID() << " on the face #"
377 << myBlock.SubMesh( ID_BOT_FACE )->GetId() );
378 prevBNode = tBotNode;
380 myShapeXYZ[ ID_BOT_FACE ] = tBotNode.GetCoords();
381 gp_XYZ botParams = tBotNode.GetParams();
383 // compute top node parameters
384 myShapeXYZ[ ID_TOP_FACE ] = gpXYZ( column.back() );
385 gp_XYZ topParams = botParams;
387 if ( column.size() > 2 ) {
388 gp_Pnt topCoords = myShapeXYZ[ ID_TOP_FACE ];
389 if ( !myBlock.ComputeParameters( topCoords, topParams, ID_TOP_FACE, topParams ))
390 return error(TCom("Can't compute normalized parameters ")
391 << "for node " << column.back()->GetID()
392 << " on the face #"<< column.back()->GetPosition()->GetShapeId() );
396 TNodeColumn::iterator columnNodes = column.begin();
397 for ( int z = 0; columnNodes != column.end(); ++columnNodes, ++z)
399 const SMDS_MeshNode* & node = *columnNodes;
400 if ( node ) continue; // skip bottom or top node
402 // params of a node to create
403 double rz = (double) z / (double) ( column.size() - 1 );
404 gp_XYZ params = botParams * ( 1 - rz ) + topParams * rz;
406 // set coords on all faces and nodes
407 const int nbSideFaces = 4;
408 int sideFaceIDs[nbSideFaces] = { SMESH_Block::ID_Fx0z,
409 SMESH_Block::ID_Fx1z,
410 SMESH_Block::ID_F0yz,
411 SMESH_Block::ID_F1yz };
412 for ( int iF = 0; iF < nbSideFaces; ++iF )
413 if ( !setFaceAndEdgesXYZ( sideFaceIDs[ iF ], params, z ))
416 // compute coords for a new node
418 if ( !SMESH_Block::ShellPoint( params, myShapeXYZ, coords ))
419 return error("Can't compute coordinates by normalized parameters");
421 SHOWYXZ("TOPFacePoint ",myShapeXYZ[ ID_TOP_FACE]);
422 SHOWYXZ("BOT Node "<< tBotNode.myNode->GetID(),gpXYZ(tBotNode.myNode));
423 SHOWYXZ("ShellPoint ",coords);
426 node = meshDS->AddNode( coords.X(), coords.Y(), coords.Z() );
427 meshDS->SetNodeInVolume( node, volumeID );
429 } // loop on bottom nodes
434 SMESHDS_SubMesh* smDS = myBlock.SubMeshDS( ID_BOT_FACE );
435 if ( !smDS ) return error(COMPERR_BAD_INPUT_MESH, "Null submesh");
437 // loop on bottom mesh faces
438 SMDS_ElemIteratorPtr faceIt = smDS->GetElements();
439 while ( faceIt->more() )
441 const SMDS_MeshElement* face = faceIt->next();
442 if ( !face || face->GetType() != SMDSAbs_Face )
444 int nbNodes = face->NbNodes();
445 if ( face->IsQuadratic() )
448 // find node columns for each node
449 vector< const TNodeColumn* > columns( nbNodes );
450 for ( int i = 0; i < nbNodes; ++i )
452 const SMDS_MeshNode* n = face->GetNode( i );
453 if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ) {
454 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.find( n );
455 if ( bot_column == myBotToColumnMap.end() )
456 return error(TCom("No nodes found above node ") << n->GetID() );
457 columns[ i ] = & bot_column->second;
460 columns[ i ] = myBlock.GetNodeColumn( n );
462 return error(TCom("No side nodes found above node ") << n->GetID() );
466 AddPrisms( columns, myHelper );
468 } // loop on bottom mesh faces
474 //=======================================================================
475 //function : Evaluate
477 //=======================================================================
479 bool StdMeshers_Prism_3D::Evaluate(SMESH_Mesh& theMesh,
480 const TopoDS_Shape& theShape,
481 MapShapeNbElems& aResMap)
483 // find face contains only triangles
484 vector < SMESH_subMesh * >meshFaces;
485 TopTools_SequenceOfShape aFaces;
486 int NumBase = 0, i = 0, NbQFs = 0;
487 for (TopExp_Explorer exp(theShape, TopAbs_FACE); exp.More(); exp.Next()) {
489 aFaces.Append(exp.Current());
490 SMESH_subMesh *aSubMesh = theMesh.GetSubMesh(exp.Current());
491 meshFaces.push_back(aSubMesh);
492 MapShapeNbElemsItr anIt = aResMap.find(meshFaces[i-1]);
493 if( anIt==aResMap.end() ) {
494 SMESH_ComputeErrorPtr& smError = aSubMesh->GetComputeError();
495 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
498 std::vector<int> aVec = (*anIt).second;
499 int nbtri = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
500 int nbqua = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
501 if( nbtri==0 && nbqua>0 ) {
510 std::vector<int> aResVec(SMDSEntity_Last);
511 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
512 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
513 aResMap.insert(std::make_pair(sm,aResVec));
514 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
515 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
519 if(NumBase==0) NumBase = 1; // only quads => set 1 faces as base
521 // find number of 1d elems for base face
523 TopTools_MapOfShape Edges1;
524 for (TopExp_Explorer exp(aFaces.Value(NumBase), TopAbs_EDGE); exp.More(); exp.Next()) {
525 Edges1.Add(exp.Current());
526 SMESH_subMesh *sm = theMesh.GetSubMesh(exp.Current());
528 MapShapeNbElemsItr anIt = aResMap.find(sm);
529 if( anIt == aResMap.end() ) continue;
530 std::vector<int> aVec = (*anIt).second;
531 nb1d += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
534 // find face opposite to base face
536 for(i=1; i<=6; i++) {
537 if(i==NumBase) continue;
538 bool IsOpposite = true;
539 for(TopExp_Explorer exp(aFaces.Value(i), TopAbs_EDGE); exp.More(); exp.Next()) {
540 if( Edges1.Contains(exp.Current()) ) {
550 // find number of 2d elems on side faces
552 for(i=1; i<=6; i++) {
553 if( i==OppNum || i==NumBase ) continue;
554 MapShapeNbElemsItr anIt = aResMap.find( meshFaces[i-1] );
555 if( anIt == aResMap.end() ) continue;
556 std::vector<int> aVec = (*anIt).second;
557 nb2d += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
560 MapShapeNbElemsItr anIt = aResMap.find( meshFaces[NumBase-1] );
561 std::vector<int> aVec = (*anIt).second;
562 bool IsQuadratic = (aVec[SMDSEntity_Quad_Triangle]>aVec[SMDSEntity_Triangle]) ||
563 (aVec[SMDSEntity_Quad_Quadrangle]>aVec[SMDSEntity_Quadrangle]);
564 int nb2d_face0_3 = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
565 int nb2d_face0_4 = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
566 int nb0d_face0 = aVec[SMDSEntity_Node];
567 int nb1d_face0_int = ( nb2d_face0_3*3 + nb2d_face0_4*4 - nb1d ) / 2;
569 std::vector<int> aResVec(SMDSEntity_Last);
570 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
572 aResVec[SMDSEntity_Quad_Penta] = nb2d_face0_3 * ( nb2d/nb1d );
573 aResVec[SMDSEntity_Quad_Hexa] = nb2d_face0_4 * ( nb2d/nb1d );
574 aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
577 aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
578 aResVec[SMDSEntity_Penta] = nb2d_face0_3 * ( nb2d/nb1d );
579 aResVec[SMDSEntity_Hexa] = nb2d_face0_4 * ( nb2d/nb1d );
581 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
582 aResMap.insert(std::make_pair(sm,aResVec));
588 //================================================================================
590 * \brief Create prisms
591 * \param columns - columns of nodes generated from nodes of a mesh face
592 * \param helper - helper initialized by mesh and shape to add prisms to
594 //================================================================================
596 void StdMeshers_Prism_3D::AddPrisms( vector<const TNodeColumn*> & columns,
597 SMESH_MesherHelper* helper)
599 SMESHDS_Mesh * meshDS = helper->GetMeshDS();
600 int shapeID = helper->GetSubShapeID();
602 int nbNodes = columns.size();
603 int nbZ = columns[0]->size();
604 if ( nbZ < 2 ) return;
606 // find out orientation
607 bool isForward = true;
608 SMDS_VolumeTool vTool;
612 const SMDS_MeshNode* botNodes[3] = { (*columns[0])[z-1],
614 (*columns[2])[z-1] };
615 const SMDS_MeshNode* topNodes[3] = { (*columns[0])[z],
618 SMDS_VolumeOfNodes tmpVol ( botNodes[0], botNodes[1], botNodes[2],
619 topNodes[0], topNodes[1], topNodes[2]);
620 vTool.Set( &tmpVol );
621 isForward = vTool.IsForward();
625 const SMDS_MeshNode* botNodes[4] = { (*columns[0])[z-1], (*columns[1])[z-1],
626 (*columns[2])[z-1], (*columns[3])[z-1] };
627 const SMDS_MeshNode* topNodes[4] = { (*columns[0])[z], (*columns[1])[z],
628 (*columns[2])[z], (*columns[3])[z] };
629 SMDS_VolumeOfNodes tmpVol ( botNodes[0], botNodes[1], botNodes[2], botNodes[3],
630 topNodes[0], topNodes[1], topNodes[2], topNodes[3]);
631 vTool.Set( &tmpVol );
632 isForward = vTool.IsForward();
637 // vertical loop on columns
638 for ( z = 1; z < nbZ; ++z )
640 SMDS_MeshElement* vol = 0;
644 const SMDS_MeshNode* botNodes[3] = { (*columns[0])[z-1],
646 (*columns[2])[z-1] };
647 const SMDS_MeshNode* topNodes[3] = { (*columns[0])[z],
651 vol = helper->AddVolume( botNodes[0], botNodes[1], botNodes[2],
652 topNodes[0], topNodes[1], topNodes[2]);
654 vol = helper->AddVolume( topNodes[0], topNodes[1], topNodes[2],
655 botNodes[0], botNodes[1], botNodes[2]);
659 const SMDS_MeshNode* botNodes[4] = { (*columns[0])[z-1], (*columns[1])[z-1],
660 (*columns[2])[z-1], (*columns[3])[z-1] };
661 const SMDS_MeshNode* topNodes[4] = { (*columns[0])[z], (*columns[1])[z],
662 (*columns[2])[z], (*columns[3])[z] };
664 vol = helper->AddVolume( botNodes[0], botNodes[1], botNodes[2], botNodes[3],
665 topNodes[0], topNodes[1], topNodes[2], topNodes[3]);
667 vol = helper->AddVolume( topNodes[0], topNodes[1], topNodes[2], topNodes[3],
668 botNodes[0], botNodes[1], botNodes[2], botNodes[3]);
673 vector<const SMDS_MeshNode*> nodes( 2*nbNodes + 4*nbNodes);
674 vector<int> quantities( 2 + nbNodes, 4 );
675 quantities[0] = quantities[1] = nbNodes;
676 columns.resize( nbNodes + 1 );
677 columns[ nbNodes ] = columns[ 0 ];
678 for ( int i = 0; i < nbNodes; ++i ) {
679 nodes[ i ] = (*columns[ i ])[z-1]; // bottom
680 nodes[ i+nbNodes ] = (*columns[ i ])[z ]; // top
682 int di = 2*nbNodes + 4*i - 1;
683 nodes[ di ] = (*columns[i ])[z-1];
684 nodes[ di+1 ] = (*columns[i+1])[z-1];
685 nodes[ di+2 ] = (*columns[i+1])[z ];
686 nodes[ di+3 ] = (*columns[i ])[z ];
688 vol = meshDS->AddPolyhedralVolume( nodes, quantities );
690 if ( vol && shapeID > 0 )
691 meshDS->SetMeshElementOnShape( vol, shapeID );
695 //================================================================================
697 * \brief Find correspondence between bottom and top nodes
698 * If elements on the bottom and top faces are topologically different,
699 * and projection is possible and allowed, perform the projection
700 * \retval bool - is a success or not
702 //================================================================================
704 bool StdMeshers_Prism_3D::assocOrProjBottom2Top()
706 SMESH_subMesh * botSM = myBlock.SubMesh( ID_BOT_FACE );
707 SMESH_subMesh * topSM = myBlock.SubMesh( ID_TOP_FACE );
709 SMESHDS_SubMesh * botSMDS = botSM->GetSubMeshDS();
710 SMESHDS_SubMesh * topSMDS = topSM->GetSubMeshDS();
712 if ( !botSMDS || botSMDS->NbElements() == 0 )
713 return error(TCom("No elememts on face #") << botSM->GetId());
715 bool needProject = false;
717 botSMDS->NbElements() != topSMDS->NbElements() ||
718 botSMDS->NbNodes() != topSMDS->NbNodes())
720 if ( myBlock.HasNotQuadElemOnTop() )
721 return error(TCom("Mesh on faces #") << botSM->GetId()
722 <<" and #"<< topSM->GetId() << " seems different" );
726 if ( 0/*needProject && !myProjectTriangles*/ )
727 return error(TCom("Mesh on faces #") << botSM->GetId()
728 <<" and #"<< topSM->GetId() << " seems different" );
729 ///RETURN_BAD_RESULT("Need to project but not allowed");
733 return projectBottomToTop();
736 TopoDS_Face botFace = TopoDS::Face( myBlock.Shape( ID_BOT_FACE ));
737 TopoDS_Face topFace = TopoDS::Face( myBlock.Shape( ID_TOP_FACE ));
738 // associate top and bottom faces
739 TAssocTool::TShapeShapeMap shape2ShapeMap;
740 if ( !TAssocTool::FindSubShapeAssociation( botFace, myBlock.Mesh(),
741 topFace, myBlock.Mesh(),
743 return error(TCom("Topology of faces #") << botSM->GetId()
744 <<" and #"<< topSM->GetId() << " seems different" );
746 // Find matching nodes of top and bottom faces
748 if ( ! TAssocTool::FindMatchingNodesOnFaces( botFace, myBlock.Mesh(),
749 topFace, myBlock.Mesh(),
750 shape2ShapeMap, n2nMap ))
751 return error(TCom("Mesh on faces #") << botSM->GetId()
752 <<" and #"<< topSM->GetId() << " seems different" );
754 // Fill myBotToColumnMap
756 int zSize = myBlock.VerticalSize();
758 TNodeNodeMap::iterator bN_tN = n2nMap.begin();
759 for ( ; bN_tN != n2nMap.end(); ++bN_tN )
761 const SMDS_MeshNode* botNode = bN_tN->first;
762 const SMDS_MeshNode* topNode = bN_tN->second;
763 if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
764 continue; // wall columns are contained in myBlock
765 // compute bottom node params
767 // if ( zSize > 2 ) {
768 // gp_XYZ paramHint(-1,-1,-1);
769 // if ( prevTNode.IsNeighbor( bN ))
770 // paramHint = prevTNode.GetParams();
771 // if ( !myBlock.ComputeParameters( bN.GetCoords(), bN.ChangeParams(),
772 // ID_BOT_FACE, paramHint ))
773 // return error(TCom("Can't compute normalized parameters for node ")
774 // << botNode->GetID() << " on the face #"<< botSM->GetId() );
777 // create node column
778 TNode2ColumnMap::iterator bN_col =
779 myBotToColumnMap.insert( make_pair ( bN, TNodeColumn() )).first;
780 TNodeColumn & column = bN_col->second;
781 column.resize( zSize );
782 column.front() = botNode;
783 column.back() = topNode;
788 //================================================================================
790 * \brief Remove quadrangles from the top face and
791 * create triangles there by projection from the bottom
792 * \retval bool - a success or not
794 //================================================================================
796 bool StdMeshers_Prism_3D::projectBottomToTop()
798 SMESH_subMesh * botSM = myBlock.SubMesh( ID_BOT_FACE );
799 SMESH_subMesh * topSM = myBlock.SubMesh( ID_TOP_FACE );
801 SMESHDS_SubMesh * botSMDS = botSM->GetSubMeshDS();
802 SMESHDS_SubMesh * topSMDS = topSM->GetSubMeshDS();
805 topSM->ComputeStateEngine( SMESH_subMesh::CLEAN );
807 SMESHDS_Mesh* meshDS = myBlock.MeshDS();
808 int shapeID = myHelper->GetSubShapeID();
809 int topFaceID = meshDS->ShapeToIndex( topSM->GetSubShape() );
811 // Fill myBotToColumnMap
813 int zSize = myBlock.VerticalSize();
815 SMDS_NodeIteratorPtr nIt = botSMDS->GetNodes();
816 while ( nIt->more() )
818 const SMDS_MeshNode* botNode = nIt->next();
819 if ( botNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE )
821 // compute bottom node params
823 gp_XYZ paramHint(-1,-1,-1);
824 if ( prevTNode.IsNeighbor( bN ))
825 paramHint = prevTNode.GetParams();
826 if ( !myBlock.ComputeParameters( bN.GetCoords(), bN.ChangeParams(),
827 ID_BOT_FACE, paramHint ))
828 return error(TCom("Can't compute normalized parameters for node ")
829 << botNode->GetID() << " on the face #"<< botSM->GetId() );
831 // compute top node coords
832 gp_XYZ topXYZ; gp_XY topUV;
833 if ( !myBlock.FacePoint( ID_TOP_FACE, bN.GetParams(), topXYZ ) ||
834 !myBlock.FaceUV ( ID_TOP_FACE, bN.GetParams(), topUV ))
835 return error(TCom("Can't compute coordinates "
836 "by normalized parameters on the face #")<< topSM->GetId() );
837 SMDS_MeshNode * topNode = meshDS->AddNode( topXYZ.X(),topXYZ.Y(),topXYZ.Z() );
838 meshDS->SetNodeOnFace( topNode, topFaceID, topUV.X(), topUV.Y() );
839 // create node column
840 TNode2ColumnMap::iterator bN_col =
841 myBotToColumnMap.insert( make_pair ( bN, TNodeColumn() )).first;
842 TNodeColumn & column = bN_col->second;
843 column.resize( zSize );
844 column.front() = botNode;
845 column.back() = topNode;
850 // loop on bottom mesh faces
851 SMDS_ElemIteratorPtr faceIt = botSMDS->GetElements();
852 while ( faceIt->more() )
854 const SMDS_MeshElement* face = faceIt->next();
855 if ( !face || face->GetType() != SMDSAbs_Face )
857 int nbNodes = face->NbNodes();
858 if ( face->IsQuadratic() )
861 // find top node in columns for each bottom node
862 vector< const SMDS_MeshNode* > nodes( nbNodes );
863 for ( int i = 0; i < nbNodes; ++i )
865 const SMDS_MeshNode* n = face->GetNode( nbNodes - i - 1 );
866 if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ) {
867 TNode2ColumnMap::iterator bot_column = myBotToColumnMap.find( n );
868 if ( bot_column == myBotToColumnMap.end() )
869 return error(TCom("No nodes found above node ") << n->GetID() );
870 nodes[ i ] = bot_column->second.back();
873 const TNodeColumn* column = myBlock.GetNodeColumn( n );
875 return error(TCom("No side nodes found above node ") << n->GetID() );
876 nodes[ i ] = column->back();
879 // create a face, with reversed orientation
880 SMDS_MeshElement* newFace = 0;
884 newFace = myHelper->AddFace(nodes[0], nodes[1], nodes[2]);
888 newFace = myHelper->AddFace( nodes[0], nodes[1], nodes[2], nodes[3] );
892 newFace = meshDS->AddPolygonalFace( nodes );
894 if ( newFace && shapeID > 0 )
895 meshDS->SetMeshElementOnShape( newFace, shapeID );
901 //================================================================================
903 * \brief Set projection coordinates of a node to a face and it's subshapes
904 * \param faceID - the face given by in-block ID
905 * \param params - node normalized parameters
906 * \retval bool - is a success
908 //================================================================================
910 bool StdMeshers_Prism_3D::setFaceAndEdgesXYZ( const int faceID, const gp_XYZ& params, int z )
912 // find base and top edges of the face
913 enum { BASE = 0, TOP, LEFT, RIGHT };
914 vector< int > edgeVec; // 0-base, 1-top
915 SMESH_Block::GetFaceEdgesIDs( faceID, edgeVec );
917 myBlock.EdgePoint( edgeVec[ BASE ], params, myShapeXYZ[ edgeVec[ BASE ]]);
918 myBlock.EdgePoint( edgeVec[ TOP ], params, myShapeXYZ[ edgeVec[ TOP ]]);
920 SHOWYXZ("\nparams ", params);
921 SHOWYXZ("TOP is " <<edgeVec[ TOP ], myShapeXYZ[ edgeVec[ TOP]]);
922 SHOWYXZ("BASE is "<<edgeVec[ BASE], myShapeXYZ[ edgeVec[ BASE]]);
924 if ( faceID == SMESH_Block::ID_Fx0z || faceID == SMESH_Block::ID_Fx1z )
926 myBlock.EdgePoint( edgeVec[ LEFT ], params, myShapeXYZ[ edgeVec[ LEFT ]]);
927 myBlock.EdgePoint( edgeVec[ RIGHT ], params, myShapeXYZ[ edgeVec[ RIGHT ]]);
929 SHOWYXZ("VER "<<edgeVec[ LEFT], myShapeXYZ[ edgeVec[ LEFT]]);
930 SHOWYXZ("VER "<<edgeVec[ RIGHT], myShapeXYZ[ edgeVec[ RIGHT]]);
932 myBlock.FacePoint( faceID, params, myShapeXYZ[ faceID ]);
933 SHOWYXZ("FacePoint "<<faceID, myShapeXYZ[ faceID]);
938 //================================================================================
940 * \brief Return true if this node and other one belong to one face
942 //================================================================================
944 bool TNode::IsNeighbor( const TNode& other ) const
946 if ( !other.myNode || !myNode ) return false;
948 SMDS_ElemIteratorPtr fIt = other.myNode->GetInverseElementIterator(SMDSAbs_Face);
949 while ( fIt->more() )
950 if ( fIt->next()->GetNodeIndex( myNode ) >= 0 )
955 //================================================================================
957 * \brief Constructor. Initialization is needed
959 //================================================================================
961 StdMeshers_PrismAsBlock::StdMeshers_PrismAsBlock()
966 StdMeshers_PrismAsBlock::~StdMeshers_PrismAsBlock()
969 delete mySide; mySide = 0;
973 //================================================================================
975 * \brief Initialization.
976 * \param helper - helper loaded with mesh and 3D shape
977 * \param shape3D - a closed shell or solid
978 * \retval bool - false if a mesh or a shape are KO
980 //================================================================================
982 bool StdMeshers_PrismAsBlock::Init(SMESH_MesherHelper* helper,
983 const TopoDS_Shape& shape3D)
986 delete mySide; mySide = 0;
988 vector< TSideFace* > sideFaces( NB_WALL_FACES, 0 );
989 vector< pair< double, double> > params ( NB_WALL_FACES );
990 mySide = new TSideFace( sideFaces, params );
993 SMESHDS_Mesh* meshDS = myHelper->GetMeshDS();
996 myShapeIDMap.Clear();
997 myShapeIndex2ColumnMap.clear();
999 int wallFaceIds[ NB_WALL_FACES ] = { // to walk around a block
1000 SMESH_Block::ID_Fx0z, SMESH_Block::ID_F1yz,
1001 SMESH_Block::ID_Fx1z, SMESH_Block::ID_F0yz
1004 myError = SMESH_ComputeError::New();
1006 // -------------------------------------------------------------
1007 // Look for top and bottom faces: not quadrangle ones or meshed
1008 // with not quadrangle elements
1009 // -------------------------------------------------------------
1011 list< SMESH_subMesh* > notQuadGeomSubMesh;
1012 list< SMESH_subMesh* > notQuadElemSubMesh;
1015 SMESH_subMesh* mainSubMesh = myHelper->GetMesh()->GetSubMeshContaining( shape3D );
1016 if ( !mainSubMesh ) return error(COMPERR_BAD_INPUT_MESH,"Null submesh of shape3D");
1018 // analyse face submeshes
1019 SMESH_subMeshIteratorPtr smIt = mainSubMesh->getDependsOnIterator(false,false);
1020 while ( smIt->more() )
1022 SMESH_subMesh* sm = smIt->next();
1023 const TopoDS_Shape& face = sm->GetSubShape();
1024 if ( face.ShapeType() != TopAbs_FACE )
1028 // is quadrangle face?
1029 list< TopoDS_Edge > orderedEdges;
1030 list< int > nbEdgesInWires;
1032 int nbWires = GetOrderedEdges( TopoDS::Face( face ),
1033 V000, orderedEdges, nbEdgesInWires );
1034 if ( nbWires != 1 || nbEdgesInWires.front() != 4 )
1035 notQuadGeomSubMesh.push_back( sm );
1037 // look for not quadrangle mesh elements
1038 if ( SMESHDS_SubMesh* smDS = sm->GetSubMeshDS() ) {
1039 bool hasNotQuad = false;
1040 SMDS_ElemIteratorPtr eIt = smDS->GetElements();
1041 while ( eIt->more() && !hasNotQuad ) {
1042 const SMDS_MeshElement* elem = eIt->next();
1043 if ( elem->GetType() == SMDSAbs_Face ) {
1044 int nbNodes = elem->NbNodes();
1045 if ( elem->IsQuadratic() )
1047 hasNotQuad = ( nbNodes != 4 );
1051 notQuadElemSubMesh.push_back( sm );
1054 return error(COMPERR_BAD_INPUT_MESH,TCom("Not meshed face #")<<sm->GetId());
1056 // check if a quadrangle face is meshed with a quadranglar grid
1057 if ( notQuadGeomSubMesh.back() != sm &&
1058 notQuadElemSubMesh.back() != sm )
1060 // count nb edges on face sides
1061 vector< int > nbEdges;
1062 nbEdges.reserve( nbEdgesInWires.front() );
1063 for ( list< TopoDS_Edge >::iterator edge = orderedEdges.begin();
1064 edge != orderedEdges.end(); ++edge )
1066 if ( SMESHDS_SubMesh* smDS = meshDS->MeshElements( *edge ))
1067 nbEdges.push_back ( smDS->NbElements() );
1069 nbEdges.push_back ( 0 );
1071 int nbQuads = sm->GetSubMeshDS()->NbElements();
1072 if ( nbEdges[0] * nbEdges[1] != nbQuads ||
1073 nbEdges[0] != nbEdges[2] ||
1074 nbEdges[1] != nbEdges[3] )
1075 notQuadElemSubMesh.push_back( sm );
1079 // ----------------------------------------------------------------------
1080 // Analyse faces mesh and topology: choose the bottom submesh.
1081 // If there are not quadrangle geom faces, they are top and bottom ones.
1082 // Not quadrangle geom faces must be only on top and bottom.
1083 // ----------------------------------------------------------------------
1085 SMESH_subMesh * botSM = 0;
1086 SMESH_subMesh * topSM = 0;
1088 int nbNotQuad = notQuadGeomSubMesh.size();
1089 int nbNotQuadMeshed = notQuadElemSubMesh.size();
1090 bool hasNotQuad = ( nbNotQuad || nbNotQuadMeshed );
1093 if ( nbNotQuad > 0 && nbNotQuad != 2 )
1094 return error(COMPERR_BAD_SHAPE,
1095 TCom("More than 2 not quadrilateral faces: ")
1097 if ( nbNotQuadMeshed > 2 )
1098 return error(COMPERR_BAD_INPUT_MESH,
1099 TCom("More than 2 faces with not quadrangle elements: ")
1102 // get found submeshes
1105 if ( nbNotQuadMeshed > 0 ) botSM = notQuadElemSubMesh.front();
1106 else botSM = notQuadGeomSubMesh.front();
1107 if ( nbNotQuadMeshed > 1 ) topSM = notQuadElemSubMesh.back();
1108 else if ( nbNotQuad > 1 ) topSM = notQuadGeomSubMesh.back();
1110 // detect other bad cases
1111 if ( nbNotQuad == 2 && nbNotQuadMeshed > 0 ) {
1113 if ( nbNotQuadMeshed == 1 )
1114 ok = ( find( notQuadGeomSubMesh.begin(),
1115 notQuadGeomSubMesh.end(), botSM ) != notQuadGeomSubMesh.end() );
1117 ok = ( notQuadGeomSubMesh == notQuadElemSubMesh );
1119 return error(COMPERR_BAD_INPUT_MESH, "Side face meshed with not quadrangle elements");
1122 myNotQuadOnTop = ( nbNotQuadMeshed > 1 );
1124 // ----------------------------------------------------------
1126 if ( nbNotQuad == 0 ) // Standard block of 6 quadrangle faces ?
1128 // SMESH_Block will perform geometry analysis, we need just to find 2
1129 // connected vertices on top and bottom
1131 TopoDS_Vertex Vbot, Vtop;
1132 if ( nbNotQuadMeshed > 0 ) // Look for vertices
1134 TopTools_IndexedMapOfShape edgeMap;
1135 TopExp::MapShapes( botSM->GetSubShape(), TopAbs_EDGE, edgeMap );
1136 // vertex 1 is any vertex of the bottom face
1137 Vbot = TopExp::FirstVertex( TopoDS::Edge( edgeMap( 1 )));
1138 // vertex 2 is end vertex of edge sharing Vbot and not belonging to the bottom face
1139 TopTools_ListIteratorOfListOfShape ancestIt = Mesh()->GetAncestors( Vbot );
1140 for ( ; Vtop.IsNull() && ancestIt.More(); ancestIt.Next() )
1142 const TopoDS_Shape & ancestor = ancestIt.Value();
1143 if ( ancestor.ShapeType() == TopAbs_EDGE && !edgeMap.FindIndex( ancestor ))
1145 TopoDS_Vertex V1, V2;
1146 TopExp::Vertices( TopoDS::Edge( ancestor ), V1, V2);
1147 if ( Vbot.IsSame ( V1 )) Vtop = V2;
1148 else if ( Vbot.IsSame ( V2 )) Vtop = V1;
1149 // check that Vtop belongs to shape3D
1150 TopExp_Explorer exp( shape3D, TopAbs_VERTEX );
1151 for ( ; exp.More(); exp.Next() )
1152 if ( Vtop.IsSame( exp.Current() ))
1159 // get shell from shape3D
1161 TopExp_Explorer exp( shape3D, TopAbs_SHELL );
1163 for ( ; exp.More(); exp.Next(), ++nbShell )
1164 shell = TopoDS::Shell( exp.Current() );
1165 // if ( nbShell != 1 )
1166 // RETURN_BAD_RESULT("There must be 1 shell in the block");
1168 // Load geometry in SMESH_Block
1169 if ( !SMESH_Block::FindBlockShapes( shell, Vbot, Vtop, myShapeIDMap )) {
1171 return error(COMPERR_BAD_SHAPE, "Can't detect top and bottom of a prism");
1174 if ( !botSM ) botSM = Mesh()->GetSubMeshContaining( myShapeIDMap( ID_BOT_FACE ));
1175 if ( !topSM ) topSM = Mesh()->GetSubMeshContaining( myShapeIDMap( ID_TOP_FACE ));
1178 } // end Standard block of 6 quadrangle faces
1179 // --------------------------------------------------------
1181 // Here the top and bottom faces are found
1182 if ( nbNotQuadMeshed == 2 ) // roughly check correspondence of horiz meshes
1184 // SMESHDS_SubMesh* topSMDS = topSM->GetSubMeshDS();
1185 // SMESHDS_SubMesh* botSMDS = botSM->GetSubMeshDS();
1186 // if ( topSMDS->NbNodes() != botSMDS->NbNodes() ||
1187 // topSMDS->NbElements() != botSMDS->NbElements() )
1188 // RETURN_BAD_RESULT("Top mesh doesn't correspond to bottom one");
1191 // ---------------------------------------------------------
1192 // If there are not quadrangle geom faces, we emulate
1193 // a block of 6 quadrangle faces.
1194 // Load SMESH_Block with faces and edges geometry
1195 // ---------------------------------------------------------
1198 // find vertex 000 - the one with smallest coordinates (for easy DEBUG :-)
1200 double minVal = DBL_MAX, minX, val;
1201 for ( TopExp_Explorer exp( botSM->GetSubShape(), TopAbs_VERTEX );
1202 exp.More(); exp.Next() )
1204 const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
1205 gp_Pnt P = BRep_Tool::Pnt( v );
1206 val = P.X() + P.Y() + P.Z();
1207 if ( val < minVal || ( val == minVal && P.X() < minX )) {
1214 // Get ordered bottom edges
1215 list< TopoDS_Edge > orderedEdges;
1216 list< int > nbVertexInWires;
1217 SMESH_Block::GetOrderedEdges( TopoDS::Face( botSM->GetSubShape().Reversed() ),
1218 V000, orderedEdges, nbVertexInWires );
1219 // if ( nbVertexInWires.size() != 1 )
1220 // RETURN_BAD_RESULT("Wrong prism geometry");
1222 // Get Wall faces corresponding to the ordered bottom edges
1223 list< TopoDS_Face > wallFaces;
1224 if ( !GetWallFaces( Mesh(), shape3D, botSM->GetSubShape(), orderedEdges, wallFaces))
1225 return error(COMPERR_BAD_SHAPE, "Can't find side faces");
1227 // Find columns of wall nodes and calculate edges' lengths
1228 // --------------------------------------------------------
1230 myParam2ColumnMaps.clear();
1231 myParam2ColumnMaps.resize( orderedEdges.size() ); // total nb edges
1233 int iE, nbEdges = nbVertexInWires.front(); // nb outer edges
1234 vector< double > edgeLength( nbEdges );
1235 map< double, int > len2edgeMap;
1237 list< TopoDS_Edge >::iterator edgeIt = orderedEdges.begin();
1238 list< TopoDS_Face >::iterator faceIt = wallFaces.begin();
1239 for ( iE = 0; iE < nbEdges; ++edgeIt, ++faceIt )
1241 TParam2ColumnMap & faceColumns = myParam2ColumnMaps[ iE ];
1242 if ( !myHelper->LoadNodeColumns( faceColumns, *faceIt, *edgeIt, meshDS ))
1243 return error(COMPERR_BAD_INPUT_MESH, TCom("Can't find regular quadrangle mesh ")
1244 << "on a side face #" << MeshDS()->ShapeToIndex( *faceIt ));
1246 SHOWYXZ("\np1 F "<<iE, gpXYZ(faceColumns.begin()->second.front() ));
1247 SHOWYXZ("p2 F "<<iE, gpXYZ(faceColumns.rbegin()->second.front() ));
1248 SHOWYXZ("V First "<<iE, BRep_Tool::Pnt( TopExp::FirstVertex(*edgeIt,true )));
1250 edgeLength[ iE ] = SMESH_Algo::EdgeLength( *edgeIt );
1252 if ( nbEdges < NB_WALL_FACES ) // fill map used to split faces
1254 SMESHDS_SubMesh* smDS = meshDS->MeshElements( *edgeIt);
1256 return error(COMPERR_BAD_INPUT_MESH, TCom("Null submesh on the edge #")
1257 << MeshDS()->ShapeToIndex( *edgeIt ));
1258 // assure length uniqueness
1259 edgeLength[ iE ] *= smDS->NbNodes() + edgeLength[ iE ] / ( 1000 + iE );
1260 len2edgeMap[ edgeLength[ iE ]] = iE;
1264 // Load columns of internal edges (forming holes)
1265 // and fill map ShapeIndex to TParam2ColumnMap for them
1266 for ( ; edgeIt != orderedEdges.end() ; ++edgeIt, ++faceIt )
1268 TParam2ColumnMap & faceColumns = myParam2ColumnMaps[ iE ];
1269 if ( !myHelper->LoadNodeColumns( faceColumns, *faceIt, *edgeIt, meshDS ))
1270 return error(COMPERR_BAD_INPUT_MESH, TCom("Can't find regular quadrangle mesh ")
1271 << "on a side face #" << MeshDS()->ShapeToIndex( *faceIt ));
1273 int id = MeshDS()->ShapeToIndex( *edgeIt );
1274 bool isForward = true; // meaningless for intenal wires
1275 myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1276 // columns for vertices
1278 const SMDS_MeshNode* n0 = faceColumns.begin()->second.front();
1279 id = n0->GetPosition()->GetShapeId();
1280 myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1282 const SMDS_MeshNode* n1 = faceColumns.rbegin()->second.front();
1283 id = n1->GetPosition()->GetShapeId();
1284 myShapeIndex2ColumnMap[ id ] = make_pair( & faceColumns, isForward );
1285 // SHOWYXZ("\np1 F "<<iE, gpXYZ(faceColumns.begin()->second.front() ));
1286 // SHOWYXZ("p2 F "<<iE, gpXYZ(faceColumns.rbegin()->second.front() ));
1287 // SHOWYXZ("V First "<<iE, BRep_Tool::Pnt( TopExp::FirstVertex(*edgeIt,true )));
1291 // Create 4 wall faces of a block
1292 // -------------------------------
1294 if ( nbEdges <= NB_WALL_FACES ) // ************* Split faces if necessary
1296 map< int, int > iE2nbSplit;
1297 if ( nbEdges != NB_WALL_FACES ) // define how to split
1299 if ( len2edgeMap.size() != nbEdges )
1300 RETURN_BAD_RESULT("Uniqueness of edge lengths not assured");
1301 map< double, int >::reverse_iterator maxLen_i = len2edgeMap.rbegin();
1302 map< double, int >::reverse_iterator midLen_i = ++len2edgeMap.rbegin();
1303 double maxLen = maxLen_i->first;
1304 double midLen = ( len2edgeMap.size() == 1 ) ? 0 : midLen_i->first;
1305 switch ( nbEdges ) {
1306 case 1: // 0-th edge is split into 4 parts
1307 iE2nbSplit.insert( make_pair( 0, 4 )); break;
1308 case 2: // either the longest edge is split into 3 parts, or both edges into halves
1309 if ( maxLen / 3 > midLen / 2 ) {
1310 iE2nbSplit.insert( make_pair( maxLen_i->second, 3 ));
1313 iE2nbSplit.insert( make_pair( maxLen_i->second, 2 ));
1314 iE2nbSplit.insert( make_pair( midLen_i->second, 2 ));
1318 // split longest into halves
1319 iE2nbSplit.insert( make_pair( maxLen_i->second, 2 ));
1322 // Create TSideFace's
1323 faceIt = wallFaces.begin();
1324 edgeIt = orderedEdges.begin();
1326 for ( iE = 0; iE < nbEdges; ++edgeIt, ++faceIt )
1329 map< int, int >::iterator i_nb = iE2nbSplit.find( iE );
1330 if ( i_nb != iE2nbSplit.end() ) {
1332 int nbSplit = i_nb->second;
1333 vector< double > params;
1334 splitParams( nbSplit, &myParam2ColumnMaps[ iE ], params );
1335 bool isForward = ( edgeIt->Orientation() == TopAbs_FORWARD );
1336 for ( int i = 0; i < nbSplit; ++i ) {
1337 double f = ( isForward ? params[ i ] : params[ nbSplit - i-1 ]);
1338 double l = ( isForward ? params[ i+1 ] : params[ nbSplit - i ]);
1339 TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1341 &myParam2ColumnMaps[ iE ], f, l );
1342 mySide->SetComponent( iSide++, comp );
1346 TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1348 &myParam2ColumnMaps[ iE ]);
1349 mySide->SetComponent( iSide++, comp );
1354 else { // **************************** Unite faces
1356 // unite first faces
1357 int nbExraFaces = nbEdges - 3;
1359 double u0 = 0, sumLen = 0;
1360 for ( iE = 0; iE < nbExraFaces; ++iE )
1361 sumLen += edgeLength[ iE ];
1363 vector< TSideFace* > components( nbExraFaces );
1364 vector< pair< double, double> > params( nbExraFaces );
1365 faceIt = wallFaces.begin();
1366 edgeIt = orderedEdges.begin();
1367 for ( iE = 0; iE < nbExraFaces; ++edgeIt, ++faceIt )
1369 components[ iE ] = new TSideFace( myHelper, wallFaceIds[ iSide ],
1371 &myParam2ColumnMaps[ iE ]);
1372 double u1 = u0 + edgeLength[ iE ] / sumLen;
1373 params[ iE ] = make_pair( u0 , u1 );
1377 mySide->SetComponent( iSide++, new TSideFace( components, params ));
1379 // fill the rest faces
1380 for ( ; iE < nbEdges; ++faceIt, ++edgeIt )
1382 TSideFace* comp = new TSideFace( myHelper, wallFaceIds[ iSide ],
1384 &myParam2ColumnMaps[ iE ]);
1385 mySide->SetComponent( iSide++, comp );
1391 // Fill geometry fields of SMESH_Block
1392 // ------------------------------------
1394 TopoDS_Face botF = TopoDS::Face( botSM->GetSubShape() );
1395 TopoDS_Face topF = TopoDS::Face( topSM->GetSubShape() );
1397 vector< int > botEdgeIdVec;
1398 SMESH_Block::GetFaceEdgesIDs( ID_BOT_FACE, botEdgeIdVec );
1400 bool isForward[NB_WALL_FACES] = { true, true, true, true };
1401 Adaptor2d_Curve2d* botPcurves[NB_WALL_FACES];
1402 Adaptor2d_Curve2d* topPcurves[NB_WALL_FACES];
1404 for ( int iF = 0; iF < NB_WALL_FACES; ++iF )
1406 TSideFace * sideFace = mySide->GetComponent( iF );
1408 RETURN_BAD_RESULT("NULL TSideFace");
1409 int fID = sideFace->FaceID();
1411 // fill myShapeIDMap
1412 if ( sideFace->InsertSubShapes( myShapeIDMap ) != 8 &&
1413 !sideFace->IsComplex())
1414 MESSAGE( ": Warning : InsertSubShapes() < 8 on side " << iF );
1416 // side faces geometry
1417 Adaptor2d_Curve2d* pcurves[NB_WALL_FACES];
1418 if ( !sideFace->GetPCurves( pcurves ))
1419 RETURN_BAD_RESULT("TSideFace::GetPCurves() failed");
1421 SMESH_Block::TFace& tFace = myFace[ fID - ID_FirstF ];
1422 tFace.Set( fID, sideFace->Surface(), pcurves, isForward );
1424 SHOWYXZ( endl<<"F "<< iF << " id " << fID << " FRW " << sideFace->IsForward(), sideFace->Value(0,0));
1425 // edges 3D geometry
1426 vector< int > edgeIdVec;
1427 SMESH_Block::GetFaceEdgesIDs( fID, edgeIdVec );
1428 for ( int isMax = 0; isMax < 2; ++isMax ) {
1430 int eID = edgeIdVec[ isMax ];
1431 SMESH_Block::TEdge& tEdge = myEdge[ eID - ID_FirstE ];
1432 tEdge.Set( eID, sideFace->HorizCurve(isMax), true);
1433 SHOWYXZ(eID<<" HOR"<<isMax<<"(0)", sideFace->HorizCurve(isMax)->Value(0));
1434 SHOWYXZ(eID<<" HOR"<<isMax<<"(1)", sideFace->HorizCurve(isMax)->Value(1));
1437 int eID = edgeIdVec[ isMax+2 ];
1438 SMESH_Block::TEdge& tEdge = myEdge[ eID - ID_FirstE ];
1439 tEdge.Set( eID, sideFace->VertiCurve(isMax), true);
1440 SHOWYXZ(eID<<" VER"<<isMax<<"(0)", sideFace->VertiCurve(isMax)->Value(0));
1441 SHOWYXZ(eID<<" VER"<<isMax<<"(1)", sideFace->VertiCurve(isMax)->Value(1));
1444 vector< int > vertexIdVec;
1445 SMESH_Block::GetEdgeVertexIDs( eID, vertexIdVec );
1446 myPnt[ vertexIdVec[0] - ID_FirstV ] = tEdge.GetCurve()->Value(0).XYZ();
1447 myPnt[ vertexIdVec[1] - ID_FirstV ] = tEdge.GetCurve()->Value(1).XYZ();
1450 // pcurves on horizontal faces
1451 for ( iE = 0; iE < NB_WALL_FACES; ++iE ) {
1452 if ( edgeIdVec[ BOTTOM_EDGE ] == botEdgeIdVec[ iE ] ) {
1453 botPcurves[ iE ] = sideFace->HorizPCurve( false, botF );
1454 topPcurves[ iE ] = sideFace->HorizPCurve( true, topF );
1458 //sideFace->dumpNodes( 4 ); // debug
1460 // horizontal faces geometry
1462 SMESH_Block::TFace& tFace = myFace[ ID_BOT_FACE - ID_FirstF ];
1463 tFace.Set( ID_BOT_FACE, new BRepAdaptor_Surface( botF ), botPcurves, isForward );
1464 SMESH_Block::Insert( botF, ID_BOT_FACE, myShapeIDMap );
1467 SMESH_Block::TFace& tFace = myFace[ ID_TOP_FACE - ID_FirstF ];
1468 tFace.Set( ID_TOP_FACE, new BRepAdaptor_Surface( topF ), topPcurves, isForward );
1469 SMESH_Block::Insert( topF, ID_TOP_FACE, myShapeIDMap );
1472 // Fill map ShapeIndex to TParam2ColumnMap
1473 // ----------------------------------------
1475 list< TSideFace* > fList;
1476 list< TSideFace* >::iterator fListIt;
1477 fList.push_back( mySide );
1478 for ( fListIt = fList.begin(); fListIt != fList.end(); ++fListIt)
1480 int nb = (*fListIt)->NbComponents();
1481 for ( int i = 0; i < nb; ++i ) {
1482 if ( TSideFace* comp = (*fListIt)->GetComponent( i ))
1483 fList.push_back( comp );
1485 if ( TParam2ColumnMap* cols = (*fListIt)->GetColumns()) {
1486 // columns for a base edge
1487 int id = MeshDS()->ShapeToIndex( (*fListIt)->BaseEdge() );
1488 bool isForward = (*fListIt)->IsForward();
1489 myShapeIndex2ColumnMap[ id ] = make_pair( cols, isForward );
1491 // columns for vertices
1492 const SMDS_MeshNode* n0 = cols->begin()->second.front();
1493 id = n0->GetPosition()->GetShapeId();
1494 myShapeIndex2ColumnMap[ id ] = make_pair( cols, isForward );
1496 const SMDS_MeshNode* n1 = cols->rbegin()->second.front();
1497 id = n1->GetPosition()->GetShapeId();
1498 myShapeIndex2ColumnMap[ id ] = make_pair( cols, !isForward );
1502 // gp_XYZ testPar(0.25, 0.25, 0), testCoord;
1503 // if ( !FacePoint( ID_BOT_FACE, testPar, testCoord ))
1504 // RETURN_BAD_RESULT("TEST FacePoint() FAILED");
1505 // SHOWYXZ("IN TEST PARAM" , testPar);
1506 // SHOWYXZ("OUT TEST CORD" , testCoord);
1507 // if ( !ComputeParameters( testCoord, testPar , ID_BOT_FACE))
1508 // RETURN_BAD_RESULT("TEST ComputeParameters() FAILED");
1509 // SHOWYXZ("OUT TEST PARAM" , testPar);
1514 //================================================================================
1516 * \brief Return pointer to column of nodes
1517 * \param node - bottom node from which the returned column goes up
1518 * \retval const TNodeColumn* - the found column
1520 //================================================================================
1522 const TNodeColumn* StdMeshers_PrismAsBlock::GetNodeColumn(const SMDS_MeshNode* node) const
1524 int sID = node->GetPosition()->GetShapeId();
1526 map<int, pair< TParam2ColumnMap*, bool > >::const_iterator col_frw =
1527 myShapeIndex2ColumnMap.find( sID );
1528 if ( col_frw != myShapeIndex2ColumnMap.end() ) {
1529 const TParam2ColumnMap* cols = col_frw->second.first;
1530 TParam2ColumnIt u_col = cols->begin();
1531 for ( ; u_col != cols->end(); ++u_col )
1532 if ( u_col->second[ 0 ] == node )
1533 return & u_col->second;
1538 //=======================================================================
1539 //function : GetLayersTransformation
1540 //purpose : Return transformations to get coordinates of nodes of each layer
1541 // by nodes of the bottom. Layer is a set of nodes at a certain step
1542 // from bottom to top.
1543 //=======================================================================
1545 bool StdMeshers_PrismAsBlock::GetLayersTransformation(vector<gp_Trsf> & trsf) const
1547 const int zSize = VerticalSize();
1548 if ( zSize < 3 ) return true;
1549 trsf.resize( zSize - 2 );
1551 // Select some node columns by which we will define coordinate system of layers
1553 vector< const TNodeColumn* > columns;
1555 const TopoDS_Shape& baseFace = Shape(ID_BOT_FACE);
1556 list< TopoDS_Edge > orderedEdges;
1557 list< int > nbEdgesInWires;
1558 GetOrderedEdges( TopoDS::Face( baseFace ), TopoDS_Vertex(), orderedEdges, nbEdgesInWires );
1560 list< TopoDS_Edge >::iterator edgeIt = orderedEdges.begin();
1561 for ( int iE = 0; iE < nbEdgesInWires.front(); ++iE, ++edgeIt )
1563 const TParam2ColumnMap& u2colMap =
1564 GetParam2ColumnMap( myHelper->GetMeshDS()->ShapeToIndex( *edgeIt ), isReverse );
1565 isReverse = ( edgeIt->Orientation() == TopAbs_REVERSED );
1566 double f = u2colMap.begin()->first, l = u2colMap.rbegin()->first;
1567 if ( isReverse ) swap ( f, l );
1568 const int nbCol = 5;
1569 for ( int i = 0; i < nbCol; ++i )
1571 double u = f + i/double(nbCol) * ( l - f );
1572 const TNodeColumn* col = & getColumn( & u2colMap, u )->second;
1573 if ( columns.empty() || col != columns.back() )
1574 columns.push_back( col );
1579 // Find tolerance to check transformations
1584 for ( int i = 0; i < columns.size(); ++i )
1585 bndBox.Add( gpXYZ( columns[i]->front() ));
1586 tol2 = bndBox.SquareExtent() * 4 * 1e-4;
1589 // Compute transformations
1592 gp_Trsf fromCsZ, toCs0;
1593 gp_Ax3 cs0 = getLayerCoordSys(0, columns, xCol );
1594 //double dist0 = cs0.Location().Distance( gpXYZ( (*columns[0])[0]));
1595 toCs0.SetTransformation( cs0 );
1596 for ( int z = 1; z < zSize-1; ++z )
1598 gp_Ax3 csZ = getLayerCoordSys(z, columns, xCol );
1599 //double distZ = csZ.Location().Distance( gpXYZ( (*columns[0])[z]));
1600 fromCsZ.SetTransformation( csZ );
1602 gp_Trsf& t = trsf[ z-1 ];
1603 t = fromCsZ * toCs0;
1604 //t.SetScaleFactor( distZ/dist0 ); - it does not work properly, wrong base point
1606 // check a transformation
1607 for ( int i = 0; i < columns.size(); ++i )
1609 gp_Pnt p0 = gpXYZ( (*columns[i])[0] );
1610 gp_Pnt pz = gpXYZ( (*columns[i])[z] );
1611 t.Transforms( p0.ChangeCoord() );
1612 if ( p0.SquareDistance( pz ) > tol2 )
1619 //================================================================================
1621 * \brief Check curve orientation of a bootom edge
1622 * \param meshDS - mesh DS
1623 * \param columnsMap - node columns map of side face
1624 * \param bottomEdge - the bootom edge
1625 * \param sideFaceID - side face in-block ID
1626 * \retval bool - true if orientation coinside with in-block froward orientation
1628 //================================================================================
1630 bool StdMeshers_PrismAsBlock::IsForwardEdge(SMESHDS_Mesh* meshDS,
1631 const TParam2ColumnMap& columnsMap,
1632 const TopoDS_Edge & bottomEdge,
1633 const int sideFaceID)
1635 bool isForward = false;
1636 if ( TAssocTool::IsClosedEdge( bottomEdge ))
1638 isForward = ( bottomEdge.Orientation() == TopAbs_FORWARD );
1642 const TNodeColumn& firstCol = columnsMap.begin()->second;
1643 const SMDS_MeshNode* bottomNode = firstCol[0];
1644 TopoDS_Shape firstVertex = SMESH_MesherHelper::GetSubShapeByNode( bottomNode, meshDS );
1645 isForward = ( firstVertex.IsSame( TopExp::FirstVertex( bottomEdge, true )));
1647 // on 2 of 4 sides first vertex is end
1648 if ( sideFaceID == ID_Fx1z || sideFaceID == ID_F0yz )
1649 isForward = !isForward;
1653 //================================================================================
1655 * \brief Find wall faces by bottom edges
1656 * \param mesh - the mesh
1657 * \param mainShape - the prism
1658 * \param bottomFace - the bottom face
1659 * \param bottomEdges - edges bounding the bottom face
1660 * \param wallFaces - faces list to fill in
1662 //================================================================================
1664 bool StdMeshers_PrismAsBlock::GetWallFaces( SMESH_Mesh* mesh,
1665 const TopoDS_Shape & mainShape,
1666 const TopoDS_Shape & bottomFace,
1667 const std::list< TopoDS_Edge >& bottomEdges,
1668 std::list< TopoDS_Face >& wallFaces)
1672 TopTools_IndexedMapOfShape faceMap;
1673 TopExp::MapShapes( mainShape, TopAbs_FACE, faceMap );
1675 list< TopoDS_Edge >::const_iterator edge = bottomEdges.begin();
1676 for ( ; edge != bottomEdges.end(); ++edge )
1678 TopTools_ListIteratorOfListOfShape ancestIt = mesh->GetAncestors( *edge );
1679 for ( ; ancestIt.More(); ancestIt.Next() )
1681 const TopoDS_Shape& ancestor = ancestIt.Value();
1682 if ( ancestor.ShapeType() == TopAbs_FACE && // face
1683 !bottomFace.IsSame( ancestor ) && // not bottom
1684 faceMap.FindIndex( ancestor )) // belongs to the prism
1686 wallFaces.push_back( TopoDS::Face( ancestor ));
1691 return ( wallFaces.size() == bottomEdges.size() );
1694 //================================================================================
1696 * \brief Constructor
1697 * \param faceID - in-block ID
1698 * \param face - geom face
1699 * \param columnsMap - map of node columns
1700 * \param first - first normalized param
1701 * \param last - last normalized param
1703 //================================================================================
1705 StdMeshers_PrismAsBlock::TSideFace::TSideFace(SMESH_MesherHelper* helper,
1707 const TopoDS_Face& face,
1708 const TopoDS_Edge& baseEdge,
1709 TParam2ColumnMap* columnsMap,
1713 myParamToColumnMap( columnsMap ),
1714 myBaseEdge( baseEdge ),
1717 mySurface.Initialize( face );
1718 myParams.resize( 1 );
1719 myParams[ 0 ] = make_pair( first, last );
1720 myIsForward = StdMeshers_PrismAsBlock::IsForwardEdge( myHelper->GetMeshDS(),
1721 *myParamToColumnMap,
1725 //================================================================================
1727 * \brief Constructor of complex side face
1729 //================================================================================
1731 StdMeshers_PrismAsBlock::TSideFace::
1732 TSideFace(const vector< TSideFace* >& components,
1733 const vector< pair< double, double> > & params)
1734 :myID( components[0] ? components[0]->myID : 0 ),
1735 myParamToColumnMap( 0 ),
1737 myIsForward( true ),
1738 myComponents( components ),
1739 myHelper( components[0] ? components[0]->myHelper : 0 )
1741 //================================================================================
1743 * \brief Copy constructor
1744 * \param other - other side
1746 //================================================================================
1748 StdMeshers_PrismAsBlock::TSideFace::TSideFace( const TSideFace& other )
1751 mySurface = other.mySurface;
1752 myBaseEdge = other.myBaseEdge;
1753 myParams = other.myParams;
1754 myIsForward = other.myIsForward;
1755 myHelper = other.myHelper;
1756 myParamToColumnMap = other.myParamToColumnMap;
1758 myComponents.resize( other.myComponents.size());
1759 for (int i = 0 ; i < myComponents.size(); ++i )
1760 myComponents[ i ] = new TSideFace( *other.myComponents[ i ]);
1763 //================================================================================
1765 * \brief Deletes myComponents
1767 //================================================================================
1769 StdMeshers_PrismAsBlock::TSideFace::~TSideFace()
1771 for (int i = 0 ; i < myComponents.size(); ++i )
1772 if ( myComponents[ i ] )
1773 delete myComponents[ i ];
1776 //================================================================================
1778 * \brief Return geometry of the vertical curve
1779 * \param isMax - true means curve located closer to (1,1,1) block point
1780 * \retval Adaptor3d_Curve* - curve adaptor
1782 //================================================================================
1784 Adaptor3d_Curve* StdMeshers_PrismAsBlock::TSideFace::VertiCurve(const bool isMax) const
1786 if ( !myComponents.empty() ) {
1788 return myComponents.back()->VertiCurve(isMax);
1790 return myComponents.front()->VertiCurve(isMax);
1792 double f = myParams[0].first, l = myParams[0].second;
1793 if ( !myIsForward ) std::swap( f, l );
1794 return new TVerticalEdgeAdaptor( myParamToColumnMap, isMax ? l : f );
1797 //================================================================================
1799 * \brief Return geometry of the top or bottom curve
1801 * \retval Adaptor3d_Curve* -
1803 //================================================================================
1805 Adaptor3d_Curve* StdMeshers_PrismAsBlock::TSideFace::HorizCurve(const bool isTop) const
1807 return new THorizontalEdgeAdaptor( this, isTop );
1810 //================================================================================
1812 * \brief Return pcurves
1813 * \param pcurv - array of 4 pcurves
1814 * \retval bool - is a success
1816 //================================================================================
1818 bool StdMeshers_PrismAsBlock::TSideFace::GetPCurves(Adaptor2d_Curve2d* pcurv[4]) const
1820 int iEdge[ 4 ] = { BOTTOM_EDGE, TOP_EDGE, V0_EDGE, V1_EDGE };
1822 for ( int i = 0 ; i < 4 ; ++i ) {
1823 Handle(Geom2d_Line) line;
1824 switch ( iEdge[ i ] ) {
1826 line = new Geom2d_Line( gp_Pnt2d( 0, 1 ), gp::DX2d() ); break;
1828 line = new Geom2d_Line( gp::Origin2d(), gp::DX2d() ); break;
1830 line = new Geom2d_Line( gp::Origin2d(), gp::DY2d() ); break;
1832 line = new Geom2d_Line( gp_Pnt2d( 1, 0 ), gp::DY2d() ); break;
1834 pcurv[ i ] = new Geom2dAdaptor_Curve( line, 0, 1 );
1839 //================================================================================
1841 * \brief Returns geometry of pcurve on a horizontal face
1842 * \param isTop - is top or bottom face
1843 * \param horFace - a horizontal face
1844 * \retval Adaptor2d_Curve2d* - curve adaptor
1846 //================================================================================
1849 StdMeshers_PrismAsBlock::TSideFace::HorizPCurve(const bool isTop,
1850 const TopoDS_Face& horFace) const
1852 return new TPCurveOnHorFaceAdaptor( this, isTop, horFace );
1855 //================================================================================
1857 * \brief Return a component corresponding to parameter
1858 * \param U - parameter along a horizontal size
1859 * \param localU - parameter along a horizontal size of a component
1860 * \retval TSideFace* - found component
1862 //================================================================================
1864 StdMeshers_PrismAsBlock::TSideFace*
1865 StdMeshers_PrismAsBlock::TSideFace::GetComponent(const double U,double & localU) const
1868 if ( myComponents.empty() )
1869 return const_cast<TSideFace*>( this );
1872 for ( i = 0; i < myComponents.size(); ++i )
1873 if ( U < myParams[ i ].second )
1875 if ( i >= myComponents.size() )
1876 i = myComponents.size() - 1;
1878 double f = myParams[ i ].first, l = myParams[ i ].second;
1879 localU = ( U - f ) / ( l - f );
1880 return myComponents[ i ];
1883 //================================================================================
1885 * \brief Find node columns for a parameter
1886 * \param U - parameter along a horizontal edge
1887 * \param col1 - the 1st found column
1888 * \param col2 - the 2nd found column
1889 * \retval r - normalized position of U between the found columns
1891 //================================================================================
1893 double StdMeshers_PrismAsBlock::TSideFace::GetColumns(const double U,
1894 TParam2ColumnIt & col1,
1895 TParam2ColumnIt & col2) const
1897 double u = U, r = 0;
1898 if ( !myComponents.empty() ) {
1899 TSideFace * comp = GetComponent(U,u);
1900 return comp->GetColumns( u, col1, col2 );
1905 double f = myParams[0].first, l = myParams[0].second;
1906 u = f + u * ( l - f );
1908 col1 = col2 = getColumn( myParamToColumnMap, u );
1909 if ( ++col2 == myParamToColumnMap->end() ) {
1914 double uf = col1->first;
1915 double ul = col2->first;
1916 r = ( u - uf ) / ( ul - uf );
1921 //================================================================================
1923 * \brief Return coordinates by normalized params
1924 * \param U - horizontal param
1925 * \param V - vertical param
1926 * \retval gp_Pnt - result point
1928 //================================================================================
1930 gp_Pnt StdMeshers_PrismAsBlock::TSideFace::Value(const Standard_Real U,
1931 const Standard_Real V) const
1933 if ( !myComponents.empty() ) {
1935 TSideFace * comp = GetComponent(U,u);
1936 return comp->Value( u, V );
1939 TParam2ColumnIt u_col1, u_col2;
1940 double vR, hR = GetColumns( U, u_col1, u_col2 );
1942 const SMDS_MeshNode* n1 = 0;
1943 const SMDS_MeshNode* n2 = 0;
1944 const SMDS_MeshNode* n3 = 0;
1945 const SMDS_MeshNode* n4 = 0;
1947 // BEGIN issue 0020680: EDF 1252 SMESH: Bad cell created by Radial prism in center of torus
1948 // Workaround for a wrongly located point returned by mySurface.Value() for
1949 // UV located near boundary of BSpline surface.
1950 // To bypass the problem, we take point from 3D curve of edge.
1951 // It solves pb of the bloc_fiss_new.py
1952 const double tol = 1e-3;
1953 if ( V < tol || V+tol >= 1. )
1955 n1 = V < tol ? u_col1->second.front() : u_col1->second.back();
1956 n3 = V < tol ? u_col2->second.front() : u_col2->second.back();
1964 TopoDS_Shape s = myHelper->GetSubShapeByNode( n1, myHelper->GetMeshDS() );
1965 if ( s.ShapeType() != TopAbs_EDGE )
1966 s = myHelper->GetSubShapeByNode( n3, myHelper->GetMeshDS() );
1967 if ( s.ShapeType() == TopAbs_EDGE )
1968 edge = TopoDS::Edge( s );
1970 if ( !edge.IsNull() )
1972 double u1 = myHelper->GetNodeU( edge, n1 );
1973 double u3 = myHelper->GetNodeU( edge, n3 );
1974 double u = u1 * ( 1 - hR ) + u3 * hR;
1975 TopLoc_Location loc; double f,l;
1976 Handle(Geom_Curve) curve = BRep_Tool::Curve( edge,loc,f,l );
1977 return curve->Value( u ).Transformed( loc );
1980 // END issue 0020680: EDF 1252 SMESH: Bad cell created by Radial prism in center of torus
1982 vR = getRAndNodes( & u_col1->second, V, n1, n2 );
1983 vR = getRAndNodes( & u_col2->second, V, n3, n4 );
1985 gp_XY uv1 = myHelper->GetNodeUV( mySurface.Face(), n1, n4);
1986 gp_XY uv2 = myHelper->GetNodeUV( mySurface.Face(), n2, n3);
1987 gp_XY uv12 = uv1 * ( 1 - vR ) + uv2 * vR;
1989 gp_XY uv3 = myHelper->GetNodeUV( mySurface.Face(), n3, n2);
1990 gp_XY uv4 = myHelper->GetNodeUV( mySurface.Face(), n4, n1);
1991 gp_XY uv34 = uv3 * ( 1 - vR ) + uv4 * vR;
1993 gp_XY uv = uv12 * ( 1 - hR ) + uv34 * hR;
1995 gp_Pnt p = mySurface.Value( uv.X(), uv.Y() );
2000 //================================================================================
2002 * \brief Return boundary edge
2003 * \param edge - edge index
2004 * \retval TopoDS_Edge - found edge
2006 //================================================================================
2008 TopoDS_Edge StdMeshers_PrismAsBlock::TSideFace::GetEdge(const int iEdge) const
2010 if ( !myComponents.empty() ) {
2012 case V0_EDGE : return myComponents.front()->GetEdge( iEdge );
2013 case V1_EDGE : return myComponents.back() ->GetEdge( iEdge );
2014 default: return TopoDS_Edge();
2018 const SMDS_MeshNode* node = 0;
2019 SMESHDS_Mesh * meshDS = myHelper->GetMesh()->GetMeshDS();
2020 TNodeColumn* column;
2025 column = & (( ++myParamToColumnMap->begin())->second );
2026 node = ( iEdge == TOP_EDGE ) ? column->back() : column->front();
2027 edge = myHelper->GetSubShapeByNode ( node, meshDS );
2028 if ( edge.ShapeType() == TopAbs_VERTEX ) {
2029 column = & ( myParamToColumnMap->begin()->second );
2030 node = ( iEdge == TOP_EDGE ) ? column->back() : column->front();
2035 bool back = ( iEdge == V1_EDGE );
2036 if ( !myIsForward ) back = !back;
2038 column = & ( myParamToColumnMap->rbegin()->second );
2040 column = & ( myParamToColumnMap->begin()->second );
2041 if ( column->size() > 0 )
2042 edge = myHelper->GetSubShapeByNode( (*column)[ 1 ], meshDS );
2043 if ( edge.IsNull() || edge.ShapeType() == TopAbs_VERTEX )
2044 node = column->front();
2049 if ( !edge.IsNull() && edge.ShapeType() == TopAbs_EDGE )
2050 return TopoDS::Edge( edge );
2052 // find edge by 2 vertices
2053 TopoDS_Shape V1 = edge;
2054 TopoDS_Shape V2 = myHelper->GetSubShapeByNode( node, meshDS );
2055 if ( V2.ShapeType() == TopAbs_VERTEX && !V2.IsSame( V1 ))
2057 TopTools_ListIteratorOfListOfShape ancestIt =
2058 myHelper->GetMesh()->GetAncestors( V1 );
2059 for ( ; ancestIt.More(); ancestIt.Next() )
2061 const TopoDS_Shape & ancestor = ancestIt.Value();
2062 if ( ancestor.ShapeType() == TopAbs_EDGE )
2063 for ( TopExp_Explorer e( ancestor, TopAbs_VERTEX ); e.More(); e.Next() )
2064 if ( V2.IsSame( e.Current() ))
2065 return TopoDS::Edge( ancestor );
2068 return TopoDS_Edge();
2071 //================================================================================
2073 * \brief Fill block subshapes
2074 * \param shapeMap - map to fill in
2075 * \retval int - nb inserted subshapes
2077 //================================================================================
2079 int StdMeshers_PrismAsBlock::TSideFace::InsertSubShapes(TBlockShapes& shapeMap) const
2084 vector< int > edgeIdVec;
2085 SMESH_Block::GetFaceEdgesIDs( myID, edgeIdVec );
2087 for ( int i = BOTTOM_EDGE; i <=V1_EDGE ; ++i ) {
2088 TopoDS_Edge e = GetEdge( i );
2089 if ( !e.IsNull() ) {
2090 nbInserted += SMESH_Block::Insert( e, edgeIdVec[ i ], shapeMap);
2094 // Insert corner vertices
2096 TParam2ColumnIt col1, col2 ;
2097 vector< int > vertIdVec;
2100 SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ V0_EDGE ], vertIdVec);
2101 GetColumns(0, col1, col2 );
2102 const SMDS_MeshNode* node0 = col1->second.front();
2103 const SMDS_MeshNode* node1 = col1->second.back();
2104 TopoDS_Shape v0 = myHelper->GetSubShapeByNode( node0, myHelper->GetMeshDS());
2105 TopoDS_Shape v1 = myHelper->GetSubShapeByNode( node1, myHelper->GetMeshDS());
2106 if ( v0.ShapeType() == TopAbs_VERTEX ) {
2107 nbInserted += SMESH_Block::Insert( v0, vertIdVec[ 0 ], shapeMap);
2109 if ( v1.ShapeType() == TopAbs_VERTEX ) {
2110 nbInserted += SMESH_Block::Insert( v1, vertIdVec[ 1 ], shapeMap);
2114 SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ V1_EDGE ], vertIdVec);
2115 GetColumns(1, col1, col2 );
2116 node0 = col2->second.front();
2117 node1 = col2->second.back();
2118 v0 = myHelper->GetSubShapeByNode( node0, myHelper->GetMeshDS());
2119 v1 = myHelper->GetSubShapeByNode( node1, myHelper->GetMeshDS());
2120 if ( v0.ShapeType() == TopAbs_VERTEX ) {
2121 nbInserted += SMESH_Block::Insert( v0, vertIdVec[ 0 ], shapeMap);
2123 if ( v1.ShapeType() == TopAbs_VERTEX ) {
2124 nbInserted += SMESH_Block::Insert( v1, vertIdVec[ 1 ], shapeMap);
2127 // TopoDS_Vertex V0, V1, Vcom;
2128 // TopExp::Vertices( myBaseEdge, V0, V1, true );
2129 // if ( !myIsForward ) std::swap( V0, V1 );
2131 // // bottom vertex IDs
2132 // SMESH_Block::GetEdgeVertexIDs( edgeIdVec[ _u0 ], vertIdVec);
2133 // SMESH_Block::Insert( V0, vertIdVec[ 0 ], shapeMap);
2134 // SMESH_Block::Insert( V1, vertIdVec[ 1 ], shapeMap);
2136 // TopoDS_Edge sideEdge = GetEdge( V0_EDGE );
2137 // if ( sideEdge.IsNull() || !TopExp::CommonVertex( botEdge, sideEdge, Vcom ))
2140 // // insert one side edge
2142 // if ( Vcom.IsSame( V0 )) edgeID = edgeIdVec[ _v0 ];
2143 // else edgeID = edgeIdVec[ _v1 ];
2144 // SMESH_Block::Insert( sideEdge, edgeID, shapeMap);
2146 // // top vertex of the side edge
2147 // SMESH_Block::GetEdgeVertexIDs( edgeID, vertIdVec);
2148 // TopoDS_Vertex Vtop = TopExp::FirstVertex( sideEdge );
2149 // if ( Vcom.IsSame( Vtop ))
2150 // Vtop = TopExp::LastVertex( sideEdge );
2151 // SMESH_Block::Insert( Vtop, vertIdVec[ 1 ], shapeMap);
2153 // // other side edge
2154 // sideEdge = GetEdge( V1_EDGE );
2155 // if ( sideEdge.IsNull() )
2157 // if ( edgeID = edgeIdVec[ _v1 ]) edgeID = edgeIdVec[ _v0 ];
2158 // else edgeID = edgeIdVec[ _v1 ];
2159 // SMESH_Block::Insert( sideEdge, edgeID, shapeMap);
2162 // TopoDS_Edge topEdge = GetEdge( TOP_EDGE );
2163 // SMESH_Block::Insert( topEdge, edgeIdVec[ _u1 ], shapeMap);
2165 // // top vertex of the other side edge
2166 // if ( !TopExp::CommonVertex( topEdge, sideEdge, Vcom ))
2168 // SMESH_Block::GetEdgeVertexIDs( edgeID, vertIdVec );
2169 // SMESH_Block::Insert( Vcom, vertIdVec[ 1 ], shapeMap);
2174 //================================================================================
2176 * \brief Dump ids of nodes of sides
2178 //================================================================================
2180 void StdMeshers_PrismAsBlock::TSideFace::dumpNodes(int nbNodes) const
2183 cout << endl << "NODES OF FACE "; SMESH_Block::DumpShapeID( myID, cout ) << endl;
2184 THorizontalEdgeAdaptor* hSize0 = (THorizontalEdgeAdaptor*) HorizCurve(0);
2185 cout << "Horiz side 0: "; hSize0->dumpNodes(nbNodes); cout << endl;
2186 THorizontalEdgeAdaptor* hSize1 = (THorizontalEdgeAdaptor*) HorizCurve(1);
2187 cout << "Horiz side 1: "; hSize1->dumpNodes(nbNodes); cout << endl;
2188 TVerticalEdgeAdaptor* vSide0 = (TVerticalEdgeAdaptor*) VertiCurve(0);
2189 cout << "Verti side 0: "; vSide0->dumpNodes(nbNodes); cout << endl;
2190 TVerticalEdgeAdaptor* vSide1 = (TVerticalEdgeAdaptor*) VertiCurve(1);
2191 cout << "Verti side 1: "; vSide1->dumpNodes(nbNodes); cout << endl;
2192 delete hSize0; delete hSize1; delete vSide0; delete vSide1;
2196 //================================================================================
2198 * \brief Creates TVerticalEdgeAdaptor
2199 * \param columnsMap - node column map
2200 * \param parameter - normalized parameter
2202 //================================================================================
2204 StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::
2205 TVerticalEdgeAdaptor( const TParam2ColumnMap* columnsMap, const double parameter)
2207 myNodeColumn = & getColumn( columnsMap, parameter )->second;
2210 //================================================================================
2212 * \brief Return coordinates for the given normalized parameter
2213 * \param U - normalized parameter
2214 * \retval gp_Pnt - coordinates
2216 //================================================================================
2218 gp_Pnt StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::Value(const Standard_Real U) const
2220 const SMDS_MeshNode* n1;
2221 const SMDS_MeshNode* n2;
2222 double r = getRAndNodes( myNodeColumn, U, n1, n2 );
2223 return gpXYZ(n1) * ( 1 - r ) + gpXYZ(n2) * r;
2226 //================================================================================
2228 * \brief Dump ids of nodes
2230 //================================================================================
2232 void StdMeshers_PrismAsBlock::TVerticalEdgeAdaptor::dumpNodes(int nbNodes) const
2235 for ( int i = 0; i < nbNodes && i < myNodeColumn->size(); ++i )
2236 cout << (*myNodeColumn)[i]->GetID() << " ";
2237 if ( nbNodes < myNodeColumn->size() )
2238 cout << myNodeColumn->back()->GetID();
2242 //================================================================================
2244 * \brief Return coordinates for the given normalized parameter
2245 * \param U - normalized parameter
2246 * \retval gp_Pnt - coordinates
2248 //================================================================================
2250 gp_Pnt StdMeshers_PrismAsBlock::THorizontalEdgeAdaptor::Value(const Standard_Real U) const
2252 return mySide->TSideFace::Value( U, myV );
2255 //================================================================================
2257 * \brief Dump ids of <nbNodes> first nodes and the last one
2259 //================================================================================
2261 void StdMeshers_PrismAsBlock::THorizontalEdgeAdaptor::dumpNodes(int nbNodes) const
2264 // Not bedugged code. Last node is sometimes incorrect
2265 const TSideFace* side = mySide;
2267 if ( mySide->IsComplex() )
2268 side = mySide->GetComponent(0,u);
2270 TParam2ColumnIt col, col2;
2271 TParam2ColumnMap* u2cols = side->GetColumns();
2272 side->GetColumns( u , col, col2 );
2274 int j, i = myV ? mySide->ColumnHeight()-1 : 0;
2276 const SMDS_MeshNode* n = 0;
2277 const SMDS_MeshNode* lastN
2278 = side->IsForward() ? u2cols->rbegin()->second[ i ] : u2cols->begin()->second[ i ];
2279 for ( j = 0; j < nbNodes && n != lastN; ++j )
2281 n = col->second[ i ];
2282 cout << n->GetID() << " ";
2283 if ( side->IsForward() )
2291 if ( mySide->IsComplex() )
2292 side = mySide->GetComponent(1,u);
2294 side->GetColumns( u , col, col2 );
2295 if ( n != col->second[ i ] )
2296 cout << col->second[ i ]->GetID();
2299 //================================================================================
2301 * \brief Return UV on pcurve for the given normalized parameter
2302 * \param U - normalized parameter
2303 * \retval gp_Pnt - coordinates
2305 //================================================================================
2307 gp_Pnt2d StdMeshers_PrismAsBlock::TPCurveOnHorFaceAdaptor::Value(const Standard_Real U) const
2309 TParam2ColumnIt u_col1, u_col2;
2310 double r = mySide->GetColumns( U, u_col1, u_col2 );
2311 gp_XY uv1 = mySide->GetNodeUV( myFace, u_col1->second[ myZ ]);
2312 gp_XY uv2 = mySide->GetNodeUV( myFace, u_col2->second[ myZ ]);
2313 return uv1 * ( 1 - r ) + uv2 * r;