Salome HOME
Update of CheckDone
[modules/smesh.git] / src / StdMeshers / StdMeshers_Projection_1D2D.cxx
1 // Copyright (C) 2007-2024  CEA, EDF, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
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, or (at your option) any later version.
10 //
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.
15 //
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
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 // File      : StdMeshers_Projection_1D2D.cxx
23 // Module    : SMESH
24 // Author    : Edward AGAPOV (eap)
25 //
26 #include "StdMeshers_Projection_1D2D.hxx"
27
28 #include "SMESHDS_Mesh.hxx"
29 #include "SMESH_Gen.hxx"
30 #include "SMESH_Mesh.hxx"
31 #include "SMESH_MesherHelper.hxx"
32 #include "SMESH_subMesh.hxx"
33 #include "SMESH_subMeshEventListener.hxx"
34 #include "StdMeshers_FaceSide.hxx"
35 #include "StdMeshers_ProjectionSource2D.hxx"
36 #include "StdMeshers_ProjectionUtils.hxx"
37
38 #include <TopoDS.hxx>
39 #include <TopExp_Explorer.hxx>
40
41 using namespace std;
42
43 namespace
44 {
45   // --------------------------------------------------------------------------------
46   /*!
47    * \brief an event listener updating submehses of EDGEs according to
48    *        events on the target FACE submesh
49    */
50   struct EventProparatorToEdges : public SMESH_subMeshEventListener
51   {
52     EventProparatorToEdges(): SMESH_subMeshEventListener(/*isDeletable=*/false,
53                                                          "Projection_1D2D::EventProparatorToEdges")
54     {}
55     static EventProparatorToEdges* Instance() { static EventProparatorToEdges E; return &E; }
56
57     static void Set(SMESH_subMesh* faceSubMesh)
58     {
59       SMESH_subMeshEventListenerData* edgeSubMeshes =
60         new SMESH_subMeshEventListenerData(/*isDeletable=*/true);
61       SMESH_Mesh* mesh = faceSubMesh->GetFather();
62       TopExp_Explorer eExp( faceSubMesh->GetSubShape(), TopAbs_EDGE );
63       for ( ; eExp.More(); eExp.Next() )
64         edgeSubMeshes->mySubMeshes.push_back( mesh->GetSubMesh( eExp.Current() ));
65
66       // set a listener
67       faceSubMesh->SetEventListener( Instance(), edgeSubMeshes, faceSubMesh );
68     }
69   };
70   // --------------------------------------------------------------------------------
71   /*!
72    * \brief Structure used to temporary remove EventProparatorToEdges from faceSubMesh
73    *  in order to prevent propagation of CLEAN event from FACE to EDGEs during 
74    *  StdMeshers_Projection_1D2D::Compute(). The CLEAN event is emitted by Pattern mapper
75    *  and causes removal of faces generated on adjacent FACEs.
76    */
77   struct UnsetterOfEventProparatorToEdges
78   {
79     SMESH_subMesh* _faceSubMesh;
80     UnsetterOfEventProparatorToEdges( SMESH_subMesh* faceSubMesh ):_faceSubMesh(faceSubMesh)
81     {
82       faceSubMesh->DeleteEventListener( EventProparatorToEdges::Instance() );
83     }
84     ~UnsetterOfEventProparatorToEdges()
85     {
86       EventProparatorToEdges::Set(_faceSubMesh);
87     }
88   };
89 }
90
91 //=======================================================================
92 //function : StdMeshers_Projection_1D2D
93 //purpose  : 
94 //=======================================================================
95
96 StdMeshers_Projection_1D2D::StdMeshers_Projection_1D2D(int hypId, SMESH_Gen* gen)
97   :StdMeshers_Projection_2D(hypId, gen)
98 {
99   _name = "Projection_1D2D";
100   _requireDiscreteBoundary = false;
101   _supportSubmeshes = true;
102 }
103
104 //=======================================================================
105 //function : Compute
106 //purpose  : 
107 //=======================================================================
108
109 bool StdMeshers_Projection_1D2D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theShape)
110 {
111   UnsetterOfEventProparatorToEdges eventBarrier( theMesh.GetSubMesh( theShape ));
112
113   // 1) Project faces
114
115   if ( !StdMeshers_Projection_2D::Compute(theMesh, theShape))
116     return false;
117
118   // 2) Create segments
119
120   SMESHDS_Mesh * meshDS = theMesh.GetMeshDS();
121
122   SMESHDS_SubMesh * faceSubMesh = meshDS->MeshElements( theShape );
123   if ( !faceSubMesh || faceSubMesh->NbElements() == 0 ) return false;
124   _quadraticMesh = faceSubMesh->GetElements()->next()->IsQuadratic();
125
126   SMESH_MesherHelper helper( theMesh );
127   helper.SetSubShape( theShape );
128
129   if ( _quadraticMesh )
130   {
131     // 2a) Move some medium nodes from FACE to EDGES. They are on FACE because
132     // EDGEs are discreteized later than FACE, in this case.
133
134     SMESH_MesherHelper posFixer( theMesh );
135     posFixer.ToFixNodeParameters( true );
136     SMDS_ElemIteratorPtr fIt = faceSubMesh->GetElements();
137     vector< const SMDS_MeshNode* > nodes;
138     double dummyU, tol = 1e-7;
139     while ( fIt->more() ) // loop on mesh faces created by StdMeshers_Projection_2D
140     {
141       const SMDS_MeshElement* f = fIt->next();
142       //if ( !f->IsQuadratic() ) continue;
143       nodes.assign( SMDS_MeshElement::iterator( f->interlacedNodesIterator() ),
144                     SMDS_MeshElement::iterator() );
145       nodes.push_back( nodes[0] );
146       for ( size_t i = 2; i < nodes.size(); i += 2 )
147       {
148         pair<int, TopAbs_ShapeEnum> idType = helper.GetMediumPos( nodes[i], nodes[i-2] );
149         if ( idType.second == TopAbs_EDGE &&
150              idType.first  != nodes[i-1]->getshapeId() )
151         {
152           faceSubMesh->RemoveNode( nodes[i-1] );
153           meshDS->SetNodeOnEdge( nodes[i-1], idType.first );
154           posFixer.SetSubShape( idType.first );
155           posFixer.CheckNodeU( TopoDS::Edge( posFixer.GetSubShape() ),
156                                nodes[i-1], dummyU=0., tol, /*force=*/true );
157         }
158       }
159     }
160   }
161   TopoDS_Face F = TopoDS::Face( theShape );
162   TError err;
163   TSideVector wires = StdMeshers_FaceSide::GetFaceWires( F, theMesh, /*ignoreMediumNodes=*/false,
164                                                          err, &helper);
165   if ( err && !err->IsOK() )
166     return error( err );
167
168   for ( size_t iWire = 0; iWire < wires.size(); ++iWire )
169   {
170     vector<const SMDS_MeshNode*> nodes = wires[ iWire ]->GetOrderedNodes();
171     if ( nodes.empty() )
172       return error("Wrong nodes on a wire");
173
174     // check that all nodes are shared by faces generated on F
175     for ( size_t i = 0; i < nodes.size(); ++i )
176     {
177       SMDS_ElemIteratorPtr fIt = nodes[i]->GetInverseElementIterator(SMDSAbs_Face);
178       bool faceFound = false;
179       while ( !faceFound && fIt->more() )
180         faceFound = ( helper.GetSubShapeID() == fIt->next()->getshapeId() );
181       if ( !faceFound )
182         return error("The existing 1D mesh mismatches the generated 2D mesh");
183     }
184
185     const bool checkExisting = ( wires[ iWire ]->NbSegments() || helper.HasSeam() );
186
187     if ( _quadraticMesh )
188     {
189       for ( size_t i = 2; i < nodes.size(); i += 2 )
190       {
191         if ( checkExisting && meshDS->FindEdge( nodes[i-2], nodes[i], nodes[i-1]))
192           continue;
193         SMDS_MeshElement* e = meshDS->AddEdge( nodes[i-2], nodes[i], nodes[i-1] );
194         meshDS->SetMeshElementOnShape( e, nodes[i-1]->getshapeId() );
195       }
196     }
197     else
198     {
199       int edgeID = meshDS->ShapeToIndex( wires[ iWire ]->Edge(0) );
200       for ( size_t i = 1; i < nodes.size(); ++i )
201       {
202         if ( checkExisting && meshDS->FindEdge( nodes[i-1], nodes[i]))
203           continue;
204         SMDS_MeshElement* e = meshDS->AddEdge( nodes[i-1], nodes[i] );
205         if ( nodes[i-1]->getshapeId() != edgeID &&
206              nodes[i  ]->getshapeId() != edgeID )
207         {
208           edgeID = helper.GetMediumPos( nodes[i-1], nodes[i] ).first;
209           if ( edgeID < 1 ) edgeID = helper.GetSubShapeID();
210         }
211         meshDS->SetMeshElementOnShape( e, edgeID );
212       }
213     }
214   }
215
216   return true;
217 }
218
219 //=======================================================================
220 //function : Evaluate
221 //purpose  : 
222 //=======================================================================
223
224 bool StdMeshers_Projection_1D2D::Evaluate(SMESH_Mesh&         theMesh,
225                                           const TopoDS_Shape& theShape,
226                                           MapShapeNbElems&    aResMap)
227 {
228   if ( !StdMeshers_Projection_2D::Evaluate(theMesh,theShape,aResMap))
229     return false;
230
231   TopoDS_Shape srcFace = _sourceHypo->GetSourceFace();
232   SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
233   if ( !srcMesh ) srcMesh = & theMesh;
234   SMESH_subMesh* srcFaceSM = srcMesh->GetSubMesh( srcFace );
235
236   namespace SPU = StdMeshers_ProjectionUtils;
237   SPU::TShapeShapeMap shape2ShapeMap;
238   SPU::InitVertexAssociation( _sourceHypo, shape2ShapeMap );
239   if ( !SPU::FindSubShapeAssociation( theShape, &theMesh, srcFace, srcMesh, shape2ShapeMap))
240     return error(COMPERR_BAD_SHAPE,"Topology of source and target faces seems different" );
241
242   MapShapeNbElems srcResMap;
243   if ( !srcFaceSM->IsMeshComputed() )
244     _gen->Evaluate( *srcMesh, srcFace, srcResMap);
245
246   SMESH_subMeshIteratorPtr smIt = srcFaceSM->getDependsOnIterator(/*includeSelf=*/false,
247                                                                   /*complexShapeFirst=*/true);
248   while ( smIt->more() )
249   {
250     SMESH_subMesh* srcSM = smIt->next();
251     TopAbs_ShapeEnum shapeType = srcSM->GetSubShape().ShapeType();
252     if ( shapeType == TopAbs_EDGE )
253     {
254       std::vector<smIdType> aVec;
255       SMESHDS_SubMesh* srcSubMeshDS = srcSM->GetSubMeshDS();
256       if ( srcSubMeshDS && srcSubMeshDS->NbElements() )
257       {
258         aVec.resize(SMDSEntity_Last, 0);
259         SMDS_ElemIteratorPtr eIt = srcSubMeshDS->GetElements();
260         _quadraticMesh = ( eIt->more() && eIt->next()->IsQuadratic() );
261
262         aVec[SMDSEntity_Node] = srcSubMeshDS->NbNodes();
263         aVec[_quadraticMesh ? SMDSEntity_Quad_Edge : SMDSEntity_Edge] = srcSubMeshDS->NbElements();
264       }
265       else
266       {
267         if ( srcResMap.empty() )
268           if ( !_gen->Evaluate( *srcMesh, srcSM->GetSubShape(), srcResMap ))
269             return error(COMPERR_BAD_INPUT_MESH,"Source mesh not evaluatable");
270         aVec = srcResMap[ srcSM ];
271         if ( aVec.empty() )
272           return error(COMPERR_BAD_INPUT_MESH,"Source mesh is wrongly evaluated");
273       }
274       TopoDS_Shape tgtEdge = shape2ShapeMap( srcSM->GetSubShape(), /*isSrc=*/true  );
275       SMESH_subMesh* tgtSM = theMesh.GetSubMesh( tgtEdge );
276       aResMap.insert(std::make_pair(tgtSM,aVec));
277     }
278     if ( shapeType == TopAbs_VERTEX ) break;
279   }
280
281   return true;
282 }
283
284 //=======================================================================
285 //function : SetEventListener
286 //purpose  : Sets a default event listener to submesh of the source face.
287 //           faceSubMesh - submesh where algo is set
288 // After being set, event listener is notified on each event of a submesh.
289 // This method is called when a submesh gets HYP_OK algo_state.
290 // Arranges that CLEAN event is translated from source submesh to
291 // the faceSubMesh submesh.
292 //=======================================================================
293
294 void StdMeshers_Projection_1D2D::SetEventListener(SMESH_subMesh* faceSubMesh)
295 {
296   // set a listener of events on a source submesh
297   StdMeshers_Projection_2D::SetEventListener(faceSubMesh);
298
299   // set a listener to the target FACE submesh in order to update submehses
300   // of EDGEs according to events on the target FACE submesh
301   EventProparatorToEdges::Set( faceSubMesh );
302 }
303