<div dir="ltr"><div>Do you have a node partitioning of your grid? If you do you really need a cell partitioning of the grid which will have multiple points a the partition overlap. <br><br></div>Can you share your data?<br></div><div class="gmail_extra"><br><div class="gmail_quote">On Fri, Jun 26, 2015 at 5:58 AM,  <span dir="ltr"><<a href="mailto:u.utku.turuncoglu@be.itu.edu.tr" target="_blank">u.utku.turuncoglu@be.itu.edu.tr</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Hi,<br>
<br>
I am trying to create structured grid in an MPI application and i am using<br>
following C++ code to implement it<br>
<br>
extern "C" void createstgrid_(double* y, double* x,<br>
                              int* nxstart, int* nxend, int* nystart, int*<br>
nyend,<br>
                              int* nx, int* ny, int* nz, int* mpiSize,<br>
int* mpiRank) {<br>
  //<br>
  // Create structured grid<br>
  //<br>
  vtkSmartPointer<vtkStructuredGrid> sg =<br>
vtkSmartPointer<vtkStructuredGrid>::New();<br>
  sg->SetExtent (*nxstart-1, *nxend-1, *nystart-1, *nyend-1, 0, 0);<br>
<br>
  //<br>
  // Create data structure to store grid coordinates<br>
  //<br>
  int pnx = ((*nxend)-(*nxstart))+1;<br>
  int pny = ((*nyend)-(*nystart))+1;<br>
  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();<br>
<br>
  //<br>
  // Insert grid coordinate data as points<br>
  //<br>
  k = 0;<br>
  for (int j = 0; j < pny; j++) {<br>
    for (int i = 0; i < pnx; i++) {<br>
      points->InsertPoint(k, x[i+j*pnx], y[i+j*pnx], 0);<br>
      k = k+1;<br>
    }<br>
  }<br>
  sg->SetPoints(points);<br>
<br>
  //<br>
  // Set grid<br>
  //<br>
  vtkCPPythonAdaptorAPI::GetCoProcessorData()->GetInputDescriptionByName("input")->SetGrid(sg);<br>
  vtkCPPythonAdaptorAPI::GetCoProcessorData()->GetInputDescriptionByName("input")->SetWholeExtent(0,<br>
*nx-1, 0, *ny-1, 0, 0);<br>
}<br>
<br>
so, when i run the code with grid writer, it creates the individual file<br>
(*.vts) for each processor and one extra for the combination (*.pvts).<br>
Anyway, the individual files looks fine in Paraview but if i open the file<br>
with pvts extension, i am getting following error,<br>
<br>
ERROR: In<br>
/Users/turuncu/Qsync/progs/paraview-4.3.1/src/VTK/IO/XML/vtkXMLPStructuredDataReader.cxx,<br>
line 461<br>
vtkXMLPStructuredGridReader (0x7f961c0910a0): No available piece provides<br>
data for the following extents:<br>
0 140  0 1  0 0<br>
0 1  1 58  0 0<br>
0 140  58 59  0 0<br>
0 1  59 114  0 0<br>
71 72  1 58  0 0<br>
71 72  59 114  0 0<br>
<br>
The UpdateExtent cannot be filled.<br>
<br>
What is the meaning of this message? It means that some of the grid data<br>
(or points) are not filled by the code. I also check the loop in the C++<br>
code and print out the number of points for each processor<br>
<br>
0 number of points = 4118<br>
1 number of points = 4047<br>
2 number of points = 4060<br>
3 number of points = 3990<br>
<br>
and the total is 16215, which is consistent with the size of the grid but<br>
if i check the individual vts files from the Paraview, the files shows<br>
different number of points for each file<br>
<br>
0 number of points = 4118<br>
1 number of points = 3976<br>
2 number of points = 4002<br>
3 number of points = 3864<br>
<br>
So, only first file matchs with the C++ code but others not. I am not sure<br>
what is happening in there and if you have any suggestion that will be<br>
great.<br>
<br>
Anyway, in this case i did not use vtkMultiBlockDataSet but after having<br>
correct representation of the grid, i am planing use SetBlock method of<br>
vtkMultiBlockDataSet to assign decomposition elements as individual<br>
blocks.<br>
<br>
Regards,<br>
<br>
--ufuk<br>
<br>
_______________________________________________<br>
Powered by <a href="http://www.kitware.com" rel="noreferrer" target="_blank">www.kitware.com</a><br>
<br>
Visit other Kitware open-source projects at <a href="http://www.kitware.com/opensource/opensource.html" rel="noreferrer" target="_blank">http://www.kitware.com/opensource/opensource.html</a><br>
<br>
Please keep messages on-topic and check the ParaView Wiki at: <a href="http://paraview.org/Wiki/ParaView" rel="noreferrer" target="_blank">http://paraview.org/Wiki/ParaView</a><br>
<br>
Search the list archives at: <a href="http://markmail.org/search/?q=ParaView" rel="noreferrer" target="_blank">http://markmail.org/search/?q=ParaView</a><br>
<br>
Follow this link to subscribe/unsubscribe:<br>
<a href="http://public.kitware.com/mailman/listinfo/paraview" rel="noreferrer" target="_blank">http://public.kitware.com/mailman/listinfo/paraview</a><br>
</blockquote></div><br></div>