Skip to content
Snippets Groups Projects
Commit 7fe30aa1 authored by Shawn Waldon's avatar Shawn Waldon Committed by Kitware Robot
Browse files

Merge topic 'time-series-to-multiblock'


4497fc65 Add a filter that groups a time series into a multiblock dataset

Acked-by: default avatarKitware Robot <kwrobot@kitware.com>
Reviewed-by: default avatarUtkarsh Ayachit <utkarsh.ayachit@kitware.com>
Merge-request: !1894
parents 2d305ac9 4497fc65
No related branches found
No related tags found
No related merge requests found
......@@ -45,6 +45,7 @@ set(Module_SRCS
vtkLinkEdgels.cxx
vtkMergeCells.cxx
vtkMultiBlockDataGroupFilter.cxx
vtkMultiBlockFromTimeSeriesFilter.cxx
vtkMultiBlockMergeFilter.cxx
vtkMultiThreshold.cxx
vtkOBBDicer.cxx
......
e9739ba863a3eb522ddee2ab32652d4e
......@@ -10,6 +10,7 @@ vtk_add_test_python(
TestDiscreteMarchingCubesAdjacentScalars.py
TestGraphLayoutFilter.py
TestMultiBlockStreamer.py
TestMultiBlockFromTimeSeries.py
TestPointConnectivityFilter.py
TestRectilinearGridToTetrahedra.py
TestSplineFilter.py
......
#!/usr/bin/python
import vtk
from vtk.test import Testing
from vtk.util.vtkAlgorithm import VTKPythonAlgorithmBase
class MovingSphereSource(VTKPythonAlgorithmBase):
def __init__(self):
VTKPythonAlgorithmBase.__init__(self,
nInputPorts=0,
nOutputPorts=1, outputType='vtkPolyData')
def RequestInformation(self, request, inInfo, outInfo):
info = outInfo.GetInformationObject(0)
t = range(0, 10)
info.Set(vtk.vtkStreamingDemandDrivenPipeline.TIME_STEPS(), t, len(t))
info.Set(vtk.vtkStreamingDemandDrivenPipeline.TIME_RANGE(), [t[0], t[-1]], 2)
return 1
def RequestData(self, request, inInfo, outInfo):
info = outInfo.GetInformationObject(0)
output = vtk.vtkPolyData.GetData(outInfo)
t = info.Get(vtk.vtkStreamingDemandDrivenPipeline.UPDATE_TIME_STEP())
sphere = vtk.vtkSphereSource()
sphere.SetCenter(t * 2, 0, 0)
sphere.Update()
output.ShallowCopy(sphere.GetOutput())
return 1
source = MovingSphereSource()
source.Update()
group = vtk.vtkMultiBlockFromTimeSeriesFilter()
group.SetInputConnection(source.GetOutputPort())
group.Update()
ren1 = vtk.vtkRenderer()
renWin = vtk.vtkRenderWindow()
renWin.AddRenderer(ren1)
iren = vtk.vtkRenderWindowInteractor()
iren.SetRenderWindow(renWin)
mapper = vtk.vtkCompositePolyDataMapper2()
mapper.SetInputConnection(group.GetOutputPort())
mapper.Update()
actor = vtk.vtkActor()
actor.SetMapper(mapper)
ren1.AddActor(actor)
ren1.ResetCamera()
renWin.SetSize(300, 300)
renWin.Render()
# render the image
#
iren.Initialize()
#iren.Start()
/*=========================================================================
Program: Visualization Toolkit
Module: vtkMultiBlockFromTimeSeriesFilter.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 "vtkMultiBlockFromTimeSeriesFilter.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkMultiBlockDataSet.h"
#include "vtkObjectFactory.h"
#include "vtkStreamingDemandDrivenPipeline.h"
vtkStandardNewMacro(vtkMultiBlockFromTimeSeriesFilter)
vtkMultiBlockFromTimeSeriesFilter::vtkMultiBlockFromTimeSeriesFilter()
{
this->UpdateTimeIndex = 0;
}
vtkMultiBlockFromTimeSeriesFilter::~vtkMultiBlockFromTimeSeriesFilter()
{
}
int vtkMultiBlockFromTimeSeriesFilter::FillInputPortInformation(
int, vtkInformation *info)
{
info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataObject");
return 1;
}
int vtkMultiBlockFromTimeSeriesFilter::RequestInformation(vtkInformation *vtkNotUsed(request),
vtkInformationVector** inInfo, vtkInformationVector* vtkNotUsed(outInfo))
{
this->UpdateTimeIndex = 0;
vtkInformation *info = inInfo[0]->GetInformationObject(0);
int len = info->Length(vtkStreamingDemandDrivenPipeline::TIME_STEPS());
double *timeSteps = info->Get(vtkStreamingDemandDrivenPipeline::TIME_STEPS());
this->TimeSteps.resize(len);
std::copy(timeSteps, timeSteps + len, this->TimeSteps.begin());
this->TempDataset = vtkSmartPointer<vtkMultiBlockDataSet>::New();
this->TempDataset->SetNumberOfBlocks(len);
return 1;
}
int vtkMultiBlockFromTimeSeriesFilter::RequestUpdateExtent(vtkInformation *vtkNotUsed(request),
vtkInformationVector** inInfo, vtkInformationVector* vtkNotUsed(outInfo))
{
if (this->TimeSteps.size() > static_cast<size_t>(this->UpdateTimeIndex))
{
vtkInformation *info = inInfo[0]->GetInformationObject(0);
info->Set(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP(),
this->TimeSteps[this->UpdateTimeIndex]);
}
return 1;
}
int vtkMultiBlockFromTimeSeriesFilter::RequestData(vtkInformation *request,
vtkInformationVector** inInfo, vtkInformationVector* outInfo)
{
vtkInformation *info = inInfo[0]->GetInformationObject(0);
vtkDataObject *data = vtkDataObject::GetData(info);
vtkSmartPointer<vtkDataObject> clone =
vtkSmartPointer<vtkDataObject>::Take(data->NewInstance());
clone->ShallowCopy(data);
this->TempDataset->SetBlock(this->UpdateTimeIndex, clone);
if (this->UpdateTimeIndex < static_cast<vtkTypeInt64>(this->TimeSteps.size()) - 1)
{
this->UpdateTimeIndex++;
request->Set(vtkStreamingDemandDrivenPipeline::CONTINUE_EXECUTING(), 1);
}
else
{
vtkMultiBlockDataSet *output = vtkMultiBlockDataSet::GetData(outInfo);
output->ShallowCopy(this->TempDataset);
for (unsigned i = 0; i < this->TempDataset->GetNumberOfBlocks(); ++i)
{
this->TempDataset->SetBlock(i, NULL);
}
request->Remove(vtkStreamingDemandDrivenPipeline::CONTINUE_EXECUTING());
}
return 1;
}
void vtkMultiBlockFromTimeSeriesFilter::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
}
/*=========================================================================
Program: Visualization Toolkit
Module: vtkMultiBlockFromTimeSeriesFilter.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 vtkMultiBlockFromTimeSeriesFilter - collects multiple inputs into one multi-group dataset
// .SECTION Description
// vtkMultiBlockFromTimeSeriesFilter is a 1 to 1 filter that merges multiple
// time steps from the input into one multiblock dataset. It will assign each
// time step from the input to one group of the multi-block dataset and will
// assign each timestep's data as a block in the multi-block datset.
#ifndef vtkMultiBlockFromTimeSeriesFilter_h
#define vtkMultiBlockFromTimeSeriesFilter_h
#include "vtkFiltersGeneralModule.h" // For export macro
#include "vtkMultiBlockDataSetAlgorithm.h"
#include "vtkSmartPointer.h" // Smart pointer
#include <vector> // Vector to hold timesteps
class vtkMultiBlockDataSet;
class VTKFILTERSGENERAL_EXPORT vtkMultiBlockFromTimeSeriesFilter : public vtkMultiBlockDataSetAlgorithm
{
public:
vtkTypeMacro(vtkMultiBlockFromTimeSeriesFilter,vtkMultiBlockDataSetAlgorithm);
void PrintSelf(ostream& os, vtkIndent indent);
static vtkMultiBlockFromTimeSeriesFilter *New();
protected:
vtkMultiBlockFromTimeSeriesFilter();
virtual ~vtkMultiBlockFromTimeSeriesFilter();
int FillInputPortInformation(int, vtkInformation *) VTK_OVERRIDE;
int RequestInformation(vtkInformation *,
vtkInformationVector **,
vtkInformationVector *) VTK_OVERRIDE;
int RequestUpdateExtent(vtkInformation *,
vtkInformationVector **,
vtkInformationVector *) VTK_OVERRIDE;
int RequestData(vtkInformation *,
vtkInformationVector **,
vtkInformationVector *) VTK_OVERRIDE;
private:
vtkMultiBlockFromTimeSeriesFilter(const vtkMultiBlockFromTimeSeriesFilter&) VTK_DELETE_FUNCTION;
void operator=(const vtkMultiBlockFromTimeSeriesFilter&) VTK_DELETE_FUNCTION;
int UpdateTimeIndex;
std::vector<double> TimeSteps;
vtkSmartPointer<vtkMultiBlockDataSet> TempDataset;
};
#endif
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