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
|
#include <iostream>
#include "FEAdaptor.h"
#include "FEDataStructures.h"
#include <vtkCellData.h>
#include <vtkCellType.h>
#include <vtkCPDataDescription.h>
#include <vtkCPInputDataDescription.h>
#include <vtkCPProcessor.h>
#include <vtkCPPythonScriptPipeline.h>
#include <vtkDoubleArray.h>
#include <vtkFloatArray.h>
#include <vtkNew.h>
#include <vtkPoints.h>
#include <vtkPointData.h>
#include <vtkUnstructuredGrid.h>
namespace
{
vtkCPProcessor* Processor = NULL;
void BuildVTKGrid(Grid& grid, vtkUnstructuredGrid* vtkgrid)
{
// create the points information
vtkNew<vtkDoubleArray> pointArray;
pointArray->SetNumberOfComponents(3);
pointArray->SetArray(grid.GetPointsArray(), static_cast<vtkIdType>(grid.GetNumberOfPoints()*3), 1);
vtkNew<vtkPoints> points;
points->SetData(pointArray.GetPointer());
vtkgrid->SetPoints(points.GetPointer());
// create the cells
size_t numCells = grid.GetNumberOfCells();
vtkgrid->Allocate(static_cast<vtkIdType>(numCells*9));
for(size_t cell=0;cell<numCells;cell++)
{
unsigned int* cellPoints = grid.GetCellPoints(cell);
vtkIdType tmp[8] = {cellPoints[0], cellPoints[1], cellPoints[2], cellPoints[3],
cellPoints[4], cellPoints[5], cellPoints[6], cellPoints[7]};
vtkgrid->InsertNextCell(VTK_HEXAHEDRON, 8, tmp);
}
}
void UpdateVTKAttributes(Grid& grid, Attributes& attributes, vtkUnstructuredGrid* vtkgrid)
{
if(vtkgrid->GetPointData()->GetNumberOfArrays() == 0)
{
// velocity array
vtkNew<vtkDoubleArray> velocity;
velocity->SetName("velocity");
velocity->SetNumberOfComponents(3);
velocity->SetNumberOfTuples(static_cast<vtkIdType>(grid.GetNumberOfPoints()));
vtkgrid->GetPointData()->AddArray(velocity.GetPointer());
}
if(vtkgrid->GetCellData()->GetNumberOfArrays() == 0)
{
// pressure array
vtkNew<vtkFloatArray> pressure;
pressure->SetName("pressure");
pressure->SetNumberOfComponents(1);
vtkgrid->GetCellData()->AddArray(pressure.GetPointer());
}
vtkDoubleArray* velocity = vtkDoubleArray::SafeDownCast(
vtkgrid->GetPointData()->GetArray("velocity"));
// The velocity array is ordered as vx0,vx1,vx2,..,vy0,vy1,vy2,..,vz0,vz1,vz2,..
// so we need to create a full copy of it with VTK's ordering of
// vx0,vy0,vz0,vx1,vy1,vz1,..
double* velocityData = attributes.GetVelocityArray();
vtkIdType numTuples = velocity->GetNumberOfTuples();
for(vtkIdType i=0;i<numTuples;i++)
{
double values[3] = {velocityData[i], velocityData[i+numTuples],
velocityData[i+2*numTuples]};
velocity->SetTypedTuple(i, values);
}
vtkFloatArray* pressure = vtkFloatArray::SafeDownCast(
vtkgrid->GetCellData()->GetArray("pressure"));
// The pressure array is a scalar array so we can reuse
// memory as long as we ordered the points properly.
float* pressureData = attributes.GetPressureArray();
pressure->SetArray(pressureData, static_cast<vtkIdType>(grid.GetNumberOfCells()), 1);
}
void BuildVTKDataStructures(Grid& grid, Attributes& attributes, vtkUnstructuredGrid* vtkgrid)
{
BuildVTKGrid(grid, vtkgrid);
UpdateVTKAttributes(grid, attributes, vtkgrid);
}
}
namespace FEAdaptor
{
void Initialize(std::vector<std::string>& scripts)
{
if(Processor == NULL)
{
Processor = vtkCPProcessor::New();
Processor->Initialize();
}
else
{
Processor->RemoveAllPipelines();
}
for(std::vector<std::string>::iterator it=scripts.begin();it!=scripts.end();it++)
{
vtkNew<vtkCPPythonScriptPipeline> pipeline;
pipeline->Initialize(it->c_str());
Processor->AddPipeline(pipeline.GetPointer());
}
}
void Finalize()
{
if(Processor)
{
Processor->Delete();
Processor = NULL;
}
}
void CoProcess(Grid& grid, Attributes& attributes, double time,
unsigned int timeStep, bool lastTimeStep)
{
vtkNew<vtkCPDataDescription> dataDescription;
dataDescription->AddInput("input");
dataDescription->SetTimeData(time, timeStep);
if(lastTimeStep == true)
{
// assume that we want to all the pipelines to execute if it
// is the last time step.
dataDescription->ForceOutputOn();
}
if(Processor->RequestDataDescription(dataDescription.GetPointer()) != 0)
{
vtkNew<vtkUnstructuredGrid> vtkgrid;
BuildVTKDataStructures(grid, attributes, vtkgrid.GetPointer());
dataDescription->GetInputDescriptionByName("input")->SetGrid(vtkgrid.GetPointer());
Processor->CoProcess(dataDescription.GetPointer());
}
}
} // end of Catalyst namespace
|