<html>
  <head>
    <meta content="text/html; charset=utf-8" http-equiv="Content-Type">
  </head>
  <body bgcolor="#FFFFFF" text="#000000">
    <div class="moz-cite-prefix">Hi,<br>
      <br>
      Actually, i am trying to create simplified case to test some ideas
      to improve the performance of my Catalyst application. In this
      case, i would like to create C++ code one for point data in
      vtkMultiBlockDataSet and other is cell data again in
      vtkMultiBlockDataSet. If i manage to create these two cases and
      ghost regions to have correct representation of the entire field,
      that would be great. As i know, for point data, i need to create
      ghost region (or overlapped regions) for point data but not for
      cell ones. Is it correct? Is there any other way to escape
      generation of ghost regions (both for point and cell data)? It is
      not easy to find any example for the cases like these. So, i am
      looking for an answer in the user lists. <br>
      <br>
      By the way, i could not see any "off-by-1" bug in the point
      coordinates. Is it possible to give more information about it.
      Thanks.<br>
      <br>
      Regards,<br>
      <br>
      --ufuk  <br>
      <br>
      On 02/06/2017 20:16, Andy Bauer wrote:<br>
    </div>
    <blockquote
cite="mid:CAMaOp+GSHL7Pmvcc2GXj_fcN3=Ke5fsU+XdtK7GW197E4H+m8Q@mail.gmail.com"
      type="cite">
      <div dir="ltr">Hi Ufuk,
        <div><br>
        </div>
        <div>If there is one vtkStructuredGrid per multi-block you don't
          need the whole extent. Using a vtkStructuredGrid in a
          multiblock though will result in each grid block to be
          considered a separate grid so if you do something like
          extracting the surface it probably won't look correct unless
          you generate your own ghost cells (the ghost cell generator
          filter in PV won't work). </div>
        <div><br>
        </div>
        <div>I only briefly looked through your code snippet but I
          suspect there's an "off-by-1" bug in the way you're setting
          your point coordinates for the structured grids but once
          that's fixed it will work the way you want it to. </div>
        <div><br>
        </div>
        <div>I probably don't know your full use case but I would
          suggest for non-Catalyst stuff to just avoid the multiblock
          here and just generate structured grid pieces.</div>
        <div><br>
        </div>
        <div>Best,</div>
        <div>Andy</div>
        <div class="gmail_extra"><br>
          <div class="gmail_quote">On Fri, Jun 2, 2017 at 8:13 AM, Ufuk
            Turuncoglu <span dir="ltr"><<a moz-do-not-send="true"
                href="mailto:ufuk.turuncoglu@itu.edu.tr" target="_blank">ufuk.turuncoglu@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
              Again,<br>
              <br>
              I create following test code to write dummy data
              (checkerboard with 0 and 1) using vtkStructuredGrid in
              multiblock format. Without extent in left and top of the
              tiles the data looks correct but there is a gap between
              pieces (see attached screenshot). To fix the issue, i add
              a extra piece of code to overlap one row and columns
              (using hasLeft and hasTop control). The problem is that
              when i run the code with this modification, i am getting
              totally wrong grid representation (screen2). I just wonder
              that do i need to set whole extent in somewhere. I tried
              to add following to set whole extent to multiblock
              dataset,<br>
              <br>
              int extent[6] = {0, N-1, 0, N-1, 0, 0};<br>
              mb->GetInformationObject(0)->S<wbr>et(vtkStreamingDemandDrivenPip<wbr>eline::WHOLE_EXTENT(),
              extent, 6);<br>
              <br>
              but it won't work as i expected. Any suggestion? Thanks.<br>
              <br>
              Regards,<br>
              <br>
              --ufuk<br>
              <br>
              <br>
              #include <iostream><br>
              #include "mpi.h"<br>
              #include "vtkStructuredGrid.h"<br>
              #include "vtkSmartPointer.h"<br>
              #include "vtkMultiBlockDataSet.h"<br>
              #include "vtkMultiPieceDataSet.h"<br>
              #include "vtkNew.h"<br>
              #include "vtkXMLPMultiBlockDataWriter.h<wbr>"<br>
              #include "vtkMPIController.h"<br>
              #include "vtkDoubleArray.h"<br>
              #include "vtkDataSet.h"<br>
              #include "vtkPointData.h"<br>
              #include "vtkInformation.h"<br>
              #include "vtkStreamingDemandDrivenPipel<wbr>ine.h"<br>
              <br>
              using namespace std;<br>
              <br>
              // global variables<br>
              int const N = 10;<br>
              int const tileX = 2;<br>
              int const tileY = 2;<br>
              <br>
              int main(int argc, char** argv) {<br>
                // initialize MPI<br>
                MPI_Init(NULL, NULL);<br>
              <br>
                // setup global controller, true means initialized
              externally<br>
                vtkNew<vtkMPIController> controller;<br>
                controller->Initialize(&argc, &argv, true);<br>
              vtkMultiProcessController::Set<wbr>GlobalController(controller.<wbr>Get());<br>
              <br>
                // get number of processor and rank<br>
                int nproc, myid;<br>
                MPI_Comm_size(MPI_COMM_WORLD, &nproc);<br>
                MPI_Comm_rank(MPI_COMM_WORLD, &myid);<br>
              <br>
                // create 2d decomposition parameters<br>
                bool hasLeft = false, hasTop = false;<br>
                if (myid%tileX == 0) hasLeft = true;<br>
                if (myid/tileY == 0) hasTop = true;<br>
              <br>
                int istr = myid%tileX*N/tileX;<br>
                int iend = istr+N/tileX-1;<br>
                if (hasLeft) iend = iend+1;<br>
                int jstr = myid/tileX*N/tileY;<br>
                int jend = jstr+N/tileY-1;<br>
                if (hasTop) jend = jend+1;<br>
                int size = (iend-istr+1)*(jend-jstr+1);<br>
                //cout << myid << " of " << nproc
              << " " << boolalpha << hasLeft <<
              " " << hasTop << endl;<br>
                //cout << myid << " " << istr <<
              " " << iend << " " << jstr << " "
              << jend << endl;<br>
              <br>
                // create structured grid<br>
                vtkSmartPointer<vtkStructuredG<wbr>rid> grid =
              vtkSmartPointer<vtkStructuredG<wbr>rid>::New();<br>
              <br>
                // define grid coordinates<br>
                vtkSmartPointer<vtkPoints> points =
              vtkSmartPointer<vtkPoints>::Ne<wbr>w();<br>
                points->SetNumberOfPoints(size<wbr>);<br>
              <br>
                int id = 0;<br>
                for (int i = istr; i <= iend; i++) {<br>
                  for (int j = jstr; j <= jend; j++) {<br>
                    points->SetPoint(id, i, j, 0);<br>
                    id = id+1;<br>
                  }<br>
                }<br>
              <br>
                grid->SetPoints(points);<br>
                grid->SetExtent(istr, iend, jstr, jend, 0, 0);<br>
                points->Delete();<br>
              <br>
                // add field, checkerboard type data with full of 0/1<br>
                vtkSmartPointer<vtkDoubleArray<wbr>> field =
              vtkSmartPointer<vtkDoubleArray<wbr>>::New();<br>
                field->SetName("dummy");<br>
                field->SetNumberOfComponents(1<wbr>);<br>
                field->SetNumberOfValues(size)<wbr>;<br>
              <br>
                id = 0;<br>
                for (int i = istr; i <= iend; i++) {<br>
                  for (int j = jstr; j <= jend; j++) {<br>
                    int val = (i+j)%2 ? 0 : 1;<br>
                    field->SetValue(id, val);<br>
                    id = id+1;<br>
                  }<br>
                }<br>
              <br>
                // create multi block grid<br>
                vtkMultiBlockDataSet* mb = vtkMultiBlockDataSet::New();<br>
                vtkNew<vtkMultiPieceDataSet> mpds;<br>
                mpds->SetNumberOfPieces(nproc)<wbr>;<br>
                mpds->SetPiece(myid, grid.GetPointer());<br>
                mb->SetNumberOfBlocks(1);<br>
                mb->SetBlock(0, mpds.GetPointer());<br>
              <br>
                // add field<br>
                vtkDataSet *ds = vtkDataSet::SafeDownCast(mpds-<wbr>>GetPiece(myid));<br>
                ds->GetPointData()->AddArray(f<wbr>ield);<br>
              <br>
                // create a writer<br>
                vtkSmartPointer<vtkXMLPMultiBl<wbr>ockDataWriter>
              mbw = vtkSmartPointer<vtkXMLPMultiBl<wbr>ockDataWriter>::New();<br>
                mbw->SetFileName("test.vtm");<br>
                mbw->SetDataModeToAscii();<br>
                mbw->SetInputData(mb);<br>
                mbw->Write();<br>
              <br>
                // finalize the MPI environment.<br>
                MPI_Finalize();<br>
              <br>
                return 0;<br>
              }<br>
              <br>
              ______________________________<wbr>_________________<br>
              Powered by <a moz-do-not-send="true"
                href="http://www.kitware.com" rel="noreferrer"
                target="_blank">www.kitware.com</a><br>
              <br>
              Visit other Kitware open-source projects at <a
                moz-do-not-send="true"
                href="http://www.kitware.com/opensource/opensource.html"
                rel="noreferrer" target="_blank">http://www.kitware.com/<wbr>opensource/opensource.html</a><br>
              <br>
              Please keep messages on-topic and check the VTK FAQ at: <a
                moz-do-not-send="true"
                href="http://www.vtk.org/Wiki/VTK_FAQ" rel="noreferrer"
                target="_blank">http://www.vtk.org/Wiki/VTK_<wbr>FAQ</a><br>
              <br>
              Search the list archives at: <a moz-do-not-send="true"
                href="http://markmail.org/search/?q=vtkusers"
                rel="noreferrer" target="_blank">http://markmail.org/search/?q=<wbr>vtkusers</a><br>
              <br>
              Follow this link to subscribe/unsubscribe:<br>
              <a moz-do-not-send="true"
                href="http://public.kitware.com/mailman/listinfo/vtkusers"
                rel="noreferrer" target="_blank">http://public.kitware.com/<wbr>mailman/listinfo/vtkusers</a><br>
              <br>
            </blockquote>
          </div>
          <br>
        </div>
      </div>
    </blockquote>
    <p><br>
    </p>
  </body>
</html>