Salome HOME
d25f7bb5a717c3ea440d6b6d6e86a1d9959cef51
[modules/paravis.git] / src / Plugins / ElevationSurface / plugin / ElevationSurfaceModule / vtkElevationSurfaceFilter.cxx
1 // Copyright (C) 2010-2021  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "vtkElevationSurfaceFilter.h"
21 #include "vtkInformation.h"
22 #include "vtkInformationVector.h"
23 #include "vtkObjectFactory.h"
24 #include "vtkElevationSurfaceFilter.h"
25 #include "vtkPolyData.h"
26 #include "vtkIdTypeArray.h"
27 #include "vtkPolyData.h"
28 #include "vtkUnstructuredGrid.h"
29 #include "vtkDemandDrivenPipeline.h"
30 #include "vtkStreamingDemandDrivenPipeline.h"
31 #include "vtkGenericCell.h"
32 #include "vtkSmartPointer.h"
33 #include "vtkPoints.h"
34 #include "vtkCellArray.h"
35 #include "vtkPointData.h"
36 #include "vtkCellData.h"
37
38 #include <math.h>
39
40 //vtkCxxRevisionMacro(vtkElevationSurfaceFilter, "$Revision$");
41 vtkStandardNewMacro(vtkElevationSurfaceFilter)
42
43 vtkElevationSurfaceFilter::vtkElevationSurfaceFilter()
44 {
45   this->SetNumberOfInputPorts(1);
46   this->SetNumberOfOutputPorts(1);
47
48   this->ScaleFactor = 0.5;
49   this->Direction[0] = 0.0;
50   this->Direction[1] = 0.0;
51   this->Direction[2] = 1.0;
52   this->AutoDetectDirection = 1;
53 }
54
55 vtkElevationSurfaceFilter::~vtkElevationSurfaceFilter()
56 {
57 }
58
59 //----------------------------------------------------------------------------
60 int vtkElevationSurfaceFilter::ProcessRequest(vtkInformation* request,
61                                          vtkInformationVector** inputVector,
62                                          vtkInformationVector* outputVector)
63 {
64   // generate the data
65   if(request->Has(vtkDemandDrivenPipeline::REQUEST_DATA()))
66     {
67     return this->RequestData(request, inputVector, outputVector);
68     }
69
70   if(request->Has(vtkStreamingDemandDrivenPipeline::REQUEST_UPDATE_EXTENT()))
71     {
72     return this->RequestUpdateExtent(request, inputVector, outputVector);
73     }
74
75   // execute information
76   if(request->Has(vtkDemandDrivenPipeline::REQUEST_INFORMATION()))
77     {
78     return this->RequestInformation(request, inputVector, outputVector);
79     }
80
81   return this->Superclass::ProcessRequest(request, inputVector, outputVector);
82 }
83
84 //----------------------------------------------------------------------------
85 int vtkElevationSurfaceFilter::FillOutputPortInformation(
86   int vtkNotUsed(port), vtkInformation* info)
87 {
88   // now add our info
89   info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkUnstructuredGrid");
90   return 1;
91 }
92
93 //----------------------------------------------------------------------------
94 int vtkElevationSurfaceFilter::FillInputPortInformation(
95   int vtkNotUsed(port), vtkInformation* info)
96 {
97   info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkPolyData");
98   return 1;
99 }
100
101 //----------------------------------------------------------------------------
102 int vtkElevationSurfaceFilter::RequestUpdateExtent(
103   vtkInformation* vtkNotUsed(request),
104   vtkInformationVector** inputVector,
105   vtkInformationVector* vtkNotUsed(outputVector))
106 {
107   int numInputPorts = this->GetNumberOfInputPorts();
108   for (int i=0; i<numInputPorts; i++)
109     {
110     int numInputConnections = this->GetNumberOfInputConnections(i);
111     for (int j=0; j<numInputConnections; j++)
112       {
113       vtkInformation* inputInfo = inputVector[i]->GetInformationObject(j);
114       inputInfo->Set(vtkStreamingDemandDrivenPipeline::EXACT_EXTENT(), 1);
115       }
116     }
117   return 1;
118 }
119
120 int vtkElevationSurfaceFilter::RequestInformation(vtkInformation * /*request*/,
121           vtkInformationVector ** /*input*/, vtkInformationVector * /*output*/)
122 {
123   return 1;
124 }
125
126 int vtkElevationSurfaceFilter::RequestData(vtkInformation * /*request*/,
127     vtkInformationVector **input, vtkInformationVector *output)
128 {
129   vtkPolyData *pdIn = vtkPolyData::SafeDownCast(
130       input[0]->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT()));
131
132   vtkUnstructuredGrid *usgOut = vtkUnstructuredGrid::SafeDownCast(
133       output->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT()));
134
135   vtkDataArray* array = this->GetInputArrayToProcess(0, input);
136
137   if(pdIn == NULL || array == NULL || usgOut == NULL
138      || array->GetNumberOfComponents() != 1)
139     {
140     vtkDebugMacro("vtkElevationSurfaceFilter no correctly configured");
141     return 1;
142     }
143
144   double dir[3];
145   if(this->AutoDetectDirection)
146     {
147     this->ComputeDirection(pdIn, dir);
148     }
149   else
150     {
151     dir[0] = this->Direction[0];
152     dir[1] = this->Direction[1];
153     dir[2] = this->Direction[2];
154     }
155
156   double len = dir[0]*dir[0] + dir[1]*dir[1] + dir[2]*dir[2];
157
158   if(len == 0)
159     {
160     dir[2] = 1;
161     len = 1.0;
162     }
163
164   len = sqrt(len);
165
166   dir[0] /= len;
167   dir[1] /= len;
168   dir[2] /= len;
169
170   dir[0] *= this->GetScaleFactor();
171   dir[1] *= this->GetScaleFactor();
172   dir[2] *= this->GetScaleFactor();
173
174   usgOut->Allocate(pdIn->GetNumberOfCells());
175
176   vtkSmartPointer<vtkPoints> newPts = vtkSmartPointer<vtkPoints>::New();
177   usgOut->SetPoints(newPts);
178
179   usgOut->GetPointData()->CopyAllocate(pdIn->GetPointData(),
180                                        2*pdIn->GetNumberOfPoints());
181   usgOut->GetCellData()->CopyAllocate(pdIn->GetCellData(),
182                                        pdIn->GetNumberOfCells());
183
184   vtkIdType ncell = pdIn->GetNumberOfCells();
185   vtkSmartPointer<vtkIdList> newIds = vtkSmartPointer<vtkIdList>::New();
186   vtkSmartPointer<vtkIdList> polyhedronIds = vtkSmartPointer<vtkIdList>::New();
187   for(vtkIdType cellId=0; cellId < ncell; cellId++)
188     {
189     vtkCell* cell = pdIn->GetCell(cellId);
190     if(cell->GetCellDimension() != 2)
191       continue;
192
193     unsigned char newCellType = VTK_EMPTY_CELL;
194     unsigned char oldCellType = (unsigned char)cell->GetCellType();
195     switch(oldCellType)
196       {
197       case VTK_TRIANGLE :
198         newCellType = VTK_WEDGE;
199         break;
200       case VTK_QUAD :
201         newCellType = VTK_HEXAHEDRON;
202         break;
203       case VTK_POLYGON :
204         newCellType = VTK_POLYHEDRON;
205       // default : add new cell types to extrude here
206       }
207     if(newCellType == VTK_EMPTY_CELL)
208       continue;
209
210     double cellScalar = array->GetTuple1(cellId);
211
212     vtkIdList* oldIds = cell->GetPointIds();
213     int oldPtsNumber = oldIds->GetNumberOfIds();
214     int newPtsNumber = oldPtsNumber * 2;
215     newIds->SetNumberOfIds(newPtsNumber);
216     double coords[VTK_CELL_SIZE*3];
217     for(int ptid = 0; ptid < oldPtsNumber; ptid++)
218       {
219       pdIn->GetPoint(oldIds->GetId(ptid), coords + 3*ptid);
220       }
221     for(int ptid = 0; ptid < oldPtsNumber; ptid++)
222       {
223       coords[(ptid+oldPtsNumber)*3+0] = coords[ptid*3+0] + cellScalar*dir[0];
224       coords[(ptid+oldPtsNumber)*3+1] = coords[ptid*3+1] + cellScalar*dir[1];
225       coords[(ptid+oldPtsNumber)*3+2] = coords[ptid*3+2] + cellScalar*dir[2];
226       }
227     for(int ptid=0; ptid<newPtsNumber; ptid++)
228       {
229       vtkIdType newId = newPts->InsertNextPoint(coords + 3*ptid);
230       newIds->SetId(ptid, newId);
231       usgOut->GetPointData()->CopyData(pdIn->GetPointData(),
232                                        oldIds->GetId(ptid % oldPtsNumber),
233                                        newIds->GetId(ptid));
234       }
235     vtkIdType newCellId;
236     if(newCellType == VTK_POLYHEDRON)
237       {
238       polyhedronIds->Initialize();
239       // in the polyhedron case, I will generate a quad for each edge
240       // of the input, and two capping faces
241       polyhedronIds->InsertNextId(2+oldPtsNumber);
242       // insert the bottom face
243       polyhedronIds->InsertNextId(oldPtsNumber);
244       for(int ptid = 0; ptid < oldPtsNumber; ptid++)
245         {
246         polyhedronIds->InsertNextId(newIds->GetId(ptid));
247         }
248       // insert the top face
249       polyhedronIds->InsertNextId(oldPtsNumber);
250       for(int ptid = oldPtsNumber; ptid < 2*oldPtsNumber; ptid++)
251         {
252         polyhedronIds->InsertNextId(newIds->GetId(ptid));
253         }
254       // insert the bording quads
255       for(int qid = 0; qid < oldPtsNumber; qid++)
256         {
257         polyhedronIds->InsertNextId(4);
258         polyhedronIds->InsertNextId(newIds->GetId(qid));
259         polyhedronIds->InsertNextId(newIds->GetId(qid+oldPtsNumber));
260         polyhedronIds->InsertNextId(newIds->GetId(qid+((oldPtsNumber+1)%oldPtsNumber)));
261         polyhedronIds->InsertNextId(newIds->GetId((qid+1)%oldPtsNumber));
262         }
263       newIds->Initialize();
264       for(int jj=0; jj<polyhedronIds->GetNumberOfIds(); jj++)
265         {
266         newIds->InsertNextId(polyhedronIds->GetId(jj));
267         }
268       }
269     newCellId = usgOut->InsertNextCell(newCellType, newIds);
270     usgOut->GetCellData()->CopyData(pdIn->GetCellData(),
271                                    cellId,
272                                    newCellId);
273     }
274
275   usgOut->GetFieldData()->ShallowCopy(pdIn->GetFieldData());
276
277   usgOut->Squeeze();
278
279   return 1;
280 }
281
282 void  vtkElevationSurfaceFilter::ComputeDirection(vtkPolyData* pdIn, double *outDir)
283 {
284   double tmp[2][3] = {{0, 0, 0}, {0, 0, 0}};
285   outDir[0] = outDir[1] = outDir[2] = 0;
286
287   vtkPoints* pts = pdIn->GetPoints();
288   vtkSmartPointer<vtkGenericCell> cell = vtkSmartPointer<vtkGenericCell>::New();
289
290   for(vtkIdType cellId = 0; cellId < pdIn->GetNumberOfCells(); cellId++)
291     {
292     pdIn->GetCell(cellId, cell);
293     if(cell->GetCellDimension() != 2)
294       continue;
295
296     vtkIdList* ptIds = cell->GetPointIds();
297     for(int i=0; i<ptIds->GetNumberOfIds(); i++)
298       {
299       vtkIdType firstId = ptIds->GetId(i);
300       vtkIdType secondId = ptIds->GetId((i+1)%ptIds->GetNumberOfIds());
301       pts->GetPoint(firstId, tmp[0]);
302       pts->GetPoint(secondId, tmp[1]);
303       outDir[0] += tmp[0][1]*tmp[1][2] - tmp[0][2]*tmp[1][1];
304       outDir[1] += tmp[0][2]*tmp[1][0] - tmp[0][0]*tmp[1][2];
305       outDir[2] += tmp[0][0]*tmp[1][1] - tmp[0][1]*tmp[1][0];
306       }
307     }
308 }
309
310 void vtkElevationSurfaceFilter::PrintSelf(ostream& os, vtkIndent indent)
311 {
312   this->Superclass::PrintSelf(os, indent);
313
314   os << indent << "ScaleFactor : " << this->ScaleFactor << endl;
315 }