Hi all, I use paraview to view results of transient combined finite discrete element analyses. My *.pvd input files look like this:
<?xml version="1.0"?> <VTKFile type="Collection" version="0.1" byte_order="LittleEndian"> <Collection> <DataSet timestep="0,0000000000e+00" group="0" part="0" file="se_relax_p_000000000.vtp"/> <DataSet timestep="0,0000000000e+00" group="0" part="1" file="se_relax_c_000000000.vtp"/> <DataSet timestep="0,0000000000e+00" group="0" part="2" file="se_relax_r_000000000.vtp"/> <DataSet timestep="0,0000000000e+00" group="0" part="3" file="se_relax_f_000000000.vtu"/> <DataSet timestep="1,0000000000e-03" group="0" part="0" file="se_relax_p_000002000.vtp"/> <DataSet timestep="1,0000000000e-03" group="0" part="1" file="se_relax_c_000002000.vtp"/> . . . The results are stored for each part (DE, contact, rigid-boundaries, FE) and time step in a single *.vtp / *.vtu file. I can load the *.pvd file into paraview apply non-temporal filters to the single parts and animate without problems. Now I wrote a temporal filter to show the trajectories of data points. It requires a time series of datasets with constant number of points and generates the trajectories of the points as a single PolyData. The filter works fine if my input does not consist of multiple parts: <?xml version="1.0"?> <VTKFile type="Collection" version="0.1" byte_order="LittleEndian"> <Collection> <DataSet timestep="0,0000000000e+00" group="0" part="0" file="se_relax_p_000000000.vtp"/> <DataSet timestep="1,0000000000e-03" group="0" part="0" file="se_relax_p_000002000.vtp"/> <DataSet timestep="2,0000000000e-03" group="0" part="0" file="se_relax_p_000004000.vtp"/> . . . But if I load my complete data with multiple parts and apply the ExtractBlock filter first to show the trajectories for one part only, ParaView complains: ERROR: In /home/wellmann/program/ParaView-3.4.0/VTK/Filtering/vtkDemandDrivenPipeline.cxx, line 822 vtkCompositeDataPipeline (0x1e6bbc0): Input for connection index 0 on input port index 0 for algorithm vtkExtractBlock(0x11f4ddc0) is of type vtkTemporalDataSet, but a vtkMultiBlockDataSet is required. In writing the filter I tried to follow the "Time Dependent Processing in a Parallel Pipeline Architecture" paper and was geared to the vtkTemporalInterpolation filter. However I am not confident about how this temporal stuff works out: In my understanding my filter tells the upstream part of the pipeline to loop over the requested time steps and join the data into a TemporalDataSet that my vtkTrajectories filter can work with. Indeed if I connect any non temporal filter to the upstream pipeline part (instead of the vtkTrajectories filter) I can animate over the time-steps so the upstream pipeline is able to provide data sets for each time step. Furthermore, as said above, the filter works fine if there is no multiblock data as source of the pipeline. Thanks for any hints, Christian (the filter files are attached) -- Christian Wellmann Institute of Continuum Mechanics Leibniz Universitaet Hannover Appelstr. 11 30167 Hannover, Germany phone: +49 511 762 2285 fax: +49 511 762 5496 email: wellm...@ikm.uni-hannover.de
<ServerManagerConfiguration> <ProxyGroup name="filters"> <SourceProxy name="Trajectories" class="vtkTrajectories" label="Trajectories"> <Documentation long_help="Displays trajectories of data points." short_help="Displays trajectories of data points."> </Documentation> <InputProperty name="Input" command="SetInputConnection"> <ProxyGroupDomain name="groups"> <Group name="sources"/> <Group name="filters"/> </ProxyGroupDomain> <DataTypeDomain name="input_type" composite_data_supported="1"> <DataType value="vtkDataObject"/> </DataTypeDomain> </InputProperty> <DoubleVectorProperty label="time begin" name="SetTimeBegin" command="SetTimeBegin" number_of_elements="1" default_values="0.0" > </DoubleVectorProperty> <DoubleVectorProperty label="time end" name="SetTimeEnd" command="SetTimeEnd" number_of_elements="1" default_values="-1.0" > </DoubleVectorProperty> <IntVectorProperty label="max. number of points" name="SetMaxNumPoints" command="SetMaxNumPoints" number_of_elements="1" default_values="0" > </IntVectorProperty> <IntVectorProperty label="minimum point ID" name="SetMinID" command="SetMinID" number_of_elements="1" default_values="0" > </IntVectorProperty> <IntVectorProperty label="maximum point ID" name="SetMaxID" command="SetMaxID" number_of_elements="1" default_values="0" > </IntVectorProperty> <IntVectorProperty label="incrment point ID" name="SetIncrID" command="SetIncrID" number_of_elements="1" default_values="1" > </IntVectorProperty> </SourceProxy> </ProxyGroup> </ServerManagerConfiguration>
#include "vtkTrajectories.h" #include "vtkTemporalDataSet.h" #include "vtkObjectFactory.h" #include "vtkDataSet.h" #include "vtkDataSetAttributes.h" #include "vtkInformation.h" #include "vtkInformationVector.h" #include "vtkStreamingDemandDrivenPipeline.h" #include "vtkCompositeDataPipeline.h" #include "vtkPoints.h" #include "vtkPolyData.h" #include "vtkFloatArray.h" #include "vtkCellArray.h" #include "vtkPointData.h" #include "vtkMath.h" #include "vtkstd/vector" vtkCxxRevisionMacro(vtkTrajectories, "$Revision: 1.0 $"); vtkStandardNewMacro(vtkTrajectories); vtkTrajectories::vtkTrajectories() { this->SetNumberOfInputPorts(1); this->SetNumberOfOutputPorts(1); this->TimeBegin = 0.0; this->TimeEnd = -1.0; this->MaxNumPoints = 0; } vtkTrajectories::~vtkTrajectories(){} void vtkTrajectories::PrintSelf(ostream& os, vtkIndent indent) { this->Superclass::PrintSelf(os,indent); os << indent << "TimeBegin: " << this->TimeBegin << "\n" << indent << "TimeEnd: " << this->TimeEnd << "\n" << indent << "MaxNumPoints: " << this->MaxNumPoints << "\n"; } int vtkTrajectories::RequestInformation(vtkInformation* vtkNotUsed(request), vtkInformationVector** inputVector, vtkInformationVector* outputVector) { // get the info objects vtkInformation *outInfo = outputVector->GetInformationObject(0); vtkInformation *inInfo = inputVector[0]->GetInformationObject(0); // this filter provides no time varying output if(outInfo->Has(vtkStreamingDemandDrivenPipeline::TIME_RANGE())) outInfo->Remove(vtkStreamingDemandDrivenPipeline::TIME_RANGE()); if(outInfo->Has(vtkStreamingDemandDrivenPipeline::TIME_STEPS())) outInfo->Remove(vtkStreamingDemandDrivenPipeline::TIME_STEPS()); // check time input double* InTimeRange(0); double* InTimeSteps(0); int NumInTimeSteps(0); if(inInfo->Has(vtkStreamingDemandDrivenPipeline::TIME_RANGE())){ InTimeRange = inInfo->Get(vtkStreamingDemandDrivenPipeline::TIME_RANGE()); }else{ vtkErrorMacro(<<"No time range to build trajectories."); return 0; } if(inInfo->Has(vtkStreamingDemandDrivenPipeline::TIME_STEPS())){ InTimeSteps = inInfo->Get(vtkStreamingDemandDrivenPipeline::TIME_STEPS()); NumInTimeSteps = inInfo->Length(vtkStreamingDemandDrivenPipeline::TIME_STEPS()); if(NumInTimeSteps < 2){ vtkErrorMacro(<<"Less than two time steps to build trajectories."); return 0; } }else{ vtkErrorMacro(<<"No discrete time steps to build trajectories."); return 0; } // if upper time bound smaller than lower use complete time range if(this->TimeEnd < this->TimeBegin){ this->TimeBegin = InTimeRange[0]; this->TimeEnd = InTimeRange[1]; }else{ // check available time bounds against trajectories time bounds if(InTimeRange[0] > this->TimeBegin){ vtkWarningMacro(<<"Lower time bound for trajectory: " << this->TimeBegin <<" smaller than available lower time bound " << InTimeRange[0] <<". Resetting trajectory time bound."); this->TimeBegin = InTimeRange[0]; } if(InTimeRange[1] < this->TimeEnd){ vtkWarningMacro(<<"Upper time bound for trajectory: " << this->TimeEnd <<" greater than available upper time bound " << InTimeRange[1] <<". Resetting trajectory time bound."); this->TimeEnd = InTimeRange[1]; } } return 1; } int vtkTrajectories::RequestUpdateExtent(vtkInformation * vtkNotUsed(request), vtkInformationVector **inputVector, vtkInformationVector * vtkNotUsed(outputVector)) { // get the info objects vtkInformation *inInfo = inputVector[0]->GetInformationObject(0); // determine the input time steps needed double* InTimeRange(0); double* InTimeSteps(0); int NumInTimeSteps(0); if(inInfo->Has(vtkStreamingDemandDrivenPipeline::TIME_RANGE())){ InTimeRange = inInfo->Get(vtkStreamingDemandDrivenPipeline::TIME_RANGE()); }else{ vtkErrorMacro(<<"No time range to build trajectories."); return 0; } if(inInfo->Has(vtkStreamingDemandDrivenPipeline::TIME_STEPS())){ InTimeSteps = inInfo->Get(vtkStreamingDemandDrivenPipeline::TIME_STEPS()); NumInTimeSteps = inInfo->Length(vtkStreamingDemandDrivenPipeline::TIME_STEPS()); if(NumInTimeSteps < 2){ vtkErrorMacro(<<"Less than two time steps to build trajectories."); return 0; } }else{ vtkErrorMacro(<<"No discrete time steps to build trajectories."); return 0; } vtkstd::vector<bool> InTimeStepsUsed(NumInTimeSteps,false); int NumInTimeStepsUsed(0); for(int i=0;i<NumInTimeSteps;i++){ if((! (InTimeSteps[i] < this->TimeBegin)) && (! (InTimeSteps[i] > this->TimeEnd))){ InTimeStepsUsed[i] = true; NumInTimeStepsUsed++; } } double* InUpdateTimes = new double[NumInTimeStepsUsed]; int k(0); for(int i=0;i<NumInTimeSteps;i++){ if(InTimeStepsUsed[i]){ InUpdateTimes[k++]=InTimeSteps[i]; } } inInfo->Set(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEPS(), InUpdateTimes,NumInTimeStepsUsed); inInfo->Set(vtkCompositeDataPipeline::REQUIRES_TIME_DOWNSTREAM(),1); delete[] InUpdateTimes; return 1; } int vtkTrajectories::RequestData(vtkInformation* vtkNotUsed(request), vtkInformationVector** inputVector, vtkInformationVector* outputVector) { vtkInformation *inInfo = inputVector[0]->GetInformationObject(0); vtkInformation *outInfo = outputVector->GetInformationObject(0); vtkTemporalDataSet *inData = vtkTemporalDataSet::SafeDownCast(inInfo->Get(vtkDataObject::DATA_OBJECT())); vtkPolyData *outData = vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())); if(! (inData && outData)) return 1; int NumTimeSteps = inData->GetNumberOfTimeSteps(); // get start points of trajectories vtkDataSet *inData_i = vtkDataSet::SafeDownCast(inData->GetTimeStep(0)); if(! inData_i){ vtkErrorMacro(<< "vtkTrajectoriesFilter not applicable to composite data sets."); return 0; } // number of points in each time-step data set const int NumPoints_i = inData_i->GetNumberOfPoints(); // which points will be used int NumPoints_i_used(NumPoints_i); bool use_point[NumPoints_i]; for(int i=0;i<NumPoints_i;i++) use_point[i]=true; if(this->MaxID > this->MinID){ if(this->MinID < 0) this->MinID = 0; else if(this->MinID >= NumPoints_i) this->MinID = NumPoints_i-1; if(this->MaxID < this->MinID) this->MaxID = this->MinID; else if(this->MaxID >= NumPoints_i) this->MaxID = NumPoints_i-1; NumPoints_i_used = (this->MaxID - this->MinID) / this->IncrID; if(NumPoints_i_used < 1) NumPoints_i_used = 1; for(int i=0;i<NumPoints_i;i++) use_point[i]=false; for(int i=0;i<NumPoints_i_used;i++) use_point[this->MinID + i * this->IncrID] = true; }else if(MaxNumPoints > 0 && MaxNumPoints < NumPoints_i){ NumPoints_i_used = MaxNumPoints; if(2 * NumPoints_i_used < MaxNumPoints){ for(int i=0;i<NumPoints_i;i++) use_point[i]=false; int k(0),ID; while(k < MaxNumPoints){ ID = static_cast<int>(NumPoints_i * vtkMath::Random()); if(! use_point[ID]){ use_point[ID]=true; k++; } } }else{ int k(0),ID,q(NumPoints_i-MaxNumPoints); while(k < q){ ID = static_cast<int>(NumPoints_i * vtkMath::Random()); if(use_point[ID]){ use_point[ID]=false; k++; } } } } int pointIDs[NumPoints_i_used],k(0); for(int j=0;j<NumPoints_i;++j){ if(use_point[j]){ pointIDs[k++]=j; } } const int NumPoints = NumPoints_i_used * NumTimeSteps; vtkPoints* newPoints = vtkPoints::New(); newPoints->SetNumberOfPoints(NumPoints); for(int i=0;i<NumTimeSteps;++i){ inData_i = vtkDataSet::SafeDownCast(inData->GetTimeStep(i)); if(! inData_i){ vtkErrorMacro(<< "vtkTrajectoriesFilter not applicable to composite data sets."); return 0; } if(inData_i->GetNumberOfPoints() != NumPoints_i){ vtkErrorMacro(<< "Change of number of points between time steps."); return 0; } int base(i*NumPoints_i_used); for(int j=0;j<NumPoints_i_used;++j){ newPoints->SetPoint(base+j,inData_i->GetPoint(pointIDs[j])); } } outData->SetPoints(newPoints); newPoints->Delete(); vtkCellArray* lines = vtkCellArray::New(); lines->Allocate(NumPoints_i_used * (NumTimeSteps + 1)); vtkIdType* LineIDs = new vtkIdType[NumTimeSteps]; for(int i=0;i<NumPoints_i_used;i++){ LineIDs[0]=i; for(int j=1;j<NumTimeSteps;j++) LineIDs[j]=LineIDs[j-1]+NumPoints_i_used; lines->InsertNextCell(NumTimeSteps,LineIDs); } delete[] LineIDs; outData->SetLines(lines); lines->Delete(); double *InTimes = inData->GetInformation()->Get(vtkDataObject::DATA_TIME_STEPS()); int NumInTimes = inData->GetInformation()->Length(vtkDataObject::DATA_TIME_STEPS()); if(NumInTimes != NumTimeSteps){ vtkErrorMacro(<< "NumInTimes != NumTimeSteps!"); return 0; } float vals[NumTimeSteps]; for(int j=0;j<NumTimeSteps;j++) vals[j] = static_cast<float>(InTimes[j]); vtkFloatArray* newScalars = vtkFloatArray::New(); newScalars->SetNumberOfValues(NumPoints); for(int j=0;j<NumTimeSteps;j++) for(int i=j*NumPoints_i_used;i<(j+1)*NumPoints_i_used;i++) newScalars->SetValue(i,vals[j]); newScalars->SetName("time"); outData->GetPointData()->SetScalars(newScalars); newScalars->Delete(); return 1; } int vtkTrajectories::FillOutputPortInformation(int vtkNotUsed(port), vtkInformation* info) { info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkPolyData"); return 1; } int vtkTrajectories::FillInputPortInformation(int vtkNotUsed(port), vtkInformation* info) { info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkTemporalDataSet"); return 1; }
/* vtkTrajectories * * vtkTrajectories converts point positions given in a time series of * data files into the point trajectories returned as PolyData. */ #ifndef __vtkTrajectories_h #define __vtkTrajectories_h #include "vtkPolyDataAlgorithm.h" class VTK_GRAPHICS_EXPORT vtkTrajectories : public vtkPolyDataAlgorithm { public: static vtkTrajectories *New(); vtkTypeRevisionMacro(vtkTrajectories,vtkPolyDataAlgorithm); void PrintSelf(ostream& os, vtkIndent indent); vtkSetMacro(TimeBegin, double); vtkGetMacro(TimeBegin, double); vtkSetMacro(TimeEnd, double); vtkGetMacro(TimeEnd, double); vtkSetMacro(MaxNumPoints, int); vtkGetMacro(MaxNumPoints, int); vtkSetMacro(MinID, int); vtkGetMacro(MinID, int); vtkSetMacro(MaxID, int); vtkGetMacro(MaxID, int); vtkSetMacro(IncrID, int); vtkGetMacro(IncrID, int); protected: vtkTrajectories(); ~vtkTrajectories(); int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *); int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *); int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *); int FillOutputPortInformation(int port, vtkInformation *info); int FillInputPortInformation(int port, vtkInformation *info); // bounds for time interval to plot trajectories double TimeBegin; double TimeEnd; // maximum number of points to use (if < 1 not considered) int MaxNumPoints; // minimum, maximum and increment of point IDs to use (if MaxID <= MinID not considered) int MinID; int MaxID; int IncrID; private: vtkTrajectories(const vtkTrajectories&); // Not implemented. void operator=(const vtkTrajectories&); // Not implemented. }; #endif
_______________________________________________ Powered by www.kitware.com Visit other Kitware open-source projects at http://www.kitware.com/opensource/opensource.html Please keep messages on-topic and check the ParaView Wiki at: http://paraview.org/Wiki/ParaView Follow this link to subscribe/unsubscribe: http://www.paraview.org/mailman/listinfo/paraview