Salome HOME
Fix test case hang-up
[modules/paravis.git] / src / Plugins / ExtraxtFieldFilter / vtkExtractFieldFilter.cxx
1 // Copyright (C) 2010-2013  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.
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 "vtkExtractFieldFilter.h"
21
22 #include <vtkInformation.h>
23 #include <vtkInformationVector.h>
24 #include <vtkObjectFactory.h>
25 #include <vtkMultiBlockDataSet.h>
26 #include <vtkDataObjectTreeIterator.h>
27 #include <vtkFieldData.h>
28 #include <vtkStringArray.h>
29 #include <vtkDataSet.h>
30 #include <vtkPointData.h>
31 #include <vtkCellData.h>
32 #include <vtkDataArray.h>
33 #include <vtkStringArray.h>
34 #include <vtkCompositeDataToUnstructuredGridFilter.h>
35
36
37 #include <string.h>
38
39 using namespace std;
40
41
42
43 bool isContainsName(const std::list<std::string>& aList, const std::string& theName)
44 {
45         std::list<std::string>::const_iterator aIt;
46         for (aIt = aList.begin(); aIt != aList.end(); aIt++) {
47                 if ((*aIt).compare(theName) == 0)
48                         return true;
49         }
50         return false;
51 }
52
53
54 void appendIfNotExists(const std::list<std::string>& theSrc, std::list<std::string>& theDest)
55 {
56         std::list<std::string>::const_iterator aIt;
57         for (aIt = theSrc.begin(); aIt != theSrc.end(); aIt++) {
58                 if (!isContainsName(theDest, *aIt))
59                         theDest.push_back(*aIt);
60         }
61 }
62
63
64
65
66 vtkStandardNewMacro(vtkExtractFieldFilter);
67
68
69
70 vtkExtractFieldFilter::vtkExtractFieldFilter()
71 :vtkMultiBlockDataSetAlgorithm()
72 {
73         this->Field = NULL;
74         this->FieldList = vtkStringArray::New();
75 }
76
77
78 vtkExtractFieldFilter::~vtkExtractFieldFilter()
79 {
80         this->FieldList->Delete();
81         if (this->Field)
82                 delete [] this->Field;
83 }
84
85 //----------------------------------------------------------------------------
86 int vtkExtractFieldFilter::RequestData(vtkInformation* vtkNotUsed(request),
87                                                                            vtkInformationVector** theInputVector,
88                                                                            vtkInformationVector* theOutputVector)
89 {
90         // get the info objects
91         vtkMultiBlockDataSet* aInput = vtkMultiBlockDataSet::GetData(theInputVector[0], 0);
92         vtkMultiBlockDataSet* aOutput = vtkMultiBlockDataSet::GetData(theOutputVector, 0);
93
94         aOutput->CopyStructure(aInput);
95
96         // Copy selected blocks over to the output.
97         int i;
98         std::list<int> toDelList;
99         for (i = 0; i < aInput->GetNumberOfBlocks(); i++) {
100                 this->CopySubTree(i, aOutput, aInput, toDelList);
101         }
102         std::list<int>::const_reverse_iterator aIt;
103         for (aIt = toDelList.rbegin(); aIt != toDelList.rend(); ++aIt)
104                 aOutput->RemoveBlock(*aIt);
105         return 1;
106 }
107
108 //----------------------------------------------------------------------------
109 int vtkExtractFieldFilter::RequestInformation(vtkInformation* reqInfo,
110                                                                                           vtkInformationVector **theInputVector,
111                                                                                           vtkInformationVector *theOutputVector)
112 {
113         // get the info objects
114         vtkMultiBlockDataSet* aInput = vtkMultiBlockDataSet::GetData(theInputVector[0], 0);
115         vtkMultiBlockDataSet* aOutput = vtkMultiBlockDataSet::GetData(theOutputVector, 0);
116
117         vtkDataObjectTreeIterator* aIter = aInput->NewTreeIterator();
118         aIter->VisitOnlyLeavesOff();
119         int i = 0;
120         std::list<std::string> aList;
121         for (aIter->InitTraversal(); !aIter->IsDoneWithTraversal(); aIter->GoToNextItem()) {
122                 std::list<std::string> aSubList = this->GetListOfFields(aIter, aInput);
123                 appendIfNotExists(aSubList, aList);
124         }
125         this->FieldList->SetNumberOfValues(aList.size());
126         std::list<std::string>::const_iterator aIt;
127         i = 0;
128         for (aIt = aList.begin(); aIt != aList.end(); aIt++) {
129                 this->FieldList->SetValue(i, *aIt);
130                 i++;
131         }
132
133         return this->Superclass::RequestInformation(reqInfo, theInputVector, theOutputVector);
134 }
135
136 //----------------------------------------------------------------------------
137 void vtkExtractFieldFilter::CopySubTree(int theLoc,
138                                                                                 vtkMultiBlockDataSet* theOutput,
139                                                                                 vtkMultiBlockDataSet* theInput,
140                                                                                 std::list<int>& toDel)
141 {
142         vtkDataObject* aInputNode = theInput->GetBlock(theLoc);
143         if (aInputNode->IsA("vtkCompositeDataSet")) {
144                 vtkMultiBlockDataSet* aCInput = vtkMultiBlockDataSet::SafeDownCast(aInputNode);
145                 vtkMultiBlockDataSet* aCOutput = vtkMultiBlockDataSet::SafeDownCast(theOutput->GetBlock(theLoc));
146                 std::list<int> toDelList;
147                 int i;
148                 for (i = 0; i < aCInput->GetNumberOfBlocks(); i++) {
149                         this->CopySubTree(i, aCOutput, aCInput, toDelList);
150                 }
151                 std::list<int>::const_reverse_iterator aIt;
152                 for (aIt = toDelList.rbegin(); aIt != toDelList.rend(); ++aIt)
153                         aCOutput->RemoveBlock(*aIt);
154                 if (aCOutput->GetNumberOfBlocks() == 0)
155                         toDel.push_back(theLoc);
156         } else {
157                 if (IsToCopy(aInputNode)) {
158                         vtkDataObject* aClone = aInputNode->NewInstance();
159                         aClone->ShallowCopy(aInputNode);
160                         theOutput->SetBlock(theLoc, aClone);
161                         aClone->Delete();
162                 } else {
163                         toDel.push_back(theLoc);
164                 }
165         }
166 }
167
168 //----------------------------------------------------------------------------
169 std::list<std::string> vtkExtractFieldFilter::GetListOfFields(vtkDataObjectTreeIterator* theLoc, vtkMultiBlockDataSet* theInput)
170 {
171         std::list<std::string> aList;
172         vtkDataObject* aInputNode = theInput->GetDataSet(theLoc);
173         if (!aInputNode->IsA("vtkCompositeDataSet")) {
174                 std::list<std::string> aSubList = this->GetListOfFields(aInputNode);
175                 appendIfNotExists(aSubList, aList);
176         } else {
177                 vtkCompositeDataSet* aCInput = vtkCompositeDataSet::SafeDownCast(aInputNode);
178                 vtkCompositeDataIterator* aIter = aCInput->NewIterator();
179                 vtkDataObjectTreeIterator* aTreeIter = vtkDataObjectTreeIterator::SafeDownCast(aIter);
180                 if (aTreeIter) {
181                         aTreeIter->VisitOnlyLeavesOff();
182                 }
183                 for (aIter->InitTraversal(); !aIter->IsDoneWithTraversal(); aIter->GoToNextItem()) {
184                         vtkDataObject* aCurNode = aIter->GetCurrentDataObject();
185                         std::list<std::string> aSubList = this->GetListOfFields(aCurNode);
186                         appendIfNotExists(aSubList, aList);
187                 }
188                 aIter->Delete();
189         }
190         return aList;
191 }
192
193 //----------------------------------------------------------------------------
194 std::list<std::string> vtkExtractFieldFilter::GetListOfFields(vtkDataObject* theObject) const
195 {
196         std::list<std::string> aList;
197
198         if (theObject->IsA("vtkDataSet")) {
199                 vtkDataSet* aDataSet = vtkDataSet::SafeDownCast(theObject);
200                 vtkPointData* aPntData = aDataSet->GetPointData();
201                 int aNbArrays = aPntData->GetNumberOfArrays();
202                 for (int i = 0; i < aNbArrays; i++) {
203                         const char* aName = aPntData->GetArrayName(i);
204                         aList.push_back(aName);
205                 }
206                 vtkCellData* aCellData = aDataSet->GetCellData();
207                 aNbArrays = aCellData->GetNumberOfArrays();
208                 for (int i = 0; i < aNbArrays; i++) {
209                         const char* aName = aCellData->GetArrayName(i);
210                         aList.push_back(aName);
211                 }
212         }
213         return aList;
214 }
215
216
217 //----------------------------------------------------------------------------
218 bool vtkExtractFieldFilter::IsToCopy(vtkDataObject* theObject) const
219 {
220         if (this->Field == NULL)
221                 return true;
222
223         std::list<std::string> aList = this->GetListOfFields(theObject);
224
225         std::list<std::string>::const_iterator aIt;
226         std::string aTestStr = this->Field;
227         for (aIt = aList.begin(); aIt != aList.end(); ++aIt)
228                 if (aTestStr.compare(*aIt) == 0)
229                         return true;
230
231         return false;
232 }
233
234
235 //----------------------------------------------------------------------------
236 void vtkExtractFieldFilter::PrintSelf(ostream& os, vtkIndent indent)
237 {
238   this->Superclass::PrintSelf(os,indent);
239   os << indent << "Field name: " << this->Field << endl;
240   this->FieldList->PrintSelf(os, indent);
241 }
242