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 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253
|
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
// This example demonstrates the use of data parallelism in VTK. The
// pipeline ( vtkMPIImageReader -> vtkContourFilter -> vtkElevationFilter )
// is created in parallel and each process is assigned 1 piece to process.
// All satellite processes send the result to the first process which
// collects and renders them.
#include <vtkMPI.h>
#include "vtkActor.h"
#include "vtkAppendPolyData.h"
#include "vtkCamera.h"
#include "vtkConeSource.h"
#include "vtkContourFilter.h"
#include "vtkDataSet.h"
#include "vtkElevationFilter.h"
#include "vtkImageData.h"
#include "vtkInformation.h"
#include "vtkMPIController.h"
#include "vtkMPIImageReader.h"
#include "vtkMath.h"
#include "vtkPolyData.h"
#include "vtkPolyDataMapper.h"
#include "vtkRegressionTestImage.h"
#include "vtkRenderWindow.h"
#include "vtkRenderWindowInteractor.h"
#include "vtkRenderer.h"
#include "vtkStreamingDemandDrivenPipeline.h"
#include "vtkTestUtilities.h"
#include "vtkWindowToImageFilter.h"
#include "vtkDebugLeaks.h"
namespace
{
const float ISO_START = 4250.0;
const float ISO_STEP = -1250.0;
const int ISO_NUM = 3;
// Just pick a tag which is available
const int ISO_VALUE_RMI_TAG = 300;
const int ISO_OUTPUT_TAG = 301;
struct ParallelIsoArgs_tmp
{
int* retVal;
int argc;
char** argv;
};
struct ParallelIsoRMIArgs_tmp
{
vtkContourFilter* ContourFilter;
vtkMultiProcessController* Controller;
vtkElevationFilter* Elevation;
};
// call back to set the iso surface value.
void SetIsoValueRMI(
void* localArg, void* vtkNotUsed(remoteArg), int vtkNotUsed(remoteArgLen), int vtkNotUsed(id))
{
ParallelIsoRMIArgs_tmp* args = (ParallelIsoRMIArgs_tmp*)localArg;
float val;
vtkMultiProcessController* controller = args->Controller;
int myid = controller->GetLocalProcessId();
int numProcs = controller->GetNumberOfProcesses();
vtkContourFilter* iso = args->ContourFilter;
val = iso->GetValue(0);
iso->SetValue(0, val + ISO_STEP);
args->Elevation->UpdatePiece(myid, numProcs, 0);
controller->Send(args->Elevation->GetOutput(), 0, ISO_OUTPUT_TAG);
}
// This will be called by all processes
void MyMain(vtkMultiProcessController* controller, void* arg)
{
vtkMPIImageReader* reader;
vtkContourFilter* iso;
vtkElevationFilter* elev;
int myid, numProcs;
float val;
ParallelIsoArgs_tmp* args = reinterpret_cast<ParallelIsoArgs_tmp*>(arg);
// Obtain the id of the running process and the total
// number of processes
myid = controller->GetLocalProcessId();
numProcs = controller->GetNumberOfProcesses();
// Create the reader, the data file name might have
// to be changed depending on where the data files are.
char* fname = vtkTestUtilities::ExpandDataFileName(args->argc, args->argv, "Data/headsq/quarter");
reader = vtkMPIImageReader::New();
reader->SetDataByteOrderToLittleEndian();
reader->SetDataExtent(0, 63, 0, 63, 1, 93);
reader->SetFilePrefix(fname);
reader->SetDataSpacing(3.2, 3.2, 1.5);
delete[] fname;
// Iso-surface.
iso = vtkContourFilter::New();
iso->SetInputConnection(reader->GetOutputPort());
iso->SetValue(0, ISO_START);
iso->ComputeScalarsOff();
iso->ComputeGradientsOff();
// Compute a different color for each process.
elev = vtkElevationFilter::New();
elev->SetInputConnection(iso->GetOutputPort());
val = (myid + 1) / static_cast<float>(numProcs);
elev->SetScalarRange(val, val + 0.001);
// Make sure all processes update at the same time.
elev->UpdatePiece(myid, numProcs, 0);
if (myid != 0)
{
// If I am not the root process
ParallelIsoRMIArgs_tmp args2;
args2.ContourFilter = iso;
args2.Controller = controller;
args2.Elevation = elev;
// Last, set up a RMI call back to change the iso surface value.
// This is done so that the root process can let this process
// know that it wants the contour value to change.
controller->AddRMI(SetIsoValueRMI, (void*)&args2, ISO_VALUE_RMI_TAG);
controller->ProcessRMIs();
}
else
{
// Create the rendering part of the pipeline
vtkAppendPolyData* app = vtkAppendPolyData::New();
vtkRenderer* ren = vtkRenderer::New();
vtkRenderWindow* renWindow = vtkRenderWindow::New();
vtkRenderWindowInteractor* iren = vtkRenderWindowInteractor::New();
vtkPolyDataMapper* mapper = vtkPolyDataMapper::New();
vtkActor* actor = vtkActor::New();
vtkCamera* cam = vtkCamera::New();
renWindow->AddRenderer(ren);
iren->SetRenderWindow(renWindow);
ren->SetBackground(0.9, 0.9, 0.9);
renWindow->SetSize(400, 400);
mapper->SetInputConnection(app->GetOutputPort());
actor->SetMapper(mapper);
ren->AddActor(actor);
cam->SetFocalPoint(100, 100, 65);
cam->SetPosition(100, 450, 65);
cam->SetViewUp(0, 0, -1);
cam->SetViewAngle(30);
cam->SetClippingRange(177.0, 536.0);
ren->SetActiveCamera(cam);
// loop through some iso surface values.
for (int j = 0; j < ISO_NUM; ++j)
{
for (int i = 1; i < numProcs; ++i)
{
// trigger the RMI to change the iso surface value.
controller->TriggerRMI(i, ISO_VALUE_RMI_TAG);
}
// set the local value
iso->SetValue(0, iso->GetValue(0) + ISO_STEP);
elev->UpdatePiece(myid, numProcs, 0);
for (int i = 1; i < numProcs; ++i)
{
vtkPolyData* pd = vtkPolyData::New();
controller->Receive(pd, i, ISO_OUTPUT_TAG);
if (j == ISO_NUM - 1)
{
app->AddInputData(pd);
}
pd->Delete();
}
}
// Tell the other processors to stop processing RMIs.
for (int i = 1; i < numProcs; ++i)
{
controller->TriggerRMI(i, vtkMultiProcessController::BREAK_RMI_TAG);
}
vtkPolyData* outputCopy = vtkPolyData::New();
outputCopy->ShallowCopy(elev->GetOutput());
app->AddInputData(outputCopy);
outputCopy->Delete();
app->Update();
outputCopy->RemoveGhostCells();
renWindow->Render();
*(args->retVal) = vtkRegressionTester::Test(args->argc, args->argv, renWindow, 10);
if (*(args->retVal) == vtkRegressionTester::DO_INTERACTOR)
{
iren->Start();
}
// Clean up
app->Delete();
ren->Delete();
renWindow->Delete();
iren->Delete();
mapper->Delete();
actor->Delete();
cam->Delete();
}
// clean up objects in all processes.
reader->Delete();
iso->Delete();
elev->Delete();
}
} // end anon namespace
int ParallelIso2(int argc, char* argv[])
{
// This is here to avoid false leak messages from vtkDebugLeaks when
// using mpich. It appears that the root process which spawns all the
// main processes waits in MPI_Init() and calls exit() when
// the others are done, causing apparent memory leaks for any objects
// created before MPI_Init().
MPI_Init(&argc, &argv);
// Note that this will create a vtkMPIController if MPI
// is configured, vtkThreadedController otherwise.
vtkMPIController* controller = vtkMPIController::New();
controller->Initialize(&argc, &argv, 1);
// Added for regression test.
// ----------------------------------------------
int retVal = 1;
ParallelIsoArgs_tmp args;
args.retVal = &retVal;
args.argc = argc;
args.argv = argv;
// ----------------------------------------------
controller->SetSingleMethod(MyMain, &args);
controller->SingleMethodExecute();
controller->Finalize();
controller->Delete();
return !retVal;
}
|