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