]> SALOME platform Git repositories - modules/smesh.git/blob - src/StdMeshers/StdMeshers_Projection_1D.cxx
Salome HOME
0021825: Error in the example of "Projection Algorithms" in the user's guide
[modules/smesh.git] / src / StdMeshers / StdMeshers_Projection_1D.cxx
1 // Copyright (C) 2007-2012  CEA/DEN, EDF R&D, 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.
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
23 //  SMESH SMESH : implementaion of SMESH idl descriptions
24 // File      : StdMeshers_Projection_1D.cxx
25 // Module    : SMESH
26 // Created   : Fri Oct 20 11:37:07 2006
27 // Author    : Edward AGAPOV (eap)
28 //
29 #include "StdMeshers_Projection_1D.hxx"
30
31 #include "StdMeshers_ProjectionSource1D.hxx"
32 #include "StdMeshers_ProjectionUtils.hxx"
33
34 #include "SMDS_MeshNode.hxx"
35 #include "SMDS_MeshElement.hxx"
36 #include "SMESHDS_Mesh.hxx"
37 #include "SMESHDS_SubMesh.hxx"
38 #include "SMESH_Mesh.hxx"
39 #include "SMESH_MesherHelper.hxx"
40 #include "SMESH_subMesh.hxx"
41 #include "SMESH_subMeshEventListener.hxx"
42 #include "SMESH_Gen.hxx"
43 #include "SMESH_Comment.hxx"
44
45 #include <BRepAdaptor_Curve.hxx>
46 #include <BRep_Tool.hxx>
47 #include <GCPnts_AbscissaPoint.hxx>
48 #include <TopExp.hxx>
49 #include <TopoDS.hxx>
50 #include <gp_Pnt.hxx>
51 #include <TopTools_ListIteratorOfListOfShape.hxx>
52
53 #include "utilities.h"
54
55
56 using namespace std;
57
58 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
59
60 typedef StdMeshers_ProjectionUtils TAssocTool;
61
62 //=======================================================================
63 //function : StdMeshers_Projection_1D
64 //purpose  : 
65 //=======================================================================
66
67 StdMeshers_Projection_1D::StdMeshers_Projection_1D(int hypId, int studyId, SMESH_Gen* gen)
68   :SMESH_1D_Algo(hypId, studyId, gen)
69 {
70   _name = "Projection_1D";
71   _shapeType = (1 << TopAbs_EDGE);      // 1 bit per shape type
72
73   _compatibleHypothesis.push_back("ProjectionSource1D");
74   _sourceHypo = 0;
75 }
76
77 //================================================================================
78 /*!
79  * \brief Destructor
80  */
81 //================================================================================
82
83 StdMeshers_Projection_1D::~StdMeshers_Projection_1D()
84 {}
85
86 //=======================================================================
87 //function : CheckHypothesis
88 //purpose  : 
89 //=======================================================================
90
91 bool StdMeshers_Projection_1D::CheckHypothesis(SMESH_Mesh&                          aMesh,
92                                                const TopoDS_Shape&                  aShape,
93                                                SMESH_Hypothesis::Hypothesis_Status& aStatus)
94 {
95   _sourceHypo = 0;
96   list <const SMESHDS_Hypothesis * >::const_iterator itl;
97
98   const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
99   if ( hyps.size() == 0 )
100   {
101     aStatus = HYP_MISSING;
102     return false;  // can't work with no hypothesis
103   }
104
105   if ( hyps.size() > 1 )
106   {
107     aStatus = HYP_ALREADY_EXIST;
108     return false;
109   }
110
111   const SMESHDS_Hypothesis *theHyp = hyps.front();
112
113   string hypName = theHyp->GetName();
114
115   aStatus = HYP_OK;
116
117   if (hypName == "ProjectionSource1D")
118   {
119     _sourceHypo = static_cast<const StdMeshers_ProjectionSource1D *>(theHyp);
120
121     // Check hypo parameters
122
123     SMESH_Mesh* srcMesh = _sourceHypo->GetSourceMesh();
124     SMESH_Mesh* tgtMesh = & aMesh;
125     if ( !srcMesh )
126       srcMesh = tgtMesh;
127
128     // check vertices
129     if ( _sourceHypo->HasVertexAssociation() )
130     {
131       // source and target vertices
132       if ( !SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSourceVertex(), srcMesh ) ||
133            !SMESH_MesherHelper::IsSubShape( _sourceHypo->GetTargetVertex(), tgtMesh ) ||
134            !SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSourceVertex(),
135                                             _sourceHypo->GetSourceEdge() ))
136       {
137         aStatus = HYP_BAD_PARAMETER;
138         SCRUTE((SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSourceVertex(), srcMesh )));
139         SCRUTE((SMESH_MesherHelper::IsSubShape( _sourceHypo->GetTargetVertex(), tgtMesh )));
140         SCRUTE((SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSourceVertex(),
141                                                 _sourceHypo->GetSourceEdge() )));
142       }
143       // PAL16202
144       else
145       {
146         bool isSub = SMESH_MesherHelper::IsSubShape( _sourceHypo->GetTargetVertex(), aShape );
147         if ( !_sourceHypo->IsCompoundSource() ) {
148           if ( !isSub ) {
149             aStatus = HYP_BAD_PARAMETER;
150             SCRUTE((SMESH_MesherHelper::IsSubShape( _sourceHypo->GetTargetVertex(), aShape)));
151           }
152         }
153         else if ( isSub ) {
154           // is Ok provided that source vertex is shared only by one edge
155           // of the source group
156           TopoDS_Shape sharingEdge;
157           TopTools_ListIteratorOfListOfShape ancestIt
158             ( aMesh.GetAncestors( _sourceHypo->GetSourceVertex() ));
159           for ( ; ancestIt.More(); ancestIt.Next() )
160           {
161             const TopoDS_Shape& ancestor = ancestIt.Value();
162             if ( ancestor.ShapeType() == TopAbs_EDGE &&
163                  SMESH_MesherHelper::IsSubShape( ancestor, _sourceHypo->GetSourceEdge() ))
164             {
165               if ( sharingEdge.IsNull() || ancestor.IsSame( sharingEdge ))
166                 sharingEdge = ancestor;
167               else {
168                 // the second encountered
169                 aStatus = HYP_BAD_PARAMETER;
170                 MESSAGE("Source vertex is shared by several edges of a group");
171                 break;
172               }
173             }
174           }
175         }
176       }
177     }
178     // check source edge
179     if ( !SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSourceEdge(), srcMesh ) ||
180          ( srcMesh == tgtMesh && aShape == _sourceHypo->GetSourceEdge() ))
181     {
182       aStatus = HYP_BAD_PARAMETER;
183       SCRUTE((SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSourceEdge(), srcMesh )));
184       SCRUTE((srcMesh == tgtMesh));
185       SCRUTE(( aShape == _sourceHypo->GetSourceEdge() ));
186     }
187   }
188   else
189   {
190     aStatus = HYP_INCOMPATIBLE;
191   }
192   return ( aStatus == HYP_OK );
193 }
194
195 //=======================================================================
196 //function : Compute
197 //purpose  : 
198 //=======================================================================
199
200 bool StdMeshers_Projection_1D::Compute(SMESH_Mesh& theMesh, const TopoDS_Shape& theShape)
201 {
202   if ( !_sourceHypo )
203     return false;
204
205   SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh(); 
206   SMESH_Mesh * tgtMesh = & theMesh;
207   if ( !srcMesh )
208     srcMesh = tgtMesh;
209
210   SMESHDS_Mesh * meshDS = theMesh.GetMeshDS();
211
212   // ---------------------------
213   // Make sub-shapes association
214   // ---------------------------
215
216   TopoDS_Edge srcEdge, tgtEdge = TopoDS::Edge( theShape.Oriented(TopAbs_FORWARD));
217   TopoDS_Shape srcShape = _sourceHypo->GetSourceEdge().Oriented(TopAbs_FORWARD);
218
219   TAssocTool::TShapeShapeMap shape2ShapeMap;
220   TAssocTool::InitVertexAssociation( _sourceHypo, shape2ShapeMap );
221   if ( !TAssocTool::FindSubShapeAssociation( tgtEdge, tgtMesh, srcShape, srcMesh,
222                                              shape2ShapeMap) ||
223        !shape2ShapeMap.IsBound( tgtEdge ))
224     return error("Vertices association failed" );
225
226   srcEdge = TopoDS::Edge( shape2ShapeMap( tgtEdge ).Oriented(TopAbs_FORWARD));
227 //   cout << " srcEdge #" << srcMesh->GetMeshDS()->ShapeToIndex( srcEdge )
228 //        << " tgtEdge #" << tgtMesh->GetMeshDS()->ShapeToIndex( tgtEdge ) << endl;
229
230   TopoDS_Vertex tgtV[2], srcV[2];
231   TopExp::Vertices( tgtEdge, tgtV[0], tgtV[1] );
232   TopExp::Vertices( srcEdge, srcV[0], srcV[1] );
233
234   // ----------------------------------------------
235   // Assure that mesh on a source edge is computed
236   // ----------------------------------------------
237
238   SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( srcEdge );
239   //SMESH_subMesh* tgtSubMesh = tgtMesh->GetSubMesh( tgtEdge );
240
241   if ( tgtMesh == srcMesh ) {
242     if ( !TAssocTool::MakeComputed( srcSubMesh ))
243       return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
244   }
245   else {
246     if ( !srcSubMesh->IsMeshComputed() )
247       return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
248   }
249   // -----------------------------------------------
250   // Find out nodes distribution on the source edge
251   // -----------------------------------------------
252
253   double srcLength = EdgeLength( srcEdge );
254   double tgtLength = EdgeLength( tgtEdge );
255   
256   vector< double > params; // sorted parameters of nodes on the source edge
257   if ( !SMESH_Algo::GetNodeParamOnEdge( srcMesh->GetMeshDS(), srcEdge, params ))
258     return error(COMPERR_BAD_INPUT_MESH,"Bad node parameters on the source edge");
259
260   int i, nbNodes = params.size();
261
262   vector< double > lengths( nbNodes - 1 ); // lengths of segments of the source edge
263   if ( srcLength > 0 )
264   {
265     BRepAdaptor_Curve curveAdaptor( srcEdge );
266     for ( i = 1; i < nbNodes; ++i )
267       lengths[ i-1 ] = GCPnts_AbscissaPoint::Length( curveAdaptor, params[i-1], params[i]);
268   }
269   else // degenerated source edge
270   {
271     for ( i = 1; i < nbNodes; ++i )
272       lengths[ i-1 ] = params[i] - params[i-1];
273     srcLength = params.back() - params[0];
274   }
275
276   bool reverse = ( srcV[0].IsSame( shape2ShapeMap( tgtV[1] )));
277   if ( tgtV[0].IsSame( tgtV[1] )) // case of closed edge
278     reverse = ( shape2ShapeMap( tgtEdge ).Orientation() == TopAbs_REVERSED );
279   if ( reverse ) // reverse lengths of segments
280     std::reverse( lengths.begin(), lengths.end() );
281
282   // ----------
283   // Make mesh
284   // ----------
285
286   // vector of target nodes
287   vector< const SMDS_MeshNode* > nodes ( nbNodes );
288
289   // Get the first and last nodes
290   nodes.front() = VertexNode( tgtV[0], meshDS );
291   nodes.back()  = VertexNode( tgtV[1], meshDS );
292   if ( !nodes.front() || !nodes.back() )
293     return error(COMPERR_BAD_INPUT_MESH,"No node on vertex");
294
295   // Compute parameters on the target edge and make internal nodes
296   // --------------------------------------------------------------
297
298   vector< double > tgtParams( nbNodes );
299
300   BRep_Tool::Range( tgtEdge, tgtParams.front(), tgtParams.back() );
301   if ( tgtLength <= 0 )
302     tgtLength = tgtParams.back() - tgtParams.front();
303   double dl = tgtLength / srcLength;
304
305   if ( tgtLength > 0 )
306   {
307     BRepAdaptor_Curve curveAdaptor( tgtEdge );
308
309     // compute params on internal nodes
310     for ( i = 1; i < nbNodes - 1; ++i )
311     {
312       // computes a point on a <curveAdaptor> at the given distance
313       // from the point at given parameter.
314       GCPnts_AbscissaPoint Discret( curveAdaptor, dl * lengths[ i-1 ], tgtParams[ i-1 ] );
315       if ( !Discret.IsDone() )
316         return error("GCPnts_AbscissaPoint failed");
317       tgtParams[ i ] = Discret.Parameter();
318     }
319     // make internal nodes 
320     for ( i = 1; i < nbNodes - 1; ++i )
321     {
322       gp_Pnt P = curveAdaptor.Value( tgtParams[ i ]);
323       SMDS_MeshNode* node = meshDS->AddNode(P.X(), P.Y(), P.Z());
324       meshDS->SetNodeOnEdge( node, tgtEdge, tgtParams[ i ]);
325       nodes[ i ] = node;
326     }
327   }
328   else // degenerated target edge
329   {
330     // compute params and make internal nodes
331     gp_Pnt P = BRep_Tool::Pnt( tgtV[0] );
332
333     for ( i = 1; i < nbNodes - 1; ++i )
334     {
335       SMDS_MeshNode* node = meshDS->AddNode(P.X(), P.Y(), P.Z());
336       tgtParams[ i ] = tgtParams[ i-1 ] + dl * lengths[ i-1 ];
337       meshDS->SetNodeOnEdge( node, tgtEdge, tgtParams[ i ]);
338       nodes[ i ] = node;
339     }
340   }
341
342   // Quadratic mesh?
343   // ----------------
344
345   bool quadratic = false;
346   SMDS_ElemIteratorPtr elemIt = srcSubMesh->GetSubMeshDS()->GetElements();
347   if ( elemIt->more() )
348     quadratic = elemIt->next()->IsQuadratic();
349   else {
350     SMDS_NodeIteratorPtr nodeIt = srcSubMesh->GetSubMeshDS()->GetNodes();
351     while ( nodeIt->more() && !quadratic )
352       quadratic = SMESH_MesherHelper::IsMedium( nodeIt->next() );
353   }
354   // enough nodes to make all edges quadratic?
355   if ( quadratic && ( nbNodes < 3 || ( nbNodes % 2 != 1 )))
356     return error(COMPERR_BAD_INPUT_MESH,
357                  SMESH_Comment("Wrong number of nodes to make quadratic mesh: ")<<nbNodes);
358
359   // Create edges
360   // -------------
361
362   SMDS_MeshElement* edge = 0;
363   int di = quadratic ? 2 : 1;
364   for ( i = di; i < nbNodes; i += di)
365   {
366     if ( quadratic )
367       edge = meshDS->AddEdge( nodes[i-2], nodes[i], nodes[i-1] );
368     else
369       edge = meshDS->AddEdge( nodes[i-1], nodes[i] );
370     meshDS->SetMeshElementOnShape(edge, tgtEdge );
371   }
372
373   return true;
374 }
375
376
377 //=======================================================================
378 //function : Evaluate
379 //purpose  : 
380 //=======================================================================
381
382 bool StdMeshers_Projection_1D::Evaluate(SMESH_Mesh& theMesh,
383                                         const TopoDS_Shape& theShape,
384                                         MapShapeNbElems& aResMap)
385 {
386   if ( !_sourceHypo )
387     return false;
388
389   SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh(); 
390   SMESH_Mesh * tgtMesh = & theMesh;
391   if ( !srcMesh )
392     srcMesh = tgtMesh;
393
394   //SMESHDS_Mesh * meshDS = theMesh.GetMeshDS();
395
396   // ---------------------------
397   // Make sub-shapes association
398   // ---------------------------
399
400   TopoDS_Edge srcEdge, tgtEdge = TopoDS::Edge( theShape.Oriented(TopAbs_FORWARD));
401   TopoDS_Shape srcShape = _sourceHypo->GetSourceEdge().Oriented(TopAbs_FORWARD);
402
403   TAssocTool::TShapeShapeMap shape2ShapeMap;
404   TAssocTool::InitVertexAssociation( _sourceHypo, shape2ShapeMap );
405   if ( !TAssocTool::FindSubShapeAssociation( tgtEdge, tgtMesh, srcShape, srcMesh,
406                                              shape2ShapeMap) ||
407        !shape2ShapeMap.IsBound( tgtEdge ))
408     return error("Vertices association failed" );
409
410   srcEdge = TopoDS::Edge( shape2ShapeMap( tgtEdge ).Oriented(TopAbs_FORWARD));
411 //   cout << " srcEdge #" << srcMesh->GetMeshDS()->ShapeToIndex( srcEdge )
412 //        << " tgtEdge #" << tgtMesh->GetMeshDS()->ShapeToIndex( tgtEdge ) << endl;
413
414   TopoDS_Vertex tgtV[2], srcV[2];
415   TopExp::Vertices( tgtEdge, tgtV[0], tgtV[1] );
416   TopExp::Vertices( srcEdge, srcV[0], srcV[1] );
417
418   // ----------------------------------------------
419   // Assure that mesh on a source edge is computed
420   // ----------------------------------------------
421
422   SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( srcEdge );
423   //SMESH_subMesh* tgtSubMesh = tgtMesh->GetSubMesh( tgtEdge );
424
425   if ( tgtMesh == srcMesh ) {
426     if ( !TAssocTool::MakeComputed( srcSubMesh ))
427       return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
428   }
429   else {
430     if ( !srcSubMesh->IsMeshComputed() )
431       return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
432   }
433   // -----------------------------------------------
434   // Find out nodes distribution on the source edge
435   // -----------------------------------------------
436
437   //double srcLength = EdgeLength( srcEdge );
438   //double tgtLength = EdgeLength( tgtEdge );
439   
440   vector< double > params; // sorted parameters of nodes on the source edge
441   if ( !SMESH_Algo::GetNodeParamOnEdge( srcMesh->GetMeshDS(), srcEdge, params ))
442     return error(COMPERR_BAD_INPUT_MESH,"Bad node parameters on the source edge");
443
444   int nbNodes = params.size();
445
446   std::vector<int> aVec(SMDSEntity_Last);
447   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i] = 0;
448
449   aVec[SMDSEntity_Node] = nbNodes;
450
451   bool quadratic = false;
452   SMDS_ElemIteratorPtr elemIt = srcSubMesh->GetSubMeshDS()->GetElements();
453   if ( elemIt->more() )
454     quadratic = elemIt->next()->IsQuadratic();
455   if(quadratic)
456     aVec[SMDSEntity_Quad_Edge] = (nbNodes-1)/2;
457   else
458     aVec[SMDSEntity_Edge] = nbNodes - 1;
459
460   SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
461   aResMap.insert(std::make_pair(sm,aVec));
462
463   return true;
464 }
465
466
467 //=============================================================================
468 /*!
469  * \brief Sets a default event listener to submesh of the source edge
470   * \param subMesh - submesh where algo is set
471  *
472  * This method is called when a submesh gets HYP_OK algo_state.
473  * After being set, event listener is notified on each event of a submesh.
474  * Arranges that CLEAN event is translated from source submesh to
475  * the submesh
476  */
477 //=============================================================================
478
479 void StdMeshers_Projection_1D::SetEventListener(SMESH_subMesh* subMesh)
480 {
481   TAssocTool::SetEventListener( subMesh,
482                                 _sourceHypo->GetSourceEdge(),
483                                 _sourceHypo->GetSourceMesh() );
484 }