1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195
|
/*=========================================================================
Program: Visualization Toolkit
Module: TestIOADIOS2VTX_VTI3DRendering.cxx
-------------------------------------------------------------------------
Copyright 2008 Sandia Corporation.
Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
the U.S. Government retains certain rights in this software.
-------------------------------------------------------------------------
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
/*
* TestIOADIOS2VTX_VTI3DRendering.cxx : simple rendering test
*
* Created on: Jun 19, 2019
* Author: William F Godoy godoywf@ornl.gov
*/
#include "vtkADIOS2VTXReader.h"
#include <numeric>
#include "vtkCamera.h"
#include "vtkCellData.h"
#include "vtkDataArray.h"
#include "vtkDataObject.h"
#include "vtkDataSetMapper.h"
#include "vtkImageData.h"
#include "vtkInformation.h"
#include "vtkLookupTable.h"
#include "vtkMPI.h"
#include "vtkMPICommunicator.h"
#include "vtkMPIController.h"
#include "vtkMultiBlockDataSet.h"
#include "vtkMultiPieceDataSet.h"
#include "vtkMultiProcessController.h"
#include "vtkNew.h"
#include "vtkRenderWindow.h"
#include "vtkRenderWindowInteractor.h"
#include "vtkRenderer.h"
#include "vtkSmartPointer.h"
#include "vtkTesting.h"
#include <adios2.h>
namespace
{
MPI_Comm MPIGetComm()
{
MPI_Comm comm = MPI_COMM_NULL;
vtkMultiProcessController* controller = vtkMultiProcessController::GetGlobalController();
vtkMPICommunicator* vtkComm = vtkMPICommunicator::SafeDownCast(controller->GetCommunicator());
if (vtkComm)
{
if (vtkComm->GetMPIComm())
{
comm = *(vtkComm->GetMPIComm()->GetHandle());
}
}
return comm;
}
int MPIGetRank()
{
MPI_Comm comm = MPIGetComm();
int rank;
MPI_Comm_rank(comm, &rank);
return rank;
}
int MPIGetSize()
{
MPI_Comm comm = MPIGetComm();
int size;
MPI_Comm_size(comm, &size);
return size;
}
std::size_t TotalElements(const std::vector<std::size_t>& dimensions) noexcept
{
return std::accumulate(dimensions.begin(), dimensions.end(), 1, std::multiplies<std::size_t>());
}
void WriteBPFile3DVars(const std::string& fileName, const adios2::Dims& shape,
const adios2::Dims& start, const adios2::Dims& count, const int rank)
{
const size_t totalElements = TotalElements(count);
const std::string extent = "0 " + std::to_string(shape[0]) + " " + "0 " +
std::to_string(shape[1]) + " " + "0 " + std::to_string(shape[2]);
const std::string imageSchema = R"( <?xml version="1.0"?>
<VTKFile type="ImageData" version="0.1" byte_order="LittleEndian">
<ImageData WholeExtent=")" +
extent +
R"(" Origin="0 0 0" Spacing="1 1 1">
<Piece Extent=")" +
extent +
R"(">
<CellData>
<DataArray Name="T" />
<DataArray Name="TIME">
time
</DataArray>
</CellData>
</Piece>
</ImageData>
</VTKFile>)";
// using adios2 C++ high-level API
std::vector<double> T(totalElements);
std::iota(T.begin(), T.end(), static_cast<double>(rank * totalElements));
adios2::fstream fw(fileName, adios2::fstream::out, MPIGetComm());
fw.write_attribute("vtk.xml", imageSchema);
fw.write("time", 0);
fw.write("T", T.data(), shape, start, count);
fw.close();
}
} // end empty namespace
int TestIOADIOS2VTX_VTI3DRendering(int argc, char* argv[])
{
vtkNew<vtkMPIController> mpiController;
mpiController->Initialize(&argc, &argv, 0);
vtkMultiProcessController::SetGlobalController(mpiController);
const int rank = MPIGetRank();
const int size = MPIGetSize();
vtkNew<vtkTesting> testing;
const std::string rootDirectory(testing->GetTempDirectory());
const std::string fileName = rootDirectory + "/heat3D_render.bp";
const adios2::Dims count{ 4, 4, 8 };
const adios2::Dims start{ static_cast<size_t>(rank) * count[0], 0, 0 };
const adios2::Dims shape{ static_cast<size_t>(size) * count[0], count[1], count[2] };
WriteBPFile3DVars(fileName, shape, start, count, rank);
vtkNew<vtkADIOS2VTXReader> adios2Reader;
adios2Reader->SetFileName(fileName.c_str());
adios2Reader->UpdateInformation();
adios2Reader->Update();
vtkMultiBlockDataSet* multiBlock = adios2Reader->GetOutput();
vtkMultiPieceDataSet* mp = vtkMultiPieceDataSet::SafeDownCast(multiBlock->GetBlock(0));
vtkImageData* imageData = vtkImageData::SafeDownCast(mp->GetPiece(rank));
// set color table
vtkSmartPointer<vtkLookupTable> lookupTable = vtkSmartPointer<vtkLookupTable>::New();
lookupTable->SetNumberOfTableValues(10);
lookupTable->SetRange(0.0, 1.0);
lookupTable->Build();
// render imageData
vtkSmartPointer<vtkDataSetMapper> mapper = vtkSmartPointer<vtkDataSetMapper>::New();
mapper->SetInputData(imageData);
mapper->SetLookupTable(lookupTable);
// mapper->SelectColorArray("T");
mapper->SetScalarModeToUseCellFieldData();
vtkSmartPointer<vtkActor> actor = vtkSmartPointer<vtkActor>::New();
actor->SetMapper(mapper);
vtkSmartPointer<vtkRenderWindow> renderWindow = vtkSmartPointer<vtkRenderWindow>::New();
vtkSmartPointer<vtkRenderer> renderer = vtkSmartPointer<vtkRenderer>::New();
// Add both renderers to the window
renderWindow->AddRenderer(renderer);
renderer->AddActor(actor);
renderer->ResetCamera();
vtkSmartPointer<vtkRenderWindowInteractor> renderWindowInteractor =
vtkSmartPointer<vtkRenderWindowInteractor>::New();
renderWindowInteractor->SetRenderWindow(renderWindow);
renderWindow->Render();
// renderWindowInteractor->Start();
mpiController->Finalize();
return EXIT_SUCCESS;
}
|