[Paraview] Reader Plugin should inherit from which class ?
Christian Richter
christian.richter at ovgu.de
Wed Nov 9 08:36:28 EST 2011
Dear List,
I've written a Plugin to read custom text files.
I'm using vtkPolyDataAlgorithm as Master-Class.
The code compiles fine and my Plugin is working in Paraview 3.12.
But if the plugin is active, it's not possible to open legecy VTK
files, because everytime my own reader class is used.
I think I should use another Superclass, but don't know which one I
should shoose.
I attached you the source code.
Thanks for help,
C. Richter
-------------- next part --------------
#include "vtkCSVImageReader.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkObjectFactory.h"
#include "vtkDataArray.h"
#include "vtkImageData.h"
#include "vtkPointData.h"
#include "vtkStreamingDemandDrivenPipeline.h"
#include "vtkTable.h"
#include "vtkVariant.h"
#include "vtkObject.h"
#include "vtkSmartPointer.h"
#include "vtkStreamingDemandDrivenPipeline.h"
#include "vtkInformation.h"
#include "vtkInformationVector.h"
#include "vtkObjectFactory.h"
#include "vtkByteSwap.h"
#include "vtkPointData.h"
#include "vtkPoints.h"
#include "vtkPolyData.h"
#include "vtkCellArray.h"
#include "vtkDataArray.h"
#include "vtkDoubleArray.h"
#include "vtkIntArray.h"
#include "vtkFloatArray.h"
#include <vtkstd/algorithm>
#include <vtkstd/vector>
#include <vtkstd/string>
#include <vtksys/ios/sstream>
vtkCxxRevisionMacro(vtkCSVImageReader, "$Revision: 0.5 $");
vtkStandardNewMacro(vtkCSVImageReader);
vtkCSVImageReader::vtkCSVImageReader()
{
this->FileName = 0;
this->File=0;
this->SetNumberOfInputPorts(0);
this->SetNumberOfOutputPorts(1);
}
vtkCSVImageReader::~vtkCSVImageReader()
{
if (this->File)
{
this->File->close();
delete this->File;
this->File = NULL;
}
this->SetFileName(0);
}
void vtkCSVImageReader::OpenFile()
{
if (!this->FileName)
{
vtkErrorMacro(<<"FileName must be specified.");
return;
}
// If the file was open close it.
if (this->File)
{
this->File->close();
delete this->File;
this->File = NULL;
}
// Open the new file.
#ifdef _WIN32
this->File = new ifstream(this->FileName, ios::in | ios::binary);
#else
this->File = new ifstream(this->FileName, ios::in);
#endif
if (! this->File || this->File->fail())
{
vtkErrorMacro(<< "Initialize: Could not open file " << this->FileName);
return;
}
}
int vtkCSVImageReader::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
{
this->OpenFile();
char line[256]; //increase to get more chars !
this->File->getline(line,sizeof(line),'\r\n');
vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
points->SetDataTypeToFloat();
points->Reset();
vtkSmartPointer<vtkFloatArray> PointID = vtkSmartPointer<vtkFloatArray>::New();
PointID->SetNumberOfComponents(1);
PointID->Reset();
PointID->SetName("Point ID");
PointID->SetComponentName(0,"Point ID");
vtkSmartPointer<vtkFloatArray> scalars = vtkSmartPointer<vtkFloatArray>::New();
scalars->SetNumberOfComponents(1);
scalars->Reset();
scalars->SetName("Scalars");
scalars->SetComponentName(0,"Scalar");
vtkSmartPointer<vtkFloatArray> colors = vtkSmartPointer<vtkFloatArray>::New();
colors->SetNumberOfComponents(4);
colors->SetComponentName(0,"Color 1");
colors->SetComponentName(1,"Color 2");
colors->SetComponentName(2,"Color 3");
colors->SetComponentName(3,"Color 4");
colors->SetName("color");
vtkSmartPointer<vtkFloatArray> v = vtkSmartPointer<vtkFloatArray>::New();
v->SetNumberOfComponents(3);
v->SetComponentName(0,"X");
v->SetComponentName(0,"Y");
v->SetComponentName(0,"Z");
v->SetName("v");
vtkSmartPointer<vtkFloatArray> w = vtkSmartPointer<vtkFloatArray>::New();
w->SetNumberOfComponents(3);
w->SetComponentName(0,"X");
w->SetComponentName(0,"Y");
w->SetComponentName(0,"Z");
w->SetName("omega");
vtkSmartPointer<vtkFloatArray> radius = vtkSmartPointer<vtkFloatArray>::New();
radius->SetNumberOfComponents(1);
radius->SetName("radius");
float val[3];
val[0]=val[1]=val[2]=0;
float scal[1];
scal[0]=0;
float vel[3];
vel[0]=vel[1]=vel[2]=0;
float omega[3];
omega[0]=omega[1]=omega[2]=0;
float c[4];
c[0]=c[1]=c[2]=c[3]=0;
int lc=0;
while ( this->File->getline(line,sizeof(line),'\r\n') ) {
if (int(line[1])!=0) {
int l = sizeof(line);
int numberOfNumbers = 0;
//Alle Leerzeichen durch '\0' = StringEnde ersetzen
for (int i = 0; i < l; i++)
{
if (line[i] == ' ')
{
line[i] = '\0';
numberOfNumbers++;
}
}
//Pointer auf den ersten Zahlenstring setzen
char* startp = line;
double id = atof(startp);
startp += strlen(startp) + 1;
PointID->InsertNextTuple1(id);
//position
for (int i = 0; i < 3; i++)
{
//Umwandlung nach double
double d = atof(startp);
val[i]=d;
//Auf den Beginn des n?sten Zahlenstrings setzen
startp += strlen(startp) + 1;
}
points->InsertNextPoint(val[0], val[1], val[2]);
//radius
double rad = atof(startp);
startp += strlen(startp) + 1;
radius->InsertNextTuple1(rad);
//velocity
for (int i = 0; i < 3; i++)
{
//Umwandlung nach double
double d = atof(startp);
vel[i]=d;
//Auf den Beginn des n?sten Zahlenstrings setzen
startp += strlen(startp) + 1;
}
v->InsertNextTuple(vel);
//omega
for (int i = 0; i < 3; i++)
{
//Umwandlung nach double
double d = atof(startp);
omega[i]=d;
//Auf den Beginn des n?sten Zahlenstrings setzen
startp += strlen(startp) + 1;
}
w->InsertNextTuple(omega);
//scalar
double sc = atof(startp);
startp += strlen(startp) + 1;
scalars->InsertNextTuple1(sc);
//colors
for (int i = 1; i <= 4; i++)
{
//Umwandlung nach double
double d = atof(startp);
c[i]=d;
//Auf den Beginn des n?sten Zahlenstrings setzen
startp += strlen(startp) + 1;
}
colors->InsertNextTuple(c);
lc++;
} //endif
}
vtkSmartPointer<vtkCellArray> vertices = vtkSmartPointer<vtkCellArray>::New();
vertices->Reset();
this->NumberOfPoints = points->GetNumberOfPoints();
for( vtkIdType j = 0; j < (vtkIdType)this->NumberOfPoints; ++j )
{
vertices->InsertNextCell( 1 );
vertices->InsertCellPoint( j );
}
// get the info object
vtkInformation *outInfo = outputVector->GetInformationObject(0);
// get the ouptut
vtkPolyData *myoutput = vtkPolyData::SafeDownCast(
outInfo->Get(vtkDataObject::DATA_OBJECT()));
myoutput->SetPoints(points);
myoutput->SetVerts(vertices);
myoutput->GetPointData()->AddArray(PointID);
myoutput->GetPointData()->AddArray(radius);
myoutput->GetPointData()->AddArray(v);
myoutput->GetPointData()->AddArray(w);
myoutput->GetPointData()->AddArray(scalars);
myoutput->GetPointData()->AddArray(colors);
myoutput->Modified();
return 1;
}
int vtkCSVImageReader::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
{
vtkInformation *outInfo = outputVector->GetInformationObject(0);
outInfo->Set(vtkStreamingDemandDrivenPipeline::MAXIMUM_NUMBER_OF_PIECES(),
-1);
return 1;
}
void vtkCSVImageReader::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent << "FileName: "
<< (this->FileName ? this->FileName : "(NULL)") << endl;
}
int vtkCSVImageReader::CanReadFile(const char *fname)
{
return 1;
}
-------------- next part --------------
#ifndef _vtkCSVImageReader_h
#define _vtkCSVImageReader_h
//#include "vtkImageAlgorithm.h"
#include <vtkPolyDataAlgorithm.h>
//#include <vtkImageAlgorithm.h>
class VTK_EXPORT vtkCSVImageReader : public vtkPolyDataAlgorithm
{
public:
static vtkCSVImageReader *New();
vtkTypeRevisionMacro(vtkCSVImageReader,vtkPolyDataAlgorithm);
void PrintSelf(ostream& os, vtkIndent indent);
// Specify file name.
vtkSetStringMacro(FileName);
vtkGetStringMacro(FileName);
int CanReadFile(const char* fname);
protected:
vtkCSVImageReader();
~vtkCSVImageReader();
int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *);
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *);
int Canreadfile(const char *fname);
char *FileName;
char *FieldDelimiterCharacters;
size_t NumberOfPoints;
void OpenFile();
ifstream *File;
private:
vtkCSVImageReader(const vtkCSVImageReader&);
void operator = (const vtkCSVImageReader&);
};
#endif
-------------- next part --------------
A non-text attachment was scrubbed...
Name: CSVImage.xml
Type: text/xml
Size: 2062 bytes
Desc: not available
URL: <http://www.paraview.org/pipermail/paraview/attachments/20111109/3f34d172/attachment-0002.bin>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: CSVImageGUI.xml
Type: text/xml
Size: 157 bytes
Desc: not available
URL: <http://www.paraview.org/pipermail/paraview/attachments/20111109/3f34d172/attachment-0003.bin>
More information about the ParaView
mailing list