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