0021096: EDF 1729 SMESH: Create a Projection 1D-2D algorithm

This commit is contained in:
eap 2011-10-13 05:35:11 +00:00
parent 9091ead47e
commit 669d7d6d0f
10 changed files with 509 additions and 90 deletions

View File

@ -55,6 +55,14 @@ Vertices</b>, which belong to the same edge of the created \b Face.
For groups of face, they must contain equal number of faces
and they must form topologically equal structures.
\n <b>Projection 1D-2D</b> algorithm differs from <b>Projection 2D</b>
algorithm in one point, namely it generates mesh segments on edges of
the face according to the projected 2D elements; thus it does not
require that edges to be meshed by any other 1D algorithm; moreover it
does not allow to mesh edges of the face using another algorithm via
definition of sub-meshes.
\n <b>Projection 3D</b> algorithm allows to define the mesh of a shape by
the projection of another already meshed shape. This algorithm works
only if all faces and edges of the target shape have been meshed as 1D-2D

View File

@ -931,6 +931,13 @@ module StdMeshers
{
};
/*!
* StdMeshers_Projection_1D2D: interface of "Projection 1D-2D" algorithm
*/
interface StdMeshers_Projection_1D2D : SMESH::SMESH_2D_Algo
{
};
/*!
* StdMeshers_Projection_1D: interface of "Projection 1D" algorithm
*/

View File

@ -242,6 +242,14 @@
output="QUAD,TRIA"
dim="2"/>
<algorithm type="Projection_1D2D"
label-id="Projection 1D-2D"
icon-id="mesh_algo_quad.png"
input=""
hypos="ProjectionSource2D"
output="QUAD,TRIA"
dim="2"/>
<algorithm type="Projection_3D"
label-id="Projection 3D"
icon-id="mesh_algo_hexa.png"

View File

@ -75,7 +75,8 @@ salomeinclude_HEADERS = \
StdMeshers_ImportSource.hxx \
StdMeshers_Import_1D.hxx \
StdMeshers_Import_1D2D.hxx \
StdMeshers_ViscousLayers.hxx
StdMeshers_ViscousLayers.hxx \
StdMeshers_Projection_1D2D.hxx
# Libraries targets
@ -129,7 +130,8 @@ dist_libStdMeshers_la_SOURCES = \
StdMeshers_ImportSource.cxx \
StdMeshers_Import_1D.cxx \
StdMeshers_Import_1D2D.cxx \
StdMeshers_ViscousLayers.cxx
StdMeshers_ViscousLayers.cxx \
StdMeshers_Projection_1D2D.cxx
# additionnal information to compil and link file

View File

@ -0,0 +1,216 @@
// Copyright (C) 2007-2011 CEA/DEN, EDF R&D, OPEN CASCADE
//
// Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
//
// This library is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 2.1 of the License.
//
// This library is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
// Lesser General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
//
// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
//
// File : StdMeshers_Projection_1D2D.cxx
// Module : SMESH
// Author : Edward AGAPOV (eap)
//
#include "StdMeshers_Projection_1D2D.hxx"
#include "SMESH_Gen.hxx"
#include "SMESH_MesherHelper.hxx"
#include "SMESH_subMeshEventListener.hxx"
#include "StdMeshers_FaceSide.hxx"
#include "StdMeshers_ProjectionSource2D.hxx"
#include "StdMeshers_ProjectionUtils.hxx"
#include <TopoDS.hxx>
#include <TopExp_Explorer.hxx>
using namespace std;
//#define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
//=======================================================================
//function : StdMeshers_Projection_1D2D
//purpose :
//=======================================================================
StdMeshers_Projection_1D2D::StdMeshers_Projection_1D2D(int hypId, int studyId, SMESH_Gen* gen)
:StdMeshers_Projection_2D(hypId, studyId, gen)
{
_name = "Projection_1D2D";
_requireDescretBoundary = false;
}
//=======================================================================
//function : Compute
//purpose :
//=======================================================================
bool StdMeshers_Projection_1D2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theShape)
{
if ( !StdMeshers_Projection_2D::Compute(theMesh, theShape))
return false;
SMESHDS_Mesh * meshDS = theMesh.GetMeshDS();
SMESHDS_SubMesh * faceSubMesh = meshDS->MeshElements( theShape );
if ( !faceSubMesh || faceSubMesh->NbElements() == 0 ) return false;
_quadraticMesh = faceSubMesh->GetElements()->next()->IsQuadratic();
SMESH_MesherHelper helper( theMesh );
helper.SetSubShape( theShape );
TopoDS_Face F = TopoDS::Face( theShape );
TError err;
TSideVector wires = StdMeshers_FaceSide::GetFaceWires( F, theMesh,
/*ignoreMediumNodes=*/false, err);
if ( err && !err->IsOK() )
return error( err );
for ( size_t iWire = 0; iWire < wires.size(); ++iWire )
{
vector<const SMDS_MeshNode*> nodes = wires[ iWire ]->GetOrderedNodes();
if ( nodes.empty() )
return error("No nodes found on a wire");
const bool checkExisting = wires[ iWire ]->NbSegments();
if ( _quadraticMesh )
{
for ( size_t i = 2; i < nodes.size(); i += 2 )
{
if ( checkExisting && meshDS->FindEdge( nodes[i-2], nodes[i], nodes[i-1]))
continue;
SMDS_MeshElement* e = meshDS->AddEdge( nodes[i-2], nodes[i], nodes[i-1] );
meshDS->SetMeshElementOnShape( e, nodes[i-1]->getshapeId() );
}
}
else
{
int edgeID = meshDS->ShapeToIndex( wires[ iWire ]->Edge(0) );
for ( size_t i = 1; i < nodes.size(); ++i )
{
if ( checkExisting && meshDS->FindEdge( nodes[i-1], nodes[i])
continue;
SMDS_MeshElement* e = meshDS->AddEdge( nodes[i-1], nodes[i] );
if ( nodes[i-1]->getshapeId() != edgeID &&
nodes[i ]->getshapeId() != edgeID )
{
edgeID = helper.GetMediumPos( nodes[i-1], nodes[i] ).first;
if ( edgeID < 1 ) edgeID = helper.GetSubShapeID();
}
meshDS->SetMeshElementOnShape( e, edgeID );
}
}
}
return true;
}
//=======================================================================
//function : Evaluate
//purpose :
//=======================================================================
bool StdMeshers_Projection_1D2D::Evaluate(SMESH_Mesh& theMesh,
const TopoDS_Shape& theShape,
MapShapeNbElems& aResMap)
{
if ( !StdMeshers_Projection_2D::Evaluate(theMesh,theShape,aResMap))
return false;
TopoDS_Shape srcFace = _sourceHypo->GetSourceFace();
SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
if ( !srcMesh ) srcMesh = & theMesh;
SMESH_subMesh* srcFaceSM = srcMesh->GetSubMesh( srcFace );
typedef StdMeshers_ProjectionUtils SPU;
SPU::TShapeShapeMap shape2ShapeMap;
SPU::InitVertexAssociation( _sourceHypo, shape2ShapeMap, theShape );
if ( !SPU::FindSubShapeAssociation( theShape, &theMesh, srcFace, srcMesh, shape2ShapeMap))
return error(COMPERR_BAD_SHAPE,"Topology of source and target faces seems different" );
MapShapeNbElems srcResMap;
if ( !srcFaceSM->IsMeshComputed() )
_gen->Evaluate( *srcMesh, srcFace, srcResMap);
SMESH_subMeshIteratorPtr smIt = srcFaceSM->getDependsOnIterator(/*includeSelf=*/false,
/*complexShapeFirst=*/true);
while ( smIt->more() )
{
SMESH_subMesh* srcSM = smIt->next();
TopAbs_ShapeEnum shapeType = srcSM->GetSubShape().ShapeType();
if ( shapeType == TopAbs_EDGE )
{
std::vector<int> aVec;
SMESHDS_SubMesh* srcSubMeshDS = srcSM->GetSubMeshDS();
if ( srcSubMeshDS && srcSubMeshDS->NbElements() )
{
aVec.resize(SMDSEntity_Last, 0);
SMDS_ElemIteratorPtr eIt = srcSubMeshDS->GetElements();
_quadraticMesh = ( eIt->more() && eIt->next()->IsQuadratic() );
aVec[SMDSEntity_Node] = srcSubMeshDS->NbNodes();
aVec[_quadraticMesh ? SMDSEntity_Quad_Edge : SMDSEntity_Edge] = srcSubMeshDS->NbElements();
}
else
{
if ( srcResMap.empty() )
if ( !_gen->Evaluate( *srcMesh, srcSM->GetSubShape(), srcResMap ))
return error(COMPERR_BAD_INPUT_MESH,"Source mesh not evaluatable");
aVec = srcResMap[ srcSM ];
if ( aVec.empty() )
return error(COMPERR_BAD_INPUT_MESH,"Source mesh is wrongly evaluated");
}
TopoDS_Shape tgtEdge = shape2ShapeMap( srcSM->GetSubShape() );
SMESH_subMesh* tgtSM = theMesh.GetSubMesh( tgtEdge );
aResMap.insert(std::make_pair(tgtSM,aVec));
}
if ( shapeType == TopAbs_VERTEX ) break;
}
return true;
}
//=======================================================================
//function : SetEventListener
//purpose : Sets a default event listener to submesh of the source face.
// whenSetToFaceSubMesh - submesh where algo is set
// After being set, event listener is notified on each event of a submesh.
// This method is called when a submesh gets HYP_OK algo_state.
// Arranges that CLEAN event is translated from source submesh to
// the whenSetToFaceSubMesh submesh.
//=======================================================================
void StdMeshers_Projection_1D2D::SetEventListener(SMESH_subMesh* whenSetToFaceSubMesh)
{
// set a listener of events on a source submesh
StdMeshers_Projection_2D::SetEventListener(whenSetToFaceSubMesh);
// set a listener to the target FACE submesh in order to update submehses
// of EDGEs according to events on the target FACE submesh
// fill a listener data with submeshes of EDGEs
SMESH_subMeshEventListenerData* data =
new SMESH_subMeshEventListenerData(/*isDeletable=*/true);
SMESH_Mesh* mesh = whenSetToFaceSubMesh->GetFather();
TopExp_Explorer eExp( whenSetToFaceSubMesh->GetSubShape(), TopAbs_EDGE );
for ( ; eExp.More(); eExp.Next() )
data->mySubMeshes.push_back( mesh->GetSubMesh( eExp.Current() ));
// set a listener
SMESH_subMeshEventListener* listener =
new SMESH_subMeshEventListener(/*isDeletable=*/true);
whenSetToFaceSubMesh->SetEventListener( listener, data, whenSetToFaceSubMesh );
}

View File

@ -0,0 +1,52 @@
// Copyright (C) 2007-2011 CEA/DEN, EDF R&D, OPEN CASCADE
//
// Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
//
// This library is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 2.1 of the License.
//
// This library is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
// Lesser General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
//
// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
//
// File : StdMeshers_Projection_1D2D.hxx
// Module : SMESH
//
#ifndef _SMESH_Projection_1D2D_HXX_
#define _SMESH_Projection_1D2D_HXX_
#include "StdMeshers_Projection_2D.hxx"
class STDMESHERS_EXPORT StdMeshers_Projection_1D2D: public StdMeshers_Projection_2D
{
public:
StdMeshers_Projection_1D2D(int hypId, int studyId, SMESH_Gen* gen);
virtual bool Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape);
virtual bool Evaluate(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape,
MapShapeNbElems& aResMap);
/*!
* \brief Sets a default event listener to submesh of the source face
* \param whenSetToSubMesh - submesh where algo is set
*
* After being set, event listener is notified on each event of a submesh.
* This method is called when a submesh gets HYP_OK algo_state.
* Arranges that CLEAN event is translated from source submesh to
* the whenSetToSubMesh submesh.
*/
virtual void SetEventListener(SMESH_subMesh* whenSetToSubMesh);
};
#endif

View File

@ -190,13 +190,22 @@ namespace {
*/
//================================================================================
bool isOldNode( const SMDS_MeshNode* node )
bool isOldNode( const SMDS_MeshNode* node/*, const bool is1DComputed*/ )
{
// old nodes are shared by edges and new ones are shared
// only by faces created by mapper
SMDS_ElemIteratorPtr invEdge = node->GetInverseElementIterator(SMDSAbs_Edge);
bool isOld = invEdge->more();
return isOld;
//if ( is1DComputed )
{
SMDS_ElemIteratorPtr invEdge = node->GetInverseElementIterator(SMDSAbs_Edge);
bool isOld = invEdge->more();
return isOld;
}
// else
// {
// SMDS_ElemIteratorPtr invFace = node->GetInverseElementIterator(SMDSAbs_Face);
// bool isNew = invFace->more();
// return !isNew;
// }
}
//================================================================================
@ -520,7 +529,7 @@ namespace {
// Make new faces
// prepare the helper adding quadratic elements if necessary
// prepare the helper to adding quadratic elements if necessary
SMESH_MesherHelper helper( *tgtMesh );
helper.SetSubShape( tgtFace );
helper.IsQuadraticSubMesh( tgtFace );
@ -582,14 +591,15 @@ namespace {
const TopoDS_Face& srcFace,
SMESH_Mesh * tgtMesh,
SMESH_Mesh * srcMesh,
const TAssocTool::TShapeShapeMap& shape2ShapeMap)
const TAssocTool::TShapeShapeMap& shape2ShapeMap,
const bool is1DComputed)
{
// 1) Preparation
// get ordered src EDGEs
TError err;
TSideVector srcWires =
StdMeshers_FaceSide::GetFaceWires( srcFace, *srcMesh,/*theIgnoreMediumNodes = */false, err);
StdMeshers_FaceSide::GetFaceWires( srcFace, *srcMesh,/*ignoreMediumNodes = */false, err);
if ( err && !err->IsOK() )
return false;
@ -605,7 +615,8 @@ namespace {
tgtWires[ iW ].reset( new StdMeshers_FaceSide( tgtFace, tgtEdges, tgtMesh,
/*theIsForward = */ true,
/*theIgnoreMediumNodes = */false));
if ( srcWires[iW]->GetUVPtStruct().size() !=
if ( is1DComputed &&
srcWires[iW]->GetUVPtStruct().size() !=
tgtWires[iW]->GetUVPtStruct().size())
return false;
}
@ -614,15 +625,11 @@ namespace {
gp_Trsf2d trsf;
{
// get ordered nodes data
const vector<UVPtStruct>& srcUVs = srcWires[0]->GetUVPtStruct();
const vector<UVPtStruct>& tgtUVs = tgtWires[0]->GetUVPtStruct();
// get 2 pairs of corresponding UVs
gp_XY srcP0( srcUVs[0].u, srcUVs[0].v );
gp_XY srcP1( srcUVs[1].u, srcUVs[1].v );
gp_XY tgtP0( tgtUVs[0].u, tgtUVs[0].v );
gp_XY tgtP1( tgtUVs[1].u, tgtUVs[1].v );
gp_Pnt2d srcP0 = srcWires[0]->Value2d(0.0);
gp_Pnt2d srcP1 = srcWires[0]->Value2d(0.333);
gp_Pnt2d tgtP0 = tgtWires[0]->Value2d(0.0);
gp_Pnt2d tgtP1 = tgtWires[0]->Value2d(0.333);
// make transformation
gp_Trsf2d fromTgtCS, toSrcCS; // from/to global CS
@ -636,12 +643,10 @@ namespace {
// check transformation
const double tol = 1e-5 * gp_Vec2d( srcP0, srcP1 ).Magnitude();
const int nbCheckPnt = Min( 10, srcUVs.size()-2 );
const int dP = ( srcUVs.size()-2 ) / nbCheckPnt;
for ( unsigned iP = 2; iP < srcUVs.size(); iP += dP )
for ( double u = 0.12; u < 1.; u += 0.1 )
{
gp_Pnt2d srcUV( srcUVs[iP].u, srcUVs[iP].v );
gp_Pnt2d tgtUV( tgtUVs[iP].u, tgtUVs[iP].v );
gp_Pnt2d srcUV = srcWires[0]->Value2d( u );
gp_Pnt2d tgtUV = tgtWires[0]->Value2d( u );
gp_Pnt2d tgtUV2 = srcUV.Transformed( trsf );
if ( tgtUV.Distance( tgtUV2 ) > tol )
return false;
@ -656,26 +661,45 @@ namespace {
// fill src2tgtNodes in with nodes on EDGEs
for ( unsigned iW = 0; iW < srcWires.size(); ++iW )
{
const vector<UVPtStruct>& srcUVs = srcWires[iW]->GetUVPtStruct();
const vector<UVPtStruct>& tgtUVs = tgtWires[iW]->GetUVPtStruct();
for ( unsigned i = 0; i < srcUVs.size(); ++i )
src2tgtNodes.insert( make_pair( srcUVs[i].node, tgtUVs[i].node ));
}
if ( is1DComputed )
{
const vector<UVPtStruct>& srcUVs = srcWires[iW]->GetUVPtStruct();
const vector<UVPtStruct>& tgtUVs = tgtWires[iW]->GetUVPtStruct();
for ( unsigned i = 0; i < srcUVs.size(); ++i )
src2tgtNodes.insert( make_pair( srcUVs[i].node, tgtUVs[i].node ));
}
else
{
for ( int iE = 0; iE < srcWires[iW]->NbEdges(); ++iE )
{
TopoDS_Vertex srcV = srcWires[iW]->FirstVertex(iE);
TopoDS_Vertex tgtV = tgtWires[iW]->FirstVertex(iE);
const SMDS_MeshNode* srcNode = SMESH_Algo::VertexNode( srcV, srcMesh->GetMeshDS() );
const SMDS_MeshNode* tgtNode = SMESH_Algo::VertexNode( tgtV, tgtMesh->GetMeshDS() );
if ( tgtNode && srcNode )
src2tgtNodes.insert( make_pair( srcNode, tgtNode ));
}
}
// make elements
SMESHDS_SubMesh* srcSubDS = srcMesh->GetMeshDS()->MeshElements( srcFace );
SMESH_MesherHelper helper( *tgtMesh );
helper.SetSubShape( tgtFace );
helper.IsQuadraticSubMesh( tgtFace );
if ( is1DComputed )
helper.IsQuadraticSubMesh( tgtFace );
else
helper.SetIsQuadratic( srcSubDS->GetElements()->next()->IsQuadratic() );
helper.SetElementsOnShape( true );
Handle(Geom_Surface) tgtSurface = BRep_Tool::Surface( tgtFace );
SMESHDS_Mesh* tgtMeshDS = tgtMesh->GetMeshDS();
SMESH_MesherHelper srcHelper( *srcMesh );
srcHelper.SetSubShape( srcFace );
const SMDS_MeshNode* nullNode = 0;
SMESHDS_SubMesh* srcSubDS = srcMesh->GetMeshDS()->MeshElements( srcFace );
SMDS_ElemIteratorPtr elemIt = srcSubDS->GetElements();
vector< const SMDS_MeshNode* > tgtNodes;
bool uvOK;
@ -693,11 +717,29 @@ namespace {
// create a new node
gp_Pnt2d srcUV = srcHelper.GetNodeUV( srcFace, srcNode,
elem->GetNode( helper.WrapIndex(i+1,nbN)), &uvOK);
gp_Pnt2d tgtUV = srcUV.Transformed( trsf );
gp_Pnt tgtP = tgtSurface->Value( tgtUV.X(), tgtUV.Y() );
SMDS_MeshNode* n = helper.AddNode( tgtP.X(), tgtP.Y(), tgtP.Z() );
n->SetPosition( new SMDS_FacePosition( tgtUV.X(), tgtUV.Y() ));
SMDS_MeshNode* n = tgtMeshDS->AddNode( tgtP.X(), tgtP.Y(), tgtP.Z() );
switch ( srcNode->GetPosition()->GetTypeOfPosition() )
{
case SMDS_TOP_FACE: {
tgtMeshDS->SetNodeOnFace( n, helper.GetSubShapeID(), tgtUV.X(), tgtUV.Y() );
break;
}
case SMDS_TOP_EDGE: {
TopoDS_Shape srcEdge = srcHelper.GetSubShapeByNode( srcNode, srcHelper.GetMeshDS() );
TopoDS_Shape tgtEdge = shape2ShapeMap( srcEdge );
double U = srcHelper.GetNodeU( TopoDS::Edge( srcEdge ), srcNode );
tgtMeshDS->SetNodeOnEdge( n, TopoDS::Edge( tgtEdge ), U);
break;
}
case SMDS_TOP_VERTEX: {
TopoDS_Shape srcV = srcHelper.GetSubShapeByNode( srcNode, srcHelper.GetMeshDS() );
TopoDS_Shape tgtV = shape2ShapeMap( srcV );
tgtMeshDS->SetNodeOnVertex( n, TopoDS::Vertex( tgtV ));
break;
}
}
srcN_tgtN->second = n;
}
tgtNodes[i] = srcN_tgtN->second;
@ -711,8 +753,7 @@ namespace {
}
return true;
} // bool projectBy2DSimilarity()
} // bool projectBy2DSimilarity(...)
} // namespace
@ -767,11 +808,26 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
}
// --------------------------
// Simple cases of projection
// --------------------------
// find out if EDGEs are meshed or not
bool is1DComputed = false;
SMESH_subMeshIteratorPtr smIt = tgtSubMesh->getDependsOnIterator(/*includeSelf=*/false,
/*complexShapeFirst=*/true);
while ( smIt->more() && !is1DComputed )
{
SMESH_subMesh* sm = smIt->next();
if ( sm->GetSubShape().ShapeType() == TopAbs_EDGE )
is1DComputed = sm->IsMeshComputed();
}
// try to project from same face with different location
if ( projectPartner( tgtFace, srcFace, tgtMesh, srcMesh, shape2ShapeMap ))
return true;
if ( projectBy2DSimilarity( tgtFace, srcFace, tgtMesh, srcMesh, shape2ShapeMap ))
if ( projectBy2DSimilarity( tgtFace, srcFace, tgtMesh, srcMesh, shape2ShapeMap, is1DComputed ))
return true;
// --------------------
@ -863,23 +919,27 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
// in failure case
MeshCleaner cleaner( tgtSubMesh );
SMESH_MeshEditor editor( tgtMesh );
// -------------------------------------------------------------------------
// mapper doesn't take care of nodes already existing on edges and vertices,
// so we must merge nodes created by it with existing ones
// -------------------------------------------------------------------------
SMESH_MeshEditor editor( tgtMesh );
SMESH_MeshEditor::TListOfListOfNodes groupsOfNodes;
// Make groups of nodes to merge
// loop on edge and vertex submeshes of a target face
SMESH_subMeshIteratorPtr smIt = tgtSubMesh->getDependsOnIterator(false,false);
smIt = tgtSubMesh->getDependsOnIterator(/*includeSelf=*/false,/*complexShapeFirst=*/false);
while ( smIt->more() )
{
SMESH_subMesh* sm = smIt->next();
SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
if ( !is1DComputed && sm->GetSubShape().ShapeType() == TopAbs_EDGE )
break;
// Sort new and old nodes of a submesh separately
bool isSeam = helper.IsRealSeam( sm->GetId() );
@ -909,7 +969,10 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
switch ( node->GetPosition()->GetTypeOfPosition() )
{
case SMDS_TOP_VERTEX: {
pos2nodes.insert( make_pair( 0, node ));
if ( !is1DComputed && !pos2nodes.empty() )
u2nodesMaps[isOld ? NEW_NODES : OLD_NODES].insert( make_pair( 0, node ));
else
pos2nodes.insert( make_pair( 0, node ));
break;
}
case SMDS_TOP_EDGE: {
@ -932,6 +995,12 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
// does not make segments/nodes on degenerated edges
continue;
if ( u2nodesMaps[ OLD_NODES ].size() == 0 &&
sm->GetSubShape().ShapeType() == TopAbs_VERTEX &&
!is1DComputed )
// old nodes are optional on vertices in the case of 1D-2D projection
continue;
RETURN_BAD_RESULT("Different nb of old and new nodes on shape #"<< sm->GetId() <<" "<<
u2nodesMaps[ OLD_NODES ].size() << " != " <<
u2nodesMaps[ NEW_NODES ].size());
@ -962,6 +1031,24 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
if ( nbFaceBeforeMerge != nbFaceAtferMerge )
return error(COMPERR_BAD_INPUT_MESH, "Probably invalid node parameters on geom faces");
// ----------------------------------------------------------------
// The mapper can't create quadratic elements, so convert if needed
// ----------------------------------------------------------------
faceIt = srcSubMesh->GetSubMeshDS()->GetElements();
bool srcIsQuad = faceIt->next()->IsQuadratic();
faceIt = tgtSubMesh->GetSubMeshDS()->GetElements();
bool tgtIsQuad = faceIt->next()->IsQuadratic();
if ( srcIsQuad && !tgtIsQuad )
{
TIDSortedElemSet tgtFaces;
faceIt = tgtSubMesh->GetSubMeshDS()->GetElements();
while ( faceIt->more() )
tgtFaces.insert( tgtFaces.end(), faceIt->next() );
editor.ConvertToQuadratic(/*theForce3d=*/false, tgtFaces);
}
// ---------------------------
// Check elements orientation
// ---------------------------
@ -1017,9 +1104,9 @@ bool StdMeshers_Projection_2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape&
//purpose :
//=======================================================================
bool StdMeshers_Projection_2D::Evaluate(SMESH_Mesh& theMesh,
bool StdMeshers_Projection_2D::Evaluate(SMESH_Mesh& theMesh,
const TopoDS_Shape& theShape,
MapShapeNbElems& aResMap)
MapShapeNbElems& aResMap)
{
if ( !_sourceHypo )
return false;
@ -1045,40 +1132,31 @@ bool StdMeshers_Projection_2D::Evaluate(SMESH_Mesh& theMesh,
TopoDS_Face srcFace = TopoDS::Face( shape2ShapeMap( tgtFace ).Oriented(TopAbs_FORWARD));
// ----------------------------------------------
// Assure that mesh on a source Face is computed
// ----------------------------------------------
// -------------------------------------------------------
// Assure that mesh on a source Face is computed/evaluated
// -------------------------------------------------------
std::vector<int> aVec;
SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( srcFace );
if ( srcSubMesh->IsMeshComputed() )
{
aVec.resize( SMDSEntity_Last, 0 );
aVec[SMDSEntity_Node] = srcSubMesh->GetSubMeshDS()->NbNodes();
if ( !srcSubMesh->IsMeshComputed() )
return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
std::vector<int> aVec(SMDSEntity_Last);
for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i] = 0;
aVec[SMDSEntity_Node] = srcSubMesh->GetSubMeshDS()->NbNodes();
//bool quadratic = false;
SMDS_ElemIteratorPtr elemIt = srcSubMesh->GetSubMeshDS()->GetElements();
while ( elemIt->more() ) {
const SMDS_MeshElement* E = elemIt->next();
if( E->NbNodes()==3 ) {
aVec[SMDSEntity_Triangle]++;
}
else if( E->NbNodes()==4 ) {
aVec[SMDSEntity_Quadrangle]++;
}
else if( E->NbNodes()==6 && E->IsQuadratic() ) {
aVec[SMDSEntity_Quad_Triangle]++;
}
else if( E->NbNodes()==8 && E->IsQuadratic() ) {
aVec[SMDSEntity_Quad_Quadrangle]++;
}
else {
aVec[SMDSEntity_Polygon]++;
}
SMDS_ElemIteratorPtr elemIt = srcSubMesh->GetSubMeshDS()->GetElements();
while ( elemIt->more() )
aVec[ elemIt->next()->GetEntityType() ]++;
}
else
{
MapShapeNbElems tmpResMap;
MapShapeNbElems& srcResMap = (srcMesh == tgtMesh) ? aResMap : tmpResMap;
if ( !_gen->Evaluate( *srcMesh, srcShape, srcResMap ))
return error(COMPERR_BAD_INPUT_MESH,"Source mesh not evaluatable");
aVec = srcResMap[ srcSubMesh ];
if ( aVec.empty() )
return error(COMPERR_BAD_INPUT_MESH,"Source mesh is wrongly evaluated");
}
SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);

View File

@ -19,13 +19,8 @@
//
// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
//
// SMESH SMESH_I : idl implementation based on 'SMESH' unit's calsses
// File : StdMeshers_Projection_3D_i.cxx
// Moved here from SMESH_Projection_3D_i.cxx
// Author : Paul RASCLE, EDF
// File : StdMeshers_Projection_1D_2D_3D_i.cxx
// Module : SMESH
// $Header$
//
#include "StdMeshers_Projection_1D_2D_3D_i.hxx"
@ -106,6 +101,40 @@ StdMeshers_Projection_2D_i::~StdMeshers_Projection_2D_i()
}
//=============================================================================
/*!
* StdMeshers_Projection_1D2D_i::StdMeshers_Projection_1D2D_i
*/
//=============================================================================
StdMeshers_Projection_1D2D_i::StdMeshers_Projection_1D2D_i( PortableServer::POA_ptr thePOA,
int theStudyId,
::SMESH_Gen* theGenImpl )
: SALOME::GenericObj_i( thePOA ),
SMESH_Hypothesis_i( thePOA ),
SMESH_Algo_i( thePOA ),
SMESH_2D_Algo_i( thePOA )
{
MESSAGE( "StdMeshers_Projection_1D2D_i::StdMeshers_Projection_1D2D_i" );
myBaseImpl = new ::StdMeshers_Projection_1D2D( theGenImpl->GetANewId(),
theStudyId,
theGenImpl );
}
//-----------------------------------------------------------------------------
StdMeshers_Projection_1D2D_i::~StdMeshers_Projection_1D2D_i()
{
MESSAGE( "StdMeshers_Projection_1D2D_i::~StdMeshers_Projection_1D2D_i" );
}
//-----------------------------------------------------------------------------
::StdMeshers_Projection_1D2D* StdMeshers_Projection_1D2D_i::GetImpl()
{
MESSAGE( "StdMeshers_Projection_1D2D_i::GetImpl" );
return ( ::StdMeshers_Projection_1D2D* )myBaseImpl;
}
//=============================================================================
/*!
* StdMeshers_Projection_1D_i::StdMeshers_Projection_1D_i

View File

@ -19,16 +19,11 @@
//
// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
//
// SMESH SMESH_I : idl implementation based on 'SMESH' unit's calsses
// File : StdMeshers_Hexa_3D_i.hxx
// Moved here from SMESH_Hexa_3D_i.hxx
// Author : Paul RASCLE, EDF
// File : StdMeshers_Projection_1D_2D_3D_i.hxx
// Module : SMESH
// $Header$
//
#ifndef _SMESH_Projection_3D_I_HXX_
#define _SMESH_Projection_3D_I_HXX_
#ifndef _SMESH_Projection_1D2D3D_I_HXX_
#define _SMESH_Projection_1D2D3D_I_HXX_
#include <SALOMEconfig.h>
#include CORBA_SERVER_HEADER(SMESH_BasicHypothesis)
@ -37,6 +32,7 @@
#include "SMESH_2D_Algo_i.hxx"
#include "SMESH_3D_Algo_i.hxx"
#include "StdMeshers_Projection_1D.hxx"
#include "StdMeshers_Projection_1D2D.hxx"
#include "StdMeshers_Projection_2D.hxx"
#include "StdMeshers_Projection_3D.hxx"
@ -64,7 +60,7 @@ public:
};
// ======================================================
// Projection 3D algorithm
// Projection 2D algorithm
// ======================================================
class StdMeshers_Projection_2D_i:
@ -85,7 +81,28 @@ public:
};
// ======================================================
// Projection 3D algorithm
// Projection 1D-2D algorithm
// ======================================================
class StdMeshers_Projection_1D2D_i:
public virtual POA_StdMeshers::StdMeshers_Projection_1D2D,
public virtual SMESH_2D_Algo_i
{
public:
// Constructor
StdMeshers_Projection_1D2D_i( PortableServer::POA_ptr thePOA,
int theStudyId,
::SMESH_Gen* theGenImpl );
// Destructor
virtual ~StdMeshers_Projection_1D2D_i();
// Get implementation
::StdMeshers_Projection_1D2D* GetImpl();
};
// ======================================================
// Projection 1D algorithm
// ======================================================
class StdMeshers_Projection_1D_i:

View File

@ -189,6 +189,8 @@ STDMESHERS_I_EXPORT
aCreator = new StdHypothesisCreator_i<StdMeshers_Hexa_3D_i>;
else if (strcmp(aHypName, "Projection_1D") == 0)
aCreator = new StdHypothesisCreator_i<StdMeshers_Projection_1D_i>;
else if (strcmp(aHypName, "Projection_1D2D") == 0)
aCreator = new StdHypothesisCreator_i<StdMeshers_Projection_1D2D_i>;
else if (strcmp(aHypName, "Projection_2D") == 0)
aCreator = new StdHypothesisCreator_i<StdMeshers_Projection_2D_i>;
else if (strcmp(aHypName, "Projection_3D") == 0)