[vtkusers] problem in writing data with vtkXMLMultiBlockDataWriter ...
Utkarsh Ayachit
utkarsh.ayachit at kitware.com
Wed May 31 10:13:04 EDT 2017
Use vtkXMLPMultiBlockDataWriter
(http://www.vtk.org/doc/nightly/html/classvtkXMLPMultiBlockDataWriter.html)
instead. For it to work correctly, you'll need to setup a global
vtkMPIController. Since you're calling MPI_Init() youself, you can do
something like following after MPI_Init().
vtkNew<vtkMPIController> controller;
controller->Initialize(argv, argc, /*initializedExternally*/ true);
vtkMultiProcessController::SetGlobalController(controller.Get());
On Wed, May 31, 2017 at 10:00 AM, Ufuk Turuncoglu
<ufuk.turuncoglu at itu.edu.tr> wrote:
> Hi,
>
> I am trying to create very simple example of creating vtkMultiBlockDataSet
> using vtkStructuredGrid to sole a problem of mine in ParaView Catalyst. By
> this way, i could test some ideas without overhead of simulation code. In
> this case, i wrote following C++ code to demonstrate it. The problem is that
> the result test.vtm file includes information only one piece of data,
>
> <?xml version="1.0"?>
> <VTKFile type="vtkMultiBlockDataSet" version="1.0" byte_order="LittleEndian"
> header_type="UInt32" compressor="vtkZLibDataCompressor">
> <vtkMultiBlockDataSet>
> <Piece index="0">
> <DataSet index="1" file="test/test_1.vts">
> </DataSet>
> </Piece>
> </vtkMultiBlockDataSet>
> </VTKFile>
>
> but if i run the code with 4 core the test/ directory contains four file. I
> also observed that the file name written to test.vtm is changing based on
> the fastest core. So, i just wonder that am i missing something in here.
>
> Thanks,
> Regards,
>
> --ufuk
>
> #include <iostream>
> #include "mpi.h"
> #include "vtkStructuredGrid.h"
> #include "vtkSmartPointer.h"
> #include "vtkMultiBlockDataSet.h"
> #include "vtkMultiPieceDataSet.h"
> #include "vtkNew.h"
> #include "vtkXMLMultiBlockDataWriter.h"
>
> using namespace std;
>
> // global variables
> int const N = 10;
> int const tileX = 2;
> int const tileY = 2;
>
> int main(int argc, char** argv) {
> // initialize MPI
> MPI_Init(NULL, NULL);
>
> // get number of processor and rank
> int nproc, myid;
> MPI_Comm_size(MPI_COMM_WORLD, &nproc);
> MPI_Comm_rank(MPI_COMM_WORLD, &myid);
>
> // create 2d decomposition parameters
> bool hasLeft = false, hasTop = false;
> if (myid%tileX == 0) hasLeft = true;
> if (myid/tileY == 0) hasTop = true;
> //cout << myid << " of " << nproc << " " << boolalpha << hasLeft << " " <<
> hasTop << endl;
>
> // create structured grid
> vtkSmartPointer<vtkStructuredGrid> grid =
> vtkSmartPointer<vtkStructuredGrid>::New();
>
> // define grid coordinates
> vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
> points->SetNumberOfPoints(N/tileX*N/tileY);
>
> int istr = myid%tileX*N/tileX;
> int iend = istr+N/tileX-1;
> int jstr = myid/tileX*N/tileY;
> int jend = jstr+N/tileY-1;
> //cout << myid << " " << istr << " " << iend << " " << jstr << " " << jend
> << endl;
>
> int id = 0;
> for (int i = istr; i < iend; i++) {
> for (int j = jstr; j < jend; j++) {
> points->SetPoint(id, i, j, 0);
> id = id+1;
> }
> }
>
> grid->SetPoints(points);
> grid->SetExtent(istr, iend, jstr, jend, 0, 0);
> points->Delete();
>
> // create multi block grid
> vtkMultiBlockDataSet* mb = vtkMultiBlockDataSet::New();
> vtkNew<vtkMultiPieceDataSet> mpds;
> mpds->SetNumberOfPieces(nproc);
> mpds->SetPiece(myid, grid.GetPointer());
> mb->SetNumberOfBlocks(1);
> mb->SetBlock(0, mpds.GetPointer());
>
> // create a writer
> vtkSmartPointer<vtkXMLMultiBlockDataWriter> mbw =
> vtkSmartPointer<vtkXMLMultiBlockDataWriter>::New();
> mbw->SetFileName("test.vtm");
> mbw->SetDataModeToAscii();
> mbw->SetInputData(mb);
> mbw->Write();
>
> // finalize the MPI environment.
> MPI_Finalize();
>
> return 0;
> }
>
> _______________________________________________
> 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 VTK FAQ at:
> http://www.vtk.org/Wiki/VTK_FAQ
>
> Search the list archives at: http://markmail.org/search/?q=vtkusers
>
> Follow this link to subscribe/unsubscribe:
> http://public.kitware.com/mailman/listinfo/vtkusers
More information about the vtkusers
mailing list