mirror of
https://git.salome-platform.org/gitpub/modules/smesh.git
synced 2024-12-26 09:20:34 +05:00
52618: Import_1D2D fails to import from a mesh group onto a group of geom faces
+ fix StdMeshers_ViscousLayers.cxx for http://www.salome-platform.org/forum/forum_14/777944993
This commit is contained in:
parent
fd1943809d
commit
066f625a46
@ -58,10 +58,9 @@ same priority.
|
|||||||
|
|
||||||
If meshing parameters are defined on sub-meshes of the same priority,
|
If meshing parameters are defined on sub-meshes of the same priority,
|
||||||
for example different 1D hypotheses are assigned to two faces sharing
|
for example different 1D hypotheses are assigned to two faces sharing
|
||||||
an edge, an arbitrary algorithm/hypothesis will be used for
|
an edge, the hypothesis assigned to a sub-shape with a lower ID will
|
||||||
meshing. This indeterminacy can be fixed by
|
be used for meshing. You can \ref submesh_order_anchor "change" mutual
|
||||||
\ref submesh_order_anchor "Changing" mutual priority of such
|
priority of such concurrent sub-meshes.
|
||||||
concurrent sub-meshes.
|
|
||||||
|
|
||||||
|
|
||||||
\n Construction of a sub-mesh consists of:
|
\n Construction of a sub-mesh consists of:
|
||||||
|
@ -876,7 +876,7 @@ module StdMeshers
|
|||||||
* interface of "Viscous Layers" hypothesis.
|
* interface of "Viscous Layers" hypothesis.
|
||||||
* This hypothesis specifies parameters of layers of prisms to build
|
* This hypothesis specifies parameters of layers of prisms to build
|
||||||
* near mesh boundary. This hypothesis can be used by several 3D algorithms:
|
* near mesh boundary. This hypothesis can be used by several 3D algorithms:
|
||||||
* NETGEN 3D, GHS3D, Hexahedron(i,j,k)
|
* NETGEN 3D, Hexahedron(i,j,k), MG_Tetra
|
||||||
*/
|
*/
|
||||||
interface StdMeshers_ViscousLayers : SMESH::SMESH_Hypothesis
|
interface StdMeshers_ViscousLayers : SMESH::SMESH_Hypothesis
|
||||||
{
|
{
|
||||||
@ -907,7 +907,7 @@ module StdMeshers
|
|||||||
short GetNumberLayers();
|
short GetNumberLayers();
|
||||||
|
|
||||||
/*!
|
/*!
|
||||||
* Set factor (>1.0) of growth of layer thickness towards inside of mesh
|
* Set factor (>=1.0) of growth of layer thickness towards inside of mesh
|
||||||
*/
|
*/
|
||||||
void SetStretchFactor(in double factor) raises (SALOME::SALOME_Exception);
|
void SetStretchFactor(in double factor) raises (SALOME::SALOME_Exception);
|
||||||
double GetStretchFactor();
|
double GetStretchFactor();
|
||||||
@ -951,7 +951,7 @@ module StdMeshers
|
|||||||
short GetNumberLayers();
|
short GetNumberLayers();
|
||||||
|
|
||||||
/*!
|
/*!
|
||||||
* Set factor (>1.0) of growth of layer thickness towards inside of mesh
|
* Set factor (>=1.0) of growth of layer thickness towards inside of mesh
|
||||||
*/
|
*/
|
||||||
void SetStretchFactor(in double factor) raises (SALOME::SALOME_Exception);
|
void SetStretchFactor(in double factor) raises (SALOME::SALOME_Exception);
|
||||||
double GetStretchFactor();
|
double GetStretchFactor();
|
||||||
|
@ -79,13 +79,15 @@ public:
|
|||||||
virtual void NotifySubMeshesHypothesisModification();
|
virtual void NotifySubMeshesHypothesisModification();
|
||||||
void SetLibName(const char* theLibName);
|
void SetLibName(const char* theLibName);
|
||||||
|
|
||||||
//void SetParameters(const char *theParameters);
|
/*!
|
||||||
//char* GetParameters() const;
|
* \brief The returned value is used by NotifySubMeshesHypothesisModification()
|
||||||
|
* to decide to call subMesh->AlgoStateEngine( MODIF_HYP, hyp ) or not
|
||||||
|
* if subMesh is ready to be computed (algo+hyp==OK) but not yet computed.
|
||||||
|
* True result is reasonable for example if EventListeners depend on
|
||||||
|
* parameters of hypothesis.
|
||||||
|
*/
|
||||||
|
virtual bool DataDependOnParams() const { return false; }
|
||||||
|
|
||||||
// void SetLastParameters(const char* theParameters);
|
|
||||||
// char* GetLastParameters() const;
|
|
||||||
// void ClearParameters();
|
|
||||||
|
|
||||||
/*!
|
/*!
|
||||||
* \brief Initialize my parameter values by the mesh built on the geometry
|
* \brief Initialize my parameter values by the mesh built on the geometry
|
||||||
* \param theMesh - the built mesh
|
* \param theMesh - the built mesh
|
||||||
|
@ -1215,7 +1215,8 @@ void SMESH_Mesh::NotifySubMeshesHypothesisModification(const SMESH_Hypothesis* h
|
|||||||
// other possible changes are not interesting. (IPAL0052457 - assigning hyp performance pb)
|
// other possible changes are not interesting. (IPAL0052457 - assigning hyp performance pb)
|
||||||
if ( aSubMesh->GetComputeState() != SMESH_subMesh::COMPUTE_OK &&
|
if ( aSubMesh->GetComputeState() != SMESH_subMesh::COMPUTE_OK &&
|
||||||
aSubMesh->GetComputeState() != SMESH_subMesh::FAILED_TO_COMPUTE &&
|
aSubMesh->GetComputeState() != SMESH_subMesh::FAILED_TO_COMPUTE &&
|
||||||
aSubMesh->GetAlgoState() != SMESH_subMesh::MISSING_HYP )
|
aSubMesh->GetAlgoState() != SMESH_subMesh::MISSING_HYP &&
|
||||||
|
!hyp->DataDependOnParams() )
|
||||||
continue;
|
continue;
|
||||||
|
|
||||||
const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
|
const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
|
||||||
|
@ -1345,15 +1345,13 @@ class StdMeshersBuilder_UseExistingElements_1D2D(Mesh_Algorithm):
|
|||||||
# @param UseExisting if ==true - searches for the existing hypothesis created with
|
# @param UseExisting if ==true - searches for the existing hypothesis created with
|
||||||
# the same parameters, else (default) - creates a new one
|
# the same parameters, else (default) - creates a new one
|
||||||
def SourceFaces(self, groups, toCopyMesh=False, toCopyGroups=False, UseExisting=False):
|
def SourceFaces(self, groups, toCopyMesh=False, toCopyGroups=False, UseExisting=False):
|
||||||
for group in groups:
|
|
||||||
from salome.smesh.smeshBuilder import AssureGeomPublished
|
|
||||||
AssureGeomPublished( self.mesh, group )
|
|
||||||
compFun = lambda hyp, args: ( hyp.GetSourceFaces() == args[0] and \
|
compFun = lambda hyp, args: ( hyp.GetSourceFaces() == args[0] and \
|
||||||
hyp.GetCopySourceMesh() == args[1], args[2] )
|
hyp.GetCopySourceMesh() == args[1], args[2] )
|
||||||
hyp = self.Hypothesis("ImportSource2D", [groups, toCopyMesh, toCopyGroups],
|
hyp = self.Hypothesis("ImportSource2D", [groups, toCopyMesh, toCopyGroups],
|
||||||
UseExisting=UseExisting, CompareMethod=compFun)
|
UseExisting=UseExisting, CompareMethod=compFun, toAdd=False)
|
||||||
hyp.SetSourceFaces(groups)
|
hyp.SetSourceFaces(groups)
|
||||||
hyp.SetCopySourceMesh(toCopyMesh, toCopyGroups)
|
hyp.SetCopySourceMesh(toCopyMesh, toCopyGroups)
|
||||||
|
self.mesh.AddHypothesis(hyp, self.geom)
|
||||||
return hyp
|
return hyp
|
||||||
|
|
||||||
pass # end of StdMeshersBuilder_UseExistingElements_1D2D class
|
pass # end of StdMeshersBuilder_UseExistingElements_1D2D class
|
||||||
|
@ -61,6 +61,7 @@ class STDMESHERS_EXPORT StdMeshers_ImportSource1D : public SMESH_Hypothesis
|
|||||||
virtual std::istream & LoadFrom(std::istream & load);
|
virtual std::istream & LoadFrom(std::istream & load);
|
||||||
virtual bool SetParametersByMesh(const SMESH_Mesh* theMesh, const TopoDS_Shape& theShape);
|
virtual bool SetParametersByMesh(const SMESH_Mesh* theMesh, const TopoDS_Shape& theShape);
|
||||||
virtual bool SetParametersByDefaults(const TDefaults& dflts, const SMESH_Mesh* theMesh=0);
|
virtual bool SetParametersByDefaults(const TDefaults& dflts, const SMESH_Mesh* theMesh=0);
|
||||||
|
virtual bool DataDependOnParams() const { return true; }
|
||||||
void RestoreGroups(const std::vector<SMESH_Group*>& groups);
|
void RestoreGroups(const std::vector<SMESH_Group*>& groups);
|
||||||
|
|
||||||
void StoreResultGroups(const std::vector<SMESH_Group*>& groups,
|
void StoreResultGroups(const std::vector<SMESH_Group*>& groups,
|
||||||
|
@ -44,15 +44,20 @@
|
|||||||
#include "Utils_SALOME_Exception.hxx"
|
#include "Utils_SALOME_Exception.hxx"
|
||||||
#include "utilities.h"
|
#include "utilities.h"
|
||||||
|
|
||||||
|
#include <BRepClass_FaceClassifier.hxx>
|
||||||
|
#include <BRepTools.hxx>
|
||||||
#include <BRep_Builder.hxx>
|
#include <BRep_Builder.hxx>
|
||||||
#include <BRep_Tool.hxx>
|
#include <BRep_Tool.hxx>
|
||||||
|
#include <Bnd_B2d.hxx>
|
||||||
|
#include <GeomAPI_ProjectPointOnSurf.hxx>
|
||||||
|
#include <GeomAdaptor_Surface.hxx>
|
||||||
|
#include <Precision.hxx>
|
||||||
#include <TopExp.hxx>
|
#include <TopExp.hxx>
|
||||||
#include <TopExp_Explorer.hxx>
|
#include <TopExp_Explorer.hxx>
|
||||||
#include <TopoDS.hxx>
|
#include <TopoDS.hxx>
|
||||||
#include <TopoDS_Compound.hxx>
|
#include <TopoDS_Compound.hxx>
|
||||||
#include <TopoDS_Edge.hxx>
|
#include <TopoDS_Edge.hxx>
|
||||||
#include <TopoDS_Vertex.hxx>
|
#include <TopoDS_Vertex.hxx>
|
||||||
#include <Precision.hxx>
|
|
||||||
|
|
||||||
#include <numeric>
|
#include <numeric>
|
||||||
|
|
||||||
@ -187,9 +192,28 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
|
|||||||
|
|
||||||
Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
|
Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
|
||||||
const bool reverse =
|
const bool reverse =
|
||||||
( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace) == TopAbs_REVERSED );
|
( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace ) == TopAbs_REVERSED );
|
||||||
gp_Pnt p; gp_Vec du, dv;
|
gp_Pnt p; gp_Vec du, dv;
|
||||||
|
|
||||||
|
BRepClass_FaceClassifier classifier;
|
||||||
|
Bnd_B2d bndBox2d;
|
||||||
|
{
|
||||||
|
Standard_Real umin,umax,vmin,vmax;
|
||||||
|
BRepTools::UVBounds(geomFace,umin,umax,vmin,vmax);
|
||||||
|
gp_XY pmin( umin,vmin ), pmax( umax,vmax );
|
||||||
|
bndBox2d.Add( pmin );
|
||||||
|
bndBox2d.Add( pmax );
|
||||||
|
if ( helper.HasSeam() )
|
||||||
|
{
|
||||||
|
const int i = helper.GetPeriodicIndex();
|
||||||
|
pmin.SetCoord( i, helper.GetOtherParam( pmin.Coord( i )));
|
||||||
|
pmax.SetCoord( i, helper.GetOtherParam( pmax.Coord( i )));
|
||||||
|
bndBox2d.Add( pmin );
|
||||||
|
bndBox2d.Add( pmax );
|
||||||
|
}
|
||||||
|
bndBox2d.Enlarge( 1e-2 * Sqrt( bndBox2d.SquareExtent() ));
|
||||||
|
}
|
||||||
|
|
||||||
set<int> subShapeIDs;
|
set<int> subShapeIDs;
|
||||||
subShapeIDs.insert( shapeID );
|
subShapeIDs.insert( shapeID );
|
||||||
|
|
||||||
@ -244,7 +268,10 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
|
|||||||
|
|
||||||
StdMeshers_Import_1D::TNodeNodeMap* n2n;
|
StdMeshers_Import_1D::TNodeNodeMap* n2n;
|
||||||
StdMeshers_Import_1D::TElemElemMap* e2e;
|
StdMeshers_Import_1D::TElemElemMap* e2e;
|
||||||
vector<const SMDS_MeshNode*> newNodes;
|
vector<TopAbs_State> nodeState;
|
||||||
|
vector<const SMDS_MeshNode*> newNodes; // of a face
|
||||||
|
set <const SMDS_MeshNode*> bndNodes; // nodes classified ON
|
||||||
|
|
||||||
for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
|
for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
|
||||||
{
|
{
|
||||||
const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
|
const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
|
||||||
@ -257,58 +284,126 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
|
|||||||
const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
|
const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
|
||||||
minGroupTol = std::min( groupTol, minGroupTol );
|
minGroupTol = std::min( groupTol, minGroupTol );
|
||||||
|
|
||||||
|
GeomAdaptor_Surface S( surface );
|
||||||
|
const double clsfTol = Min( S.UResolution( 0.1 * groupTol ),
|
||||||
|
S.VResolution( 0.1 * groupTol ));
|
||||||
|
|
||||||
|
StdMeshers_Import_1D::TNodeNodeMap::iterator n2nIt;
|
||||||
|
pair< StdMeshers_Import_1D::TNodeNodeMap::iterator, bool > it_isnew;
|
||||||
|
|
||||||
SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
|
SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
|
||||||
SMDS_MeshNode *tmpNode = helper.AddNode(0,0,0);
|
|
||||||
gp_XY uv( Precision::Infinite(), Precision::Infinite() );
|
|
||||||
while ( srcElems->more() ) // loop on group contents
|
while ( srcElems->more() ) // loop on group contents
|
||||||
{
|
{
|
||||||
const SMDS_MeshElement* face = srcElems->next();
|
const SMDS_MeshElement* face = srcElems->next();
|
||||||
|
|
||||||
// find or create nodes of a new face
|
// find or create nodes of a new face
|
||||||
newNodes.resize( face->NbNodes() );
|
nodeState.resize( face->NbNodes() );
|
||||||
|
newNodes.resize( nodeState.size() );
|
||||||
newNodes.back() = 0;
|
newNodes.back() = 0;
|
||||||
int nbCreatedNodes = 0;
|
int nbCreatedNodes = 0;
|
||||||
|
bool isOut = false, isIn = false; // if at least one node isIn - do not classify other nodes
|
||||||
SMDS_MeshElement::iterator node = face->begin_nodes();
|
SMDS_MeshElement::iterator node = face->begin_nodes();
|
||||||
for ( size_t i = 0; i < newNodes.size(); ++i, ++node )
|
for ( size_t i = 0; i < newNodes.size(); ++i, ++node )
|
||||||
{
|
{
|
||||||
StdMeshers_Import_1D::TNodeNodeMap::iterator n2nIt =
|
SMESH_TNodeXYZ nXYZ = *node;
|
||||||
n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 )).first;
|
nodeState[ i ] = TopAbs_UNKNOWN;
|
||||||
if ( n2nIt->second )
|
|
||||||
|
it_isnew = n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 ));
|
||||||
|
n2nIt = it_isnew.first;
|
||||||
|
|
||||||
|
const SMDS_MeshNode* & newNode = n2nIt->second;
|
||||||
|
if ( !it_isnew.second && !newNode )
|
||||||
|
break; // a node is mapped to NULL - it is OUT of the FACE
|
||||||
|
|
||||||
|
if ( newNode )
|
||||||
{
|
{
|
||||||
if ( !subShapeIDs.count( n2nIt->second->getshapeId() )) // node already on an EDGE
|
if ( !subShapeIDs.count( newNode->getshapeId() ))
|
||||||
break;
|
break; // node is Imported onto other FACE
|
||||||
|
if ( !isIn && bndNodes.count( *node ))
|
||||||
|
nodeState[ i ] = TopAbs_ON;
|
||||||
}
|
}
|
||||||
else
|
else
|
||||||
{
|
{
|
||||||
// find a pre-existing node
|
// find a node pre-existing on EDGE or VERTEX
|
||||||
dist2foundNodes.clear();
|
dist2foundNodes.clear();
|
||||||
existingNodeOcTr.NodesAround( SMESH_TNodeXYZ( *node ), dist2foundNodes, groupTol );
|
existingNodeOcTr.NodesAround( nXYZ, dist2foundNodes, groupTol );
|
||||||
if ( !dist2foundNodes.empty() )
|
if ( !dist2foundNodes.empty() )
|
||||||
(*n2nIt).second = dist2foundNodes.begin()->second;
|
|
||||||
}
|
|
||||||
if ( !n2nIt->second )
|
|
||||||
{
|
|
||||||
// find out if node lies on theShape
|
|
||||||
tmpNode->setXYZ( (*node)->X(), (*node)->Y(), (*node)->Z());
|
|
||||||
uv.SetCoord( Precision::Infinite(), Precision::Infinite() );
|
|
||||||
if ( helper.CheckNodeUV( geomFace, tmpNode, uv, groupTol, /*force=*/true ))
|
|
||||||
{
|
{
|
||||||
SMDS_MeshNode* newNode = tgtMesh->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
|
newNode = dist2foundNodes.begin()->second;
|
||||||
n2nIt->second = newNode;
|
nodeState[ i ] = TopAbs_ON;
|
||||||
tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
|
|
||||||
nbCreatedNodes++;
|
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
if ( !(newNodes[i] = n2nIt->second ))
|
|
||||||
|
if ( !newNode )
|
||||||
|
{
|
||||||
|
// find out if node lies on the surface of theShape
|
||||||
|
gp_XY uv( Precision::Infinite(), 0 );
|
||||||
|
isOut = ( !helper.CheckNodeUV( geomFace, *node, uv, groupTol, /*force=*/true ) ||
|
||||||
|
bndBox2d.IsOut( uv ));
|
||||||
|
if ( !isOut && !isIn ) // classify
|
||||||
|
{
|
||||||
|
classifier.Perform( geomFace, uv, clsfTol );
|
||||||
|
nodeState[i] = classifier.State();
|
||||||
|
isOut = ( nodeState[i] == TopAbs_OUT );
|
||||||
|
}
|
||||||
|
if ( !isOut ) // create a new node
|
||||||
|
{
|
||||||
|
newNode = tgtMesh->AddNode( nXYZ.X(), nXYZ.Y(), nXYZ.Z());
|
||||||
|
tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
|
||||||
|
nbCreatedNodes++;
|
||||||
|
if ( nodeState[i] == TopAbs_ON )
|
||||||
|
bndNodes.insert( *node );
|
||||||
|
else
|
||||||
|
isIn = true;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
if ( !(newNodes[i] = newNode ) || isOut )
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
|
|
||||||
if ( !newNodes.back() )
|
if ( !newNodes.back() )
|
||||||
continue; // not all nodes of the face lie on theShape
|
continue; // not all nodes of the face lie on theShape
|
||||||
|
|
||||||
|
if ( !isIn ) // if all nodes are on FACE boundary, a mesh face can be OUT
|
||||||
|
{
|
||||||
|
// check state of nodes created for other faces
|
||||||
|
for ( size_t i = 0; i < nodeState.size() && !isIn; ++i )
|
||||||
|
{
|
||||||
|
if ( nodeState[i] != TopAbs_UNKNOWN ) continue;
|
||||||
|
gp_XY uv = helper.GetNodeUV( geomFace, newNodes[i] );
|
||||||
|
classifier.Perform( geomFace, uv, clsfTol );
|
||||||
|
nodeState[i] = classifier.State();
|
||||||
|
isIn = ( nodeState[i] == TopAbs_IN );
|
||||||
|
}
|
||||||
|
if ( !isIn ) // classify face center
|
||||||
|
{
|
||||||
|
gp_XYZ gc( 0., 0., 0 );
|
||||||
|
for ( size_t i = 0; i < newNodes.size(); ++i )
|
||||||
|
gc += SMESH_TNodeXYZ( newNodes[i] );
|
||||||
|
gc /= newNodes.size();
|
||||||
|
|
||||||
|
TopLoc_Location loc;
|
||||||
|
GeomAPI_ProjectPointOnSurf& proj = helper.GetProjector( geomFace,
|
||||||
|
loc,
|
||||||
|
helper.MaxTolerance( geomFace ));
|
||||||
|
if ( !loc.IsIdentity() ) loc.Transformation().Inverted().Transforms( gc );
|
||||||
|
proj.Perform( gc );
|
||||||
|
if ( !proj.IsDone() || proj.NbPoints() < 1 )
|
||||||
|
continue;
|
||||||
|
Quantity_Parameter U,V;
|
||||||
|
proj.LowerDistanceParameters(U,V);
|
||||||
|
gp_XY uv( U,V );
|
||||||
|
classifier.Perform( geomFace, uv, clsfTol );
|
||||||
|
if ( classifier.State() != TopAbs_IN )
|
||||||
|
continue;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
// try to find already created face
|
// try to find already created face
|
||||||
SMDS_MeshElement * newFace = 0;
|
SMDS_MeshElement * newFace = 0;
|
||||||
if ( nbCreatedNodes == 0 &&
|
if ( nbCreatedNodes == 0 &&
|
||||||
tgtMesh->FindElement(newNodes, SMDSAbs_Face, /*noMedium=*/false))
|
tgtMesh->FindElement(newNodes, SMDSAbs_Face, /*noMedium=*/false))
|
||||||
continue; // repeated face in source groups already created
|
continue; // repeated face in source groups already created
|
||||||
|
|
||||||
// check future face orientation
|
// check future face orientation
|
||||||
const int nbCorners = face->NbCornerNodes();
|
const int nbCorners = face->NbCornerNodes();
|
||||||
@ -319,7 +414,7 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
|
|||||||
gp_Vec geomNorm;
|
gp_Vec geomNorm;
|
||||||
do
|
do
|
||||||
{
|
{
|
||||||
uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
|
gp_XY uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
|
||||||
surface->D1( uv.X(),uv.Y(), p, du,dv );
|
surface->D1( uv.X(),uv.Y(), p, du,dv );
|
||||||
geomNorm = reverse ? dv^du : du^dv;
|
geomNorm = reverse ? dv^du : du^dv;
|
||||||
}
|
}
|
||||||
@ -384,9 +479,15 @@ bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape &
|
|||||||
// }
|
// }
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
helper.GetMeshDS()->RemoveNode(tmpNode);
|
// Remove OUT nodes from n2n map
|
||||||
|
for ( n2nIt = n2n->begin(); n2nIt != n2n->end(); )
|
||||||
|
if ( !n2nIt->second )
|
||||||
|
n2n->erase( n2nIt++ );
|
||||||
|
else
|
||||||
|
++n2nIt;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
// ==========================================================
|
// ==========================================================
|
||||||
// Put nodes on geom edges and create edges on them;
|
// Put nodes on geom edges and create edges on them;
|
||||||
// check if the whole geom face is covered by imported faces
|
// check if the whole geom face is covered by imported faces
|
||||||
|
@ -1590,7 +1590,7 @@ SMESH_ComputeErrorPtr _ViscousBuilder::CheckHypotheses( SMESH_Mesh& mesh
|
|||||||
|
|
||||||
|
|
||||||
findSolidsWithLayers();
|
findSolidsWithLayers();
|
||||||
bool ok = findFacesWithLayers();
|
bool ok = findFacesWithLayers( true );
|
||||||
|
|
||||||
// remove _MeshOfSolid's of _SolidData's
|
// remove _MeshOfSolid's of _SolidData's
|
||||||
for ( size_t i = 0; i < _sdVec.size(); ++i )
|
for ( size_t i = 0; i < _sdVec.size(); ++i )
|
||||||
|
Loading…
Reference in New Issue
Block a user