Skip to content
Snippets Groups Projects
Commit 7a18dbd8 authored by Sujin Philip's avatar Sujin Philip Committed by Kitware Robot
Browse files

Merge topic 'add-vtkExtractTimeSteps-filter'


eff36261 Add vtkExtractTimeSteps filter and test

Acked-by: default avatarKitware Robot <kwrobot@kitware.com>
Merge-request: !1939
parents 0b8bbcb4 eff36261
No related branches found
No related tags found
No related merge requests found
......@@ -20,6 +20,7 @@ set(Module_SRCS
vtkExtractSelection.cxx
vtkExtractTemporalFieldData.cxx
vtkExtractTensorComponents.cxx
vtkExtractTimeSteps.cxx
vtkExtractUnstructuredGrid.cxx
vtkExtractVectorComponents.cxx
......
......@@ -3,5 +3,6 @@ vtk_add_test_cxx(${vtk-module}CxxTests tests
TestExtractSelection.cxx
TestExtraction.cxx
TestExtractRectilinearGrid.cxx,NO_VALID,NO_DATA
TestExtractTimeSteps.cxx,NO_VALID
)
vtk_test_cxx_executable(${vtk-module}CxxTests tests)
/*=========================================================================
Program: Visualization Toolkit
Module: TestExtractTimeSteps.cxx
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
#include "vtkExtractTimeSteps.h"
#include "vtkExodusIIReader.h"
#include "vtkInformation.h"
#include "vtkNew.h"
#include "vtkStreamingDemandDrivenPipeline.h"
#include "vtkTestUtilities.h"
#include <cmath>
enum
{
TEST_PASSED_RETVAL = 0,
TEST_FAILED_RETVAL = 1
};
const double e = 1e-5;
int TestExtractTimeSteps(int argc, char *argv[])
{
const char *fname = vtkTestUtilities::ExpandDataFileName(argc, argv, "Data/can.ex2");
vtkNew<vtkExodusIIReader> reader;
reader->SetFileName(fname);
vtkNew<vtkExtractTimeSteps> extracter;
extracter->SetInputConnection(reader->GetOutputPort());
extracter->GenerateTimeStepIndices(0, 30, 5);
extracter->AddTimeStepIndex(30);
extracter->AddTimeStepIndex(35);
extracter->AddTimeStepIndex(30);
extracter->AddTimeStepIndex(40);
extracter->AddTimeStepIndex(43);
int numSteps = extracter->GetNumberOfTimeSteps();
if (numSteps != 10)
{
std::cout << "vtkExtractTimeSteps add time-steps failed" << std::endl;
return TEST_FAILED_RETVAL;
}
int tsteps[10];
extracter->GetTimeStepIndices(tsteps);
extracter->ClearTimeStepIndices();
extracter->SetTimeStepIndices(numSteps, tsteps);
extracter->Update();
double expected[10] = { 0.0000, 0.0005, 0.0010, 0.0015, 0.0020, 0.0025,
0.0030, 0.0035, 0.0040, 0.0043 };
double *result = NULL;
vtkInformation *info = extracter->GetOutputInformation(0);
if (info->Has(vtkStreamingDemandDrivenPipeline::TIME_STEPS()))
{
if (info->Length(vtkStreamingDemandDrivenPipeline::TIME_STEPS()) != 10)
{
std::cout << "got incorrect number of time steps" << std::endl;
return TEST_FAILED_RETVAL;
}
result = info->Get(vtkStreamingDemandDrivenPipeline::TIME_STEPS());
}
if (!result)
{
std::cout << "result has no time steps" << std::endl;
return TEST_FAILED_RETVAL;
}
else
{
for (int i = 0; i < 10; ++i)
{
if (std::abs(expected[i] - result[i]) > e)
{
std::cout << "extracted time steps values do not match" << std::endl;
return TEST_FAILED_RETVAL;
}
}
}
return TEST_PASSED_RETVAL;
}
......@@ -4,6 +4,7 @@ vtk_module(vtkFiltersExtraction
TEST_DEPENDS
vtkIOLegacy
vtkIOXML
vtkIOExodus
vtkRendering${VTK_RENDERING_BACKEND}
vtkTestingRendering
vtkInteractionStyle
......
/*=========================================================================
Program: Visualization Toolkit
Module: vtkExtractTimeSteps.cxx
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
#include "vtkExtractTimeSteps.h"
#include "vtkDataObject.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkObjectFactory.h"
#include "vtkStreamingDemandDrivenPipeline.h"
#include <algorithm>
#include <vector>
vtkStandardNewMacro(vtkExtractTimeSteps);
//----------------------------------------------------------------------------
void vtkExtractTimeSteps::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
int count = static_cast<int>(this->TimeStepIndices.size());
os << indent << "Number of Time Steps: " << count << std::endl;
if (count > 0)
{
std::set<int>::iterator it = this->TimeStepIndices.begin();
os << indent << "Time Step Indices: " << *it++;
for (int i = 1; i < std::min(count, 4); ++i)
{
os << ", " << *it++;
}
if (count > 9)
{
std::advance(it, count - 8);
os << ", ... ";
}
while (it != this->TimeStepIndices.end())
{
os << ", " << *it++;
}
os << std::endl;
}
}
//----------------------------------------------------------------------------
void vtkExtractTimeSteps::AddTimeStepIndex(int timeStepIndex)
{
if (this->TimeStepIndices.insert(timeStepIndex).second)
{
this->Modified();
}
}
void vtkExtractTimeSteps::SetTimeStepIndices(int count, const int *timeStepIndices)
{
this->TimeStepIndices.clear();
this->TimeStepIndices.insert(timeStepIndices, timeStepIndices + count);
this->Modified();
}
void vtkExtractTimeSteps::GetTimeStepIndices(int *timeStepIndices) const
{
std::copy(this->TimeStepIndices.begin(), this->TimeStepIndices.end(), timeStepIndices);
}
void vtkExtractTimeSteps::GenerateTimeStepIndices(int begin, int end, int step)
{
if (step != 0)
{
this->TimeStepIndices.clear();
for (int i = begin; i < end; i += step)
{
this->TimeStepIndices.insert(i);
}
this->Modified();
}
}
//----------------------------------------------------------------------------
int vtkExtractTimeSteps::RequestInformation(vtkInformation*,
vtkInformationVector **inputVector,
vtkInformationVector *outputVector)
{
// get the info objects
vtkInformation *outInfo = outputVector->GetInformationObject(0);
vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
if (!this->TimeStepIndices.empty() &&
inInfo->Has(vtkStreamingDemandDrivenPipeline::TIME_STEPS()))
{
double *inTimes =
inInfo->Get(vtkStreamingDemandDrivenPipeline::TIME_STEPS());
int numTimes =
inInfo->Length(vtkStreamingDemandDrivenPipeline::TIME_STEPS());
std::vector<double> outTimes;
for (std::set<int>::iterator it = this->TimeStepIndices.begin();
it != this->TimeStepIndices.end(); ++it)
{
if (*it >= 0 && *it < numTimes)
{
outTimes.push_back(inTimes[*it]);
}
}
if (!outTimes.empty())
{
outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(), &outTimes[0],
static_cast<int>(outTimes.size()));
double range[2] = { outTimes.front(), outTimes.back() };
outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(), range, 2);
}
}
return 1;
}
//----------------------------------------------------------------------------
int vtkExtractTimeSteps::RequestData(vtkInformation *,
vtkInformationVector **inputVector,
vtkInformationVector *outputVector)
{
vtkDataObject* inData = vtkDataObject::GetData(inputVector[0], 0);
vtkDataObject* outData = vtkDataObject::GetData(outputVector, 0);
if (inData && outData)
{
outData->ShallowCopy(inData);
}
return 1;
}
/*=========================================================================
Program: Visualization Toolkit
Module: vtkExtractTimeSteps.h
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
// .NAME vtkExtractTimeSteps - extract specific time-steps from dataset
// .SECTION Description
// vtkExtractTimeSteps extracts the specified time steps from the input dataset.
// The timesteps to be extracted are specified by their indices. If no
// time step is specified, all of the input time steps are extracted.
// This filter is useful when one wants to work with only a sub-set of the input
// time steps.
#ifndef vtkExtractTimeSteps_h
#define vtkExtractTimeSteps_h
#include "vtkFiltersExtractionModule.h" // for export macro
#include "vtkPassInputTypeAlgorithm.h"
#include <set> // for time step indices
class VTKFILTERSEXTRACTION_EXPORT vtkExtractTimeSteps : public vtkPassInputTypeAlgorithm
{
public:
vtkTypeMacro(vtkExtractTimeSteps, vtkPassInputTypeAlgorithm);
void PrintSelf(ostream& os, vtkIndent indent);
static vtkExtractTimeSteps *New();
// Description:
// Get the number of time steps that will be extracted
int GetNumberOfTimeSteps() const
{
return static_cast<int>(this->TimeStepIndices.size());
}
// Description:
// Add a time step index. Not added if the index already exists.
void AddTimeStepIndex(int timeStepIndex);
// Description:
// Get/Set an array of time step indices. For the Get function,
// timeStepIndices should be big enough for GetNumberOfTimeSteps() values.
void SetTimeStepIndices(int count, const int *timeStepIndices);
void GetTimeStepIndices(int *timeStepIndices) const;
// Description:
// Generate a range of indices in [begin, end) with a step size of 'step'
void GenerateTimeStepIndices(int begin, int end, int step);
// Description:
// Clear the time step indices
void ClearTimeStepIndices()
{
this->TimeStepIndices.clear();
this->Modified();
}
protected:
vtkExtractTimeSteps() {};
~vtkExtractTimeSteps() {};
virtual int RequestData(vtkInformation *, vtkInformationVector **,
vtkInformationVector *);
virtual int RequestInformation(vtkInformation *, vtkInformationVector **,
vtkInformationVector *);
std::set<int> TimeStepIndices;
private:
vtkExtractTimeSteps(const vtkExtractTimeSteps&) VTK_DELETE_FUNCTION;
void operator=(const vtkExtractTimeSteps&) VTK_DELETE_FUNCTION;
};
#endif // vtkExtractTimeSteps_h
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment