Salome HOME
Join modifications from BR_Dev_For_4_0 tag V4_1_1.
[modules/smesh.git] / src / StdMeshers / StdMeshers_Projection_1D.cxx
1 //  SMESH SMESH : implementaion of SMESH idl descriptions
2 //
3 //  Copyright (C) 2003  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 //
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
30 #include "StdMeshers_Projection_1D.hxx"
31
32 #include "StdMeshers_ProjectionSource1D.hxx"
33 #include "StdMeshers_ProjectionUtils.hxx"
34
35 #include "SMDS_MeshNode.hxx"
36 #include "SMDS_MeshElement.hxx"
37 #include "SMESHDS_Mesh.hxx"
38 #include "SMESHDS_SubMesh.hxx"
39 #include "SMESH_Mesh.hxx"
40 #include "SMESH_MesherHelper.hxx"
41 #include "SMESH_subMesh.hxx"
42 #include "SMESH_subMeshEventListener.hxx"
43 #include "SMESH_Gen.hxx"
44 #include "SMESH_Comment.hxx"
45
46 #include <BRepAdaptor_Curve.hxx>
47 #include <BRep_Tool.hxx>
48 #include <GCPnts_AbscissaPoint.hxx>
49 #include <TopExp.hxx>
50 #include <TopoDS.hxx>
51 #include <gp_Pnt.hxx>
52 #include <TopTools_ListIteratorOfListOfShape.hxx>
53
54 #include "utilities.h"
55
56
57 using namespace std;
58
59 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
60
61 typedef StdMeshers_ProjectionUtils TAssocTool;
62
63 //=======================================================================
64 //function : StdMeshers_Projection_1D
65 //purpose  : 
66 //=======================================================================
67
68 StdMeshers_Projection_1D::StdMeshers_Projection_1D(int hypId, int studyId, SMESH_Gen* gen)
69   :SMESH_1D_Algo(hypId, studyId, 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 ( !TAssocTool::IsSubShape( _sourceHypo->GetSourceVertex(), srcMesh ) ||
134            !TAssocTool::IsSubShape( _sourceHypo->GetTargetVertex(), tgtMesh ) ||
135            !TAssocTool::IsSubShape( _sourceHypo->GetSourceVertex(),
136                                     _sourceHypo->GetSourceEdge() ))
137       {
138         aStatus = HYP_BAD_PARAMETER;
139         SCRUTE((TAssocTool::IsSubShape( _sourceHypo->GetSourceVertex(), srcMesh )));
140         SCRUTE((TAssocTool::IsSubShape( _sourceHypo->GetTargetVertex(), tgtMesh )));
141         SCRUTE((TAssocTool::IsSubShape( _sourceHypo->GetSourceVertex(),
142                                         _sourceHypo->GetSourceEdge() )));
143       }
144       // PAL16202
145       else 
146       {
147         bool isSub = TAssocTool::IsSubShape( _sourceHypo->GetTargetVertex(), aShape );
148         if ( !_sourceHypo->IsCompoundSource() ) {
149           if ( !isSub ) {
150             aStatus = HYP_BAD_PARAMETER;
151             SCRUTE((TAssocTool::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                  TAssocTool::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 ( !TAssocTool::IsSubShape( _sourceHypo->GetSourceEdge(), srcMesh ) ||
181          ( srcMesh == tgtMesh && aShape == _sourceHypo->GetSourceEdge() ))
182     {
183       aStatus = HYP_BAD_PARAMETER;
184       SCRUTE((TAssocTool::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 subshapes 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, tgtEdge );
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   if ( tgtMesh == srcMesh ) {
243     if ( !TAssocTool::MakeComputed( srcSubMesh ))
244       return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
245   }
246   else {
247     if ( !srcSubMesh->IsMeshComputed() )
248       return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
249   }
250   // -----------------------------------------------
251   // Find out nodes distribution on the source edge
252   // -----------------------------------------------
253
254   double srcLength = EdgeLength( srcEdge );
255   double tgtLength = EdgeLength( tgtEdge );
256   
257   vector< double > params; // sorted parameters of nodes on the source edge
258   if ( !SMESH_Algo::GetNodeParamOnEdge( srcMesh->GetMeshDS(), srcEdge, params ))
259     return error(COMPERR_BAD_INPUT_MESH,"Bad node parameters on the source edge");
260
261   int i, nbNodes = params.size();
262
263   vector< double > lengths( nbNodes - 1 ); // lengths of segments of the source edge
264   if ( srcLength > 0 )
265   {
266     BRepAdaptor_Curve curveAdaptor( srcEdge );
267     for ( i = 1; i < nbNodes; ++i )
268       lengths[ i-1 ] = GCPnts_AbscissaPoint::Length( curveAdaptor, params[i-1], params[i]);
269   }
270   else // degenerated source edge
271   {
272     for ( i = 1; i < nbNodes; ++i )
273       lengths[ i-1 ] = params[i] - params[i-1];
274     srcLength = params.back() - params[0];
275   }
276
277   bool reverse = ( srcV[0].IsSame( shape2ShapeMap( tgtV[1] )));
278   if ( tgtV[0].IsSame( tgtV[1] )) // case of closed edge
279     reverse = ( shape2ShapeMap( tgtEdge ).Orientation() == TopAbs_REVERSED );
280   if ( reverse ) // reverse lengths of segments
281     std::reverse( lengths.begin(), lengths.end() );
282
283   // ----------
284   // Make mesh
285   // ----------
286
287   // vector of target nodes
288   vector< const SMDS_MeshNode* > nodes ( nbNodes );
289
290   // Get the first and last nodes
291   nodes.front() = VertexNode( tgtV[0], meshDS );
292   nodes.back()  = VertexNode( tgtV[1], meshDS );
293   if ( !nodes.front() || !nodes.back() )
294     return error(COMPERR_BAD_INPUT_MESH,"No node on vertex");
295
296   // Compute parameters on the target edge and make internal nodes
297   // --------------------------------------------------------------
298
299   vector< double > tgtParams( nbNodes );
300
301   BRep_Tool::Range( tgtEdge, tgtParams.front(), tgtParams.back() );
302   if ( tgtLength <= 0 )
303     tgtLength = tgtParams.back() - tgtParams.front();
304   double dl = tgtLength / srcLength;
305
306   if ( tgtLength > 0 )
307   {
308     BRepAdaptor_Curve curveAdaptor( tgtEdge );
309
310     // compute params on internal nodes
311     for ( i = 1; i < nbNodes - 1; ++i )
312     {
313       // computes a point on a <curveAdaptor> at the given distance
314       // from the point at given parameter.
315       GCPnts_AbscissaPoint Discret( curveAdaptor, dl * lengths[ i-1 ], tgtParams[ i-1 ] );
316       if ( !Discret.IsDone() )
317         return error("GCPnts_AbscissaPoint failed");
318       tgtParams[ i ] = Discret.Parameter();
319     }
320     // make internal nodes 
321     for ( i = 1; i < nbNodes - 1; ++i )
322     {
323       gp_Pnt P = curveAdaptor.Value( tgtParams[ i ]);
324       SMDS_MeshNode* node = meshDS->AddNode(P.X(), P.Y(), P.Z());
325       meshDS->SetNodeOnEdge( node, tgtEdge, tgtParams[ i ]);
326       nodes[ i ] = node;
327     }
328   }
329   else // degenerated target edge
330   {
331     // compute params and make internal nodes
332     gp_Pnt P = BRep_Tool::Pnt( tgtV[0] );
333
334     for ( i = 1; i < nbNodes - 1; ++i )
335     {
336       SMDS_MeshNode* node = meshDS->AddNode(P.X(), P.Y(), P.Z());
337       tgtParams[ i ] = tgtParams[ i-1 ] + dl * lengths[ i-1 ];
338       meshDS->SetNodeOnEdge( node, tgtEdge, tgtParams[ i ]);
339       nodes[ i ] = node;
340     }
341   }
342
343   // Quadratic mesh?
344   // ----------------
345
346   bool quadratic = false;
347   SMDS_ElemIteratorPtr elemIt = srcSubMesh->GetSubMeshDS()->GetElements();
348   if ( elemIt->more() )
349     quadratic = elemIt->next()->IsQuadratic();
350   else {
351     SMDS_NodeIteratorPtr nodeIt = srcSubMesh->GetSubMeshDS()->GetNodes();
352     while ( nodeIt->more() && !quadratic )
353       quadratic = SMESH_MesherHelper::IsMedium( nodeIt->next() );
354   }
355   // enough nodes to make all edges quadratic?
356   if ( quadratic && ( nbNodes < 3 || ( nbNodes % 2 != 1 )))
357     return error(COMPERR_BAD_INPUT_MESH,
358                  SMESH_Comment("Wrong number of nodes to make quadratic mesh: ")<<nbNodes);
359
360   // Create edges
361   // -------------
362
363   SMDS_MeshElement* edge = 0;
364   int di = quadratic ? 2 : 1;
365   for ( i = di; i < nbNodes; i += di)
366   {
367     if ( quadratic )
368       edge = meshDS->AddEdge( nodes[i-2], nodes[i], nodes[i-1] );
369     else
370       edge = meshDS->AddEdge( nodes[i-1], nodes[i] );
371     meshDS->SetMeshElementOnShape(edge, tgtEdge );
372   }
373
374   return true;
375 }
376
377 //=============================================================================
378 /*!
379  * \brief Sets a default event listener to submesh of the source edge
380   * \param subMesh - submesh where algo is set
381  *
382  * This method is called when a submesh gets HYP_OK algo_state.
383  * After being set, event listener is notified on each event of a submesh.
384  * Arranges that CLEAN event is translated from source submesh to
385  * the submesh
386  */
387 //=============================================================================
388
389 void StdMeshers_Projection_1D::SetEventListener(SMESH_subMesh* subMesh)
390 {
391   TAssocTool::SetEventListener( subMesh,
392                                 _sourceHypo->GetSourceEdge(),
393                                 _sourceHypo->GetSourceMesh() );
394 }