1 //=============================================================================
2 //
3 //  Copyright (c) Kitware, Inc.
4 //  All rights reserved.
5 //  See LICENSE.txt for details.
6 //
7 //  This software is distributed WITHOUT ANY WARRANTY; without even
8 //  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
9 //  PURPOSE.  See the above copyright notice for more information.
10 //
11 //  Copyright 2012 Sandia Corporation.
12 //  Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
13 //  the U.S. Government retains certain rights in this software.
14 //
15 //=============================================================================
16 #include "vtkmProbe.h"
17 
18 #include "vtkCellData.h"
19 #include "vtkDataSet.h"
20 #include "vtkExecutive.h"
21 #include "vtkImageData.h"
22 #include "vtkInformation.h"
23 #include "vtkInformationVector.h"
24 #include "vtkObjectFactory.h"
25 #include "vtkStreamingDemandDrivenPipeline.h"
26 #include "vtkPointData.h"
27 
28 #include "vtkmlib/ArrayConverters.h"
29 #include "vtkmlib/DataSetConverters.h"
30 #include "vtkmlib/Storage.h"
31 
32 #include "vtkmCellSetExplicit.h"
33 #include "vtkmCellSetSingleType.h"
34 #include "vtkmFilterPolicy.h"
35 
36 #include "vtkm/filter/Probe.h"
37 
vtkStandardNewMacro(vtkmProbe)38 vtkStandardNewMacro(vtkmProbe)
39 
40 //------------------------------------------------------------------------------
41 vtkmProbe::vtkmProbe()
42 {
43   this->SetNumberOfInputPorts(2);
44   this->PassCellArrays = false;
45   this->PassPointArrays = false;
46   this->PassFieldArrays = true;
47   this->ValidPointMaskArrayName = "vtkValidPointMask";
48   this->ValidCellMaskArrayName = "vtkValidCellMask";
49 }
50 
51 //------------------------------------------------------------------------------
SetSourceData(vtkDataObject * input)52 void vtkmProbe::SetSourceData(vtkDataObject* input)
53 {
54   this->SetInputData(1, input);
55 }
56 
57 //------------------------------------------------------------------------------
GetSource()58 vtkDataObject* vtkmProbe::GetSource()
59 {
60   if (this->GetNumberOfInputConnections(1) < 1)
61   {
62      return nullptr;
63   }
64   return this->GetExecutive()->GetInputData(1, 0);
65 }
66 
67 //------------------------------------------------------------------------------
SetSourceConnection(vtkAlgorithmOutput * algOutput)68 void vtkmProbe::SetSourceConnection(vtkAlgorithmOutput *algOutput)
69 {
70   this->SetInputConnection(1, algOutput);
71 }
72 
73 //------------------------------------------------------------------------------
RequestData(vtkInformation * vtkNotUsed (request),vtkInformationVector ** inputVector,vtkInformationVector * outputVector)74 int vtkmProbe::RequestData(vtkInformation* vtkNotUsed(request),
75                                  vtkInformationVector** inputVector,
76                                  vtkInformationVector* outputVector)
77 {
78   // Get the info objects
79   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
80   vtkInformation *sourceInfo = inputVector[1]->GetInformationObject(0);
81   vtkInformation *outInfo = outputVector->GetInformationObject(0);
82 
83   // Get the input and output
84   vtkDataSet* input = vtkDataSet::SafeDownCast(
85                  inInfo->Get(vtkDataSet::DATA_OBJECT()));
86   vtkDataSet* source = vtkDataSet::SafeDownCast(
87                  sourceInfo->Get(vtkDataSet::DATA_OBJECT()));
88   vtkDataSet* output = vtkDataSet::SafeDownCast(
89                  outInfo->Get(vtkDataSet::DATA_OBJECT()));
90 
91   // Copy the input to the output as a starting point
92   output->CopyStructure(input);
93 
94   try
95   {
96     // Convert the input dataset to a vtkm::cont::DataSet
97     vtkm::cont::DataSet in = tovtkm::Convert(input);
98     // VTK-m's probe filter requires the source to have at least a cellSet.
99     vtkm::cont::DataSet so = tovtkm::Convert(source, tovtkm::FieldsFlag::PointsAndCells);
100     if (!so.GetNumberOfCellSets())
101     {
102       vtkErrorMacro(<< "The source geometry does not have any cell set,"
103                        "aborting vtkmProbe filter");
104       return 0;
105     }
106 
107     vtkmInputFilterPolicy policy;
108     vtkm::filter::Probe probe;
109     // The input in VTK is the geometry in VTKM and the source in VTK is the input
110     // in VTKM.
111     probe.SetGeometry(in);
112 
113     auto result = probe.Execute(so, policy);
114     for (vtkm::Id i=0; i < result.GetNumberOfFields(); i++)
115     {
116       const vtkm::cont::Field& field = result.GetField(i);
117       vtkDataArray* fieldArray = fromvtkm::Convert(field);
118       if (field.GetAssociation() == vtkm::cont::Field::Association::POINTS)
119       {
120         if (strcmp(fieldArray->GetName(), "HIDDEN") == 0)
121         {
122           fieldArray->SetName(this->ValidPointMaskArrayName.c_str());
123         }
124         output->GetPointData()->AddArray(fieldArray);
125       }
126       else if (field.GetAssociation() == vtkm::cont::Field::Association::CELL_SET)
127       {
128         if (strcmp(fieldArray->GetName(), "HIDDEN") == 0)
129         {
130           fieldArray->SetName(this->ValidCellMaskArrayName.c_str());
131         }
132         output->GetCellData()->AddArray(fieldArray);
133       }
134       fieldArray->FastDelete();
135     }
136   }
137   catch (const vtkm::cont::Error& e)
138   {
139     vtkErrorMacro(<< "VTK-m error: " << e.GetMessage());
140     return 0;
141   }
142 
143   this->PassAttributeData( input, source, output);
144 
145   return 1;
146 }
147 
148 //------------------------------------------------------------------------------
RequestInformation(vtkInformation * vtkNotUsed (request),vtkInformationVector ** inputVector,vtkInformationVector * outputVector)149 int vtkmProbe::RequestInformation(vtkInformation* vtkNotUsed(request),
150                                   vtkInformationVector** inputVector,
151                                   vtkInformationVector* outputVector)
152 {
153   // Update the whole extent in the output
154   vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
155   vtkInformation* sourceInfo = inputVector[1]->GetInformationObject(0);
156   vtkInformation* outInfo = outputVector->GetInformationObject(0);
157   int wholeExtent[6];
158   if (inInfo && outInfo)
159   {
160     outInfo->CopyEntry(sourceInfo,
161                        vtkStreamingDemandDrivenPipeline::TIME_STEPS());
162     outInfo->CopyEntry(sourceInfo,
163                        vtkStreamingDemandDrivenPipeline::TIME_RANGE());
164     inInfo->Get(vtkStreamingDemandDrivenPipeline::WHOLE_EXTENT(), wholeExtent);
165     outInfo->Set(vtkStreamingDemandDrivenPipeline::WHOLE_EXTENT(), wholeExtent, 6);
166 
167   // Make sure that the scalar type and number of components
168   // are propagated from the source not the input.
169   if (vtkImageData::HasScalarType(sourceInfo))
170   {
171     vtkImageData::SetScalarType(vtkImageData::GetScalarType(sourceInfo),
172                                 outInfo);
173   }
174   if (vtkImageData::HasNumberOfScalarComponents(sourceInfo))
175   {
176     vtkImageData::SetNumberOfScalarComponents(
177       vtkImageData::GetNumberOfScalarComponents(sourceInfo),
178       outInfo);
179   }
180     return 1;
181   }
182   vtkErrorMacro("Missing input or output info!");
183   return 0;
184 }
185 
186 //------------------------------------------------------------------------------
RequestUpdateExtent(vtkInformation * vtkNotUsed (request),vtkInformationVector ** inputVector,vtkInformationVector * outputVector)187 int vtkmProbe::RequestUpdateExtent(vtkInformation* vtkNotUsed(request),
188                                    vtkInformationVector** inputVector,
189                                    vtkInformationVector* outputVector)
190 {
191   vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
192   vtkInformation* sourceInfo = inputVector[1]->GetInformationObject(0);
193   vtkInformation* outInfo = outputVector->GetInformationObject(0);
194 
195   if (inInfo && outInfo)
196   { // Source's update exetent should be independent of the resampling extent
197     inInfo->Set(vtkStreamingDemandDrivenPipeline::EXACT_EXTENT(), 1);
198     sourceInfo->Remove(vtkStreamingDemandDrivenPipeline::UPDATE_EXTENT());
199     if (sourceInfo->Has(vtkStreamingDemandDrivenPipeline::WHOLE_EXTENT()))
200     {
201       sourceInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_EXTENT(),
202         sourceInfo->Get(vtkStreamingDemandDrivenPipeline::WHOLE_EXTENT()), 6);
203     }
204     return 1;
205   }
206   vtkErrorMacro("Missing input or output info!");
207   return 0;
208 
209 }
210 
211 //------------------------------------------------------------------------------
PassAttributeData(vtkDataSet * input,vtkDataObject * vtkNotUsed (source),vtkDataSet * output)212 void vtkmProbe::PassAttributeData(vtkDataSet* input,
213                                   vtkDataObject* vtkNotUsed(source),
214                                   vtkDataSet* output)
215 {
216   if (this->PassPointArrays)
217   { // Copy point data arrays
218     int numPtArrays = input->GetPointData()->GetNumberOfArrays();
219     for (int i=0; i < numPtArrays; i++)
220     {
221       vtkDataArray* da = input->GetPointData()->GetArray(i);
222       if (da && !output->GetPointData()->HasArray(da->GetName()))
223       {
224         output->GetPointData()->AddArray(da);
225       }
226     }
227 
228     // Set active attributes in the output to the active attributes in the input
229     for (int i = 0; i < vtkDataSetAttributes::NUM_ATTRIBUTES; ++i)
230     {
231       vtkAbstractArray* da = input->GetPointData()->GetAttribute(i);
232       if (da && da->GetName() && !output->GetPointData()->GetAttribute(i))
233       {
234         output->GetPointData()->SetAttribute(da, i);
235       }
236     }
237   }
238 
239   // copy cell data arrays
240   if (this->PassCellArrays)
241   {
242     int numCellArrays = input->GetCellData()->GetNumberOfArrays();
243     for (int i=0; i<numCellArrays; ++i)
244     {
245       vtkDataArray *da = input->GetCellData()->GetArray(i);
246       if (!output->GetCellData()->HasArray(da->GetName()))
247       {
248         output->GetCellData()->AddArray(da);
249       }
250     }
251 
252     // Set active attributes in the output to the active attributes in the input
253     for (int i = 0; i < vtkDataSetAttributes::NUM_ATTRIBUTES; ++i)
254     {
255       vtkAbstractArray* da = input->GetCellData()->GetAttribute(i);
256       if (da && da->GetName() && !output->GetCellData()->GetAttribute(i))
257       {
258         output->GetCellData()->SetAttribute(da, i);
259       }
260     }
261   }
262 
263   if (this->PassFieldArrays)
264   {
265     // nothing to do, vtkDemandDrivenPipeline takes care of that.
266   }
267   else
268   {
269     output->GetFieldData()->Initialize();
270   }
271 }
272 
273 //------------------------------------------------------------------------------
PrintSelf(ostream & os,vtkIndent indent)274 void vtkmProbe::PrintSelf(ostream& os, vtkIndent indent)
275 {
276   this->Superclass::PrintSelf(os, indent);
277   os << indent << "PassPointArrays: " << this->PassPointArrays << "\n";
278   os << indent << "PassCellArrays: " << this->PassCellArrays << "\n";
279   os << indent << "PassFieldArray: " << this->PassFieldArrays << "\n";
280 }
281