Salome HOME
Merge branch 'V9_10_BR'
[modules/paravis.git] / src / Plugins / StaticMesh / plugin / StaticMeshModule / vtkStaticDataSetSurfaceFilter.cxx
1 /*=========================================================================
2
3   Program:   Visualization Toolkit
4   Module:    vtkStaticDataSetSurfaceFilter.cxx
5
6   Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7   All rights reserved.
8   See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10      This software is distributed WITHOUT ANY WARRANTY; without even
11      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12      PURPOSE.  See the above copyright notice for more information.
13
14 =========================================================================*/
15 #include "vtkStaticDataSetSurfaceFilter.h"
16
17 #include <vtkCellData.h>
18 #include <vtkInformation.h>
19 #include <vtkInformationVector.h>
20 #include <vtkObjectFactory.h>
21 #include <vtkPointData.h>
22 #include <vtkUnstructuredGrid.h>
23 #include <vtksys/SystemTools.hxx>
24
25 vtkStandardNewMacro(vtkStaticDataSetSurfaceFilter);
26
27 //-----------------------------------------------------------------------------
28 int vtkStaticDataSetSurfaceFilter::UnstructuredGridExecute(vtkDataSet* input, vtkPolyData* output)
29 {
30   vtkUnstructuredGrid* inputUG = vtkUnstructuredGrid::SafeDownCast(input);
31   if (!inputUG)
32   {
33     // Rely on superclass for any input which is not a vtkUnstructuredGrid
34     return this->Superclass::UnstructuredGridExecute(input, output);
35   }
36
37   // Check if cache is still valid
38   if (this->InputMeshTime == inputUG->GetMeshMTime() && this->FilterMTime == this->GetMTime())
39   {
40     if (vtksys::SystemTools::HasEnv("VTK_DEBUG_STATIC_MESH"))
41     {
42       vtkWarningMacro("Using static mesh cache");
43     }
44
45     // Use cache as base
46     output->ShallowCopy(this->Cache.Get());
47
48     // Recover original ids
49     vtkPointData* outPD = output->GetPointData();
50     vtkCellData* outCD = output->GetCellData();
51     vtkIdTypeArray* origPointArray =
52       vtkIdTypeArray::SafeDownCast(outPD->GetArray(this->GetOriginalPointIdsName()));
53     vtkIdTypeArray* origCellArray =
54       vtkIdTypeArray::SafeDownCast(outCD->GetArray(this->GetOriginalCellIdsName()));
55     if (!origPointArray || !origCellArray)
56     {
57       vtkErrorMacro(
58         "OriginalPointIds or OriginalCellIds are missing, cannot use static mesh cache");
59       return this->Superclass::UnstructuredGridExecute(input, output);
60     }
61
62     // Recover input point and cell data
63     vtkPointData* inPD = input->GetPointData();
64     vtkCellData* inCD = input->GetCellData();
65
66     // Update output point data
67     vtkNew<vtkIdList> pointIds;
68     const vtkIdType origPANbTuples = origPointArray->GetNumberOfTuples();
69     pointIds->SetNumberOfIds(origPANbTuples);
70     for (vtkIdType i = 0; i < origPANbTuples; i++)
71     {
72       pointIds->SetId(i, origPointArray->GetTuple1(i));
73     }
74
75     // Remove array that have disappeared from input
76     for (int iArr = outPD->GetNumberOfArrays() - 1; iArr >= 0; iArr--)
77     {
78       outPD->RemoveArray(iArr);
79     }
80
81     // Update or create arrays present in input
82     for (int iArr = 0; iArr < inPD->GetNumberOfArrays(); iArr++)
83     {
84       vtkAbstractArray* outArr = outPD->GetAbstractArray(inPD->GetArrayName(iArr));
85       if (outArr)
86       {
87         inPD->GetAbstractArray(iArr)->GetTuples(pointIds.Get(), outArr);
88       }
89       else
90       {
91         // New array in input, create it in output
92         vtkAbstractArray* inArr = inPD->GetAbstractArray(iArr);
93         outArr = inArr->NewInstance();
94         outArr->SetName(inArr->GetName());
95         outArr->SetNumberOfComponents(inArr->GetNumberOfComponents());
96         outArr->SetNumberOfTuples(output->GetNumberOfPoints());
97         inArr->GetTuples(pointIds.Get(), outArr);
98         outPD->AddArray(outArr);
99         outArr->Delete();
100       }
101     }
102
103     // Update output cell data
104     vtkNew<vtkIdList> cellIds;
105     const vtkIdType origCANbTuples = origCellArray->GetNumberOfTuples();
106     cellIds->SetNumberOfIds(origCANbTuples);
107     for (vtkIdType i = 0; i < origCANbTuples; i++)
108     {
109       cellIds->SetId(i, origCellArray->GetTuple1(i));
110     }
111
112     // EDF26573 : Initially At the begining Remove array that have disappeared from input
113     // But after investigation ShallowCopy of Cache leads to a same array on input/output and leads to conflicts
114     // on inArr->GetTuples(cellIds.Get(), outArr); operation. So here a deepCopy of arrays is done inconditionnaly
115     for (int iArr = outCD->GetNumberOfArrays() - 1; iArr >= 0; iArr--)
116     {
117       outCD->RemoveArray(iArr);
118     }
119
120     for (int iArr = 0; iArr < inCD->GetNumberOfArrays(); iArr++)
121     {
122       vtkAbstractArray* outArr = outCD->GetAbstractArray(inCD->GetArrayName(iArr));
123       if(outArr)
124       {// EDF26573 : normally this condition is never fetched
125         inCD->GetAbstractArray(iArr)->GetTuples(cellIds.Get(), outArr);
126       }
127       else
128       {
129         // New array in input, create it in output
130         vtkAbstractArray* inArr = inCD->GetAbstractArray(iArr);
131         outArr = inArr->NewInstance();
132         outArr->SetName(inArr->GetName());
133         outArr->SetNumberOfComponents(inArr->GetNumberOfComponents());
134         outArr->SetNumberOfTuples(output->GetNumberOfCells());
135         inArr->GetTuples(cellIds.Get(), outArr);
136         outCD->AddArray(outArr);
137         outArr->Delete();
138       }
139     }
140
141     // Update output field data
142     output->GetFieldData()->ShallowCopy(input->GetFieldData());
143     return 1;
144   }
145   else
146   {
147     // Cache is not valid, Execute supercall algorithm
148     if (vtksys::SystemTools::HasEnv("VTK_DEBUG_STATIC_MESH"))
149     {
150       vtkWarningMacro("Building static mesh cache");
151     }
152
153     int ret = this->Superclass::UnstructuredGridExecute(input, output);
154
155     // Update the cache with superclass output
156     this->Cache->ShallowCopy(output);
157     this->InputMeshTime = inputUG->GetMeshMTime();
158     this->FilterMTime = this->GetMTime();
159     return ret;
160   }
161 }
162
163 //----------------------------------------------------------------------------
164 void vtkStaticDataSetSurfaceFilter::PrintSelf(ostream& os, vtkIndent indent)
165 {
166   this->Superclass::PrintSelf(os, indent);
167   os << indent << "Cache: " << this->Cache << endl;
168   os << indent << "Input Mesh Time: " << this->InputMeshTime << endl;
169   os << indent << "Filter mTime: " << this->FilterMTime << endl;
170 }