[Paraview] Writing a Paraview Reader - Only part of the data available

Dan Lipsa dan.lipsa at kitware.com
Mon Aug 24 10:31:05 EDT 2015


James,
You are missing image->SetDimensions() - that may be a reason why your
scalar does not have all values.

vtkRectilinear grid is like an image data with variable extents and warped
grid.


On Fri, Aug 21, 2015 at 6:55 AM, James Furness <pcxjf1 at nottingham.ac.uk>
wrote:

> Hello,
>
> I have a program that will produce an output either in csv, or a binary
> format (which I have written and can alter).
>
> The data points are constructed by taking 3 vectors, and stepping at set
> intervals along these vectors to construct a 3D set of data. Ultimately
> there shouldn’t be a restriction on the three vectors, so the data may not
> be in a rectilinear array. But I’m leaving this for now and assuming the
> data is constructed from 3 orthogonal Axis aligned vectors to create a
> rectilinear structured grid. When I have a better grasp of preview readers
> I’ll update the reader to account for this. I guess it would require an
> unstructured data type to handle non-orthogonal grids?
>
> Also, the user is free to select what data they save from the program.
> This could be many fields, both vector and scalar quantities. What has been
> selected is stored in a header to the file, and known on reading. For my
> toy example I save only one scalar quantity.
>
> I’ve been trying to write a reader for paraview to import the binary files
> output. The csv format works ok, but the files can become large and it is a
> pain to have to construct the vector quantities from 3 scalars in each
> case, I want to avoid the user having to do this work. Another reason for
> writing a custom reader is to include additional information into the
> binary file for visualisation (e.g. the position of atoms, not a field
> quantity).
>
> With that background and reasoning I have managed to follow the
> documentation enough to read a basic test file, get the origin and point
> spacing, and read a single scalar quantity into a vtkDataArray owned by a
> vtkImageData object. Paraview picks this up seemingly correctly, with one
> flaw:
>
> It only takes 2 elements from each dimension, 0 and 1 displaying 8
> elements in total. I am confident it has read the other values correctly as
> the Information tab reports “X extent 0 to 3” “X range: -1 to 1” (as
> expected) and similar for the other dimensions. If I view the data in
> spreadsheet layout the 8 values it has seem correct, but the others are
> missing.
>
> The code I am using in the reader’s RequestData(  ) function is added
> below.
>
> Does anyone know why this may be happening and how to fix it? Also, any
> advice on how I should structure the reader to handle the non-orthogonal
> data?
>
> Thanks for your time, and thank you for such a stellar program. The
> results paraview produces for this data is brilliant, hence my want to make
> it convenient for our users!
>
> Regards,
> James Furness
>
>
> CODE for RequestData( ) Method:
> ———————————————————————
>
> int LondonReader::RequestData(
>         vtkInformation*,
>         vtkInformationVector**,
>         vtkInformationVector *outputVector)
> {
>     vtkWarningMacro("Requesting the data!");
>
>
>     ifstream finp;
>     finp.open(this->FileName, ios::in | ios::binary);
>
>     if (finp.is_open()) {
>         cerr << "File is open without problem!" << endl;
>     }   else    {
>         cerr << "File failed to open :(" << endl;
>         return 0;
>     }
>
>     // size of real numbers may not be 8. Check for this from file header
>     int realSize;
>     finp.read((char*)&realSize, sizeof(int));
>     if(realSize != 8)   {
>         cerr << "Not implimented yet!" << endl;
>         return 0;
>     }
>
>     // number of data fields
>     int nFields;
>     finp.read((char*)&nFields, sizeof(int));
>
>     vtkImageData* image = vtkImageData::GetData(outputVector);
>
>     // Read the dimensions of the grid and set the extent accordingly
>     int gridDim[3];
>     finp.read((char*)&gridDim, 3*sizeof(int));
>     int extent[6] = {0, gridDim[0]-1, 0, gridDim[1]-1, 0, gridDim[2]-1};
>     image->SetExtent(extent);
>
>
>     // Read the field names from the file
>     std::vector<std::string> fields;
>     std::string strBuf;
>     for (int i = 0; i < nFields; i++)   {
>         std::getline( finp, strBuf, '\0');
>         fields.push_back(strBuf);
>         cerr << "Printing Fields (" << i << "): " << fields[i] << endl;
>     }
>
>     // setup image for only one field for test case
>     image->AllocateScalars(VTK_FLOAT, 1);
>     vtkDataArray* scalars = image->GetPointData()->GetScalars();
>
>     // currently there is only one field 'rho'
>     scalars->SetName(fields[3].c_str());
>
>     double x, y, z, rho;
>     double oX, oY, oZ;      //origin coordinates
>     double sX, sY, sZ;      //spacing of points
>
>     for (vtkIdType itx = 0; itx < gridDim[0]; itx++)    {
>         for (vtkIdType ity = 0; ity < gridDim[1]; ity++)    {
>             for (vtkIdType itz = 0; itz < gridDim[2]; itz++)    {
>                 finp.read((char*)&x, realSize);
>                 finp.read((char*)&y, realSize);
>                 finp.read((char*)&z, realSize);
>                 finp.read((char*)&rho, realSize);
>
>                 // Find and set the origin and spacing
>                 if (itx == 0 && ity == 0 && itz == 0)   {
>                     image->SetOrigin(x, y, z);
>                     oX = x; oY = y; oZ = z;
>                 }   else if (itx == 1 && ity == 0 && itz == 0)  {
>                     sX = x - oX;
>                 }   else if (itx == 0 && ity == 1 && itz == 0)  {
>                     sY = y - oY;
>                 }   else if (itx == 0 && ity == 0 && itz == 1)  {
>                     sZ = z - oZ;
>                 }
>                 //check correct read.
>                 cerr << x << "," << y << "," << z << "," << rho << ", at "
> << itx*(gridDim[1]*gridDim[2]) + ity*gridDim[2] + itz << endl;
>
>                 //add value
>                 scalars->SetTuple1(itx*gridDim[1]*gridDim[2] +
> ity*gridDim[2] + itz,
>                         rho);
>             }
>         }
>     }
>
>     image->SetSpacing(sX, sY, sZ);
>
>     image->Print(cerr);
>
>     return 1;
> }
>
>
>
>
>
> This message and any attachment are intended solely for the addressee
> and may contain confidential information. If you have received this
> message in error, please send it back to me, and immediately delete it.
>
> Please do not use, copy or disclose the information contained in this
> message or in any attachment.  Any views or opinions expressed by the
> author of this email do not necessarily reflect the views of the
> University of Nottingham.
>
> This message has been checked for viruses but the contents of an
> attachment may still contain software viruses which could damage your
> computer system, you are advised to perform your own checks. Email
> communications with the University of Nottingham may be monitored as
> permitted by UK legislation.
>
> _______________________________________________
> 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
>
> Search the list archives at: http://markmail.org/search/?q=ParaView
>
> Follow this link to subscribe/unsubscribe:
> http://public.kitware.com/mailman/listinfo/paraview
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://public.kitware.com/pipermail/paraview/attachments/20150824/99916dd7/attachment.html>


More information about the ParaView mailing list