vtk-m2/examples/particle_advection/ParticleAdvection.cxx

85 lines
2.8 KiB
C++
Raw Normal View History

//============================================================================
// Copyright (c) Kitware, Inc.
// All rights reserved.
// See LICENSE.txt for details.
2019-04-15 23:24:21 +00:00
//
// 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.
//============================================================================
#include <vtkm/cont/DataSet.h>
2019-08-28 11:48:20 +00:00
#include <vtkm/filter/Streamline.h>
#include <vtkm/io/reader/VTKDataSetReader.h>
#include <vtkm/io/writer/VTKDataSetWriter.h>
2019-08-28 15:14:30 +00:00
// Example computing streamlines.
// An example vector field is available in the vtk-m data directory: magField.vtk
// Example usage:
// this will advect 200 particles 50 steps using a step size of 0.01
//
// Particle_Advection <path-to-data-dir>/magField.vtk vec 200 50 0.01 output.vtk
//
2019-08-28 11:48:20 +00:00
int main(int argc, char** argv)
{
2019-08-28 11:48:20 +00:00
if (argc != 7)
{
2019-08-28 11:48:20 +00:00
std::cerr << "Usage: " << argv[0] << " dataFile varName numSeeds numSteps stepSize outputFile"
<< std::endl;
return -1;
}
2019-08-28 11:48:20 +00:00
std::string dataFile = argv[1];
std::string varName = argv[2];
vtkm::Id numSeeds = std::stoi(argv[3]);
vtkm::Id numSteps = std::stoi(argv[4]);
vtkm::FloatDefault stepSize = std::stof(argv[5]);
std::string outputFile = argv[6];
2019-08-28 11:48:20 +00:00
vtkm::cont::DataSet ds;
2019-08-28 15:14:30 +00:00
if (dataFile.find(".vtk") != std::string::npos)
{
2019-08-28 11:48:20 +00:00
vtkm::io::reader::VTKDataSetReader rdr(dataFile);
ds = rdr.ReadDataSet();
}
2019-08-28 11:48:20 +00:00
else
{
2019-08-28 11:48:20 +00:00
std::cerr << "Unsupported data file: " << dataFile << std::endl;
return -1;
}
2019-08-28 15:14:30 +00:00
//create seeds randomly placed withing the bounding box of the data.
2019-08-28 11:48:20 +00:00
vtkm::Bounds bounds = ds.GetCoordinateSystem().GetBounds();
std::vector<vtkm::Vec3f> seeds;
2019-08-28 11:48:20 +00:00
for (int i = 0; i < numSeeds; i++)
{
2019-08-28 11:48:20 +00:00
vtkm::Vec3f p;
vtkm::FloatDefault rx = (vtkm::FloatDefault)rand() / (vtkm::FloatDefault)RAND_MAX;
vtkm::FloatDefault ry = (vtkm::FloatDefault)rand() / (vtkm::FloatDefault)RAND_MAX;
vtkm::FloatDefault rz = (vtkm::FloatDefault)rand() / (vtkm::FloatDefault)RAND_MAX;
p[0] = static_cast<vtkm::FloatDefault>(bounds.X.Min + rx * bounds.X.Length());
p[1] = static_cast<vtkm::FloatDefault>(bounds.Y.Min + ry * bounds.Y.Length());
p[2] = static_cast<vtkm::FloatDefault>(bounds.Z.Min + rz * bounds.Z.Length());
seeds.push_back(p);
}
2019-08-28 11:48:20 +00:00
vtkm::cont::ArrayHandle<vtkm::Vec3f> seedArray = vtkm::cont::make_ArrayHandle(seeds);
2019-08-28 11:48:20 +00:00
//compute streamlines
vtkm::filter::Streamline streamline;
2019-08-28 11:48:20 +00:00
streamline.SetStepSize(stepSize);
streamline.SetNumberOfSteps(numSteps);
streamline.SetSeeds(seedArray);
2019-08-28 11:48:20 +00:00
streamline.SetActiveField(varName);
auto output = streamline.Execute(ds);
2019-08-28 11:48:20 +00:00
vtkm::io::writer::VTKDataSetWriter wrt(outputFile);
wrt.WriteDataSet(output);
return 0;
}