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