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
|
#include "vtkCPVTKPipeline.h"
#include <vtkCommunicator.h>
#include <vtkCompleteArrays.h>
#include <vtkCPDataDescription.h>
#include <vtkDataArray.h>
#include <vtkCPInputDataDescription.h>
#include <vtkMultiProcessController.h>
#include <vtkNew.h>
#include <vtkObjectFactory.h>
#include <vtkPointData.h>
#include <vtkPVArrayCalculator.h>
#include <vtkPVTrivialProducer.h>
#include <vtkSMProxyManager.h>
#include <vtkThreshold.h>
#include <vtkUnstructuredGrid.h>
#include <vtkXMLPUnstructuredGridWriter.h>
#include <string>
#include <sstream>
vtkStandardNewMacro(vtkCPVTKPipeline);
//----------------------------------------------------------------------------
vtkCPVTKPipeline::vtkCPVTKPipeline()
{
this->OutputFrequency = 0;
}
//----------------------------------------------------------------------------
vtkCPVTKPipeline::~vtkCPVTKPipeline()
{
}
//----------------------------------------------------------------------------
void vtkCPVTKPipeline::Initialize(int outputFrequency, std::string& fileName)
{
this->OutputFrequency = outputFrequency;
this->FileName = fileName;
}
//----------------------------------------------------------------------------
int vtkCPVTKPipeline::RequestDataDescription(
vtkCPDataDescription* dataDescription)
{
if(!dataDescription)
{
vtkWarningMacro("dataDescription is NULL.");
return 0;
}
if(this->FileName.empty())
{
vtkWarningMacro("No output file name given to output results to.");
return 0;
}
if(dataDescription->GetForceOutput() == true ||
(this->OutputFrequency != 0 && dataDescription->GetTimeStep() % this->OutputFrequency == 0) )
{
dataDescription->GetInputDescriptionByName("input")->AllFieldsOn();
dataDescription->GetInputDescriptionByName("input")->GenerateMeshOn();
return 1;
}
return 0;
}
//----------------------------------------------------------------------------
int vtkCPVTKPipeline::CoProcess(
vtkCPDataDescription* dataDescription)
{
if(!dataDescription)
{
vtkWarningMacro("DataDescription is NULL");
return 0;
}
vtkUnstructuredGrid* grid = vtkUnstructuredGrid::SafeDownCast(
dataDescription->GetInputDescriptionByName("input")->GetGrid());
if(grid == NULL)
{
vtkWarningMacro("DataDescription is missing input unstructured grid.");
return 0;
}
if(this->RequestDataDescription(dataDescription) == 0)
{
return 1;
}
vtkNew<vtkPVTrivialProducer> producer;
producer->SetOutput(grid);
vtkNew<vtkPVArrayCalculator> calculator;
calculator->SetInputConnection(producer->GetOutputPort());
calculator->SetAttributeMode(1);
calculator->SetResultArrayName("velocity magnitude");
calculator->SetFunction("mag(velocity)");
// update now so that we can get the global data bounds of
// the velocity magnitude for thresholding
calculator->Update();
double range[2];
vtkUnstructuredGrid::SafeDownCast(calculator->GetOutput())->GetPointData()
->GetArray("velocity magnitude")->GetRange(range, 0);
double globalRange[2];
vtkMultiProcessController::GetGlobalController()->AllReduce(
range+1, globalRange+1, 1, vtkCommunicator::MAX_OP);
vtkNew<vtkThreshold> threshold;
threshold->SetInputConnection(calculator->GetOutputPort());
threshold->SetInputArrayToProcess(
0, 0, 0, "vtkDataObject::FIELD_ASSOCIATION_POINTS", "velocity magnitude");
threshold->ThresholdBetween(0.9*globalRange[1], globalRange[1]);
// If process 0 doesn't have any points or cells, the writer may
// have problems in parallel so we use completeArrays to fill in
// the missing information.
vtkNew<vtkCompleteArrays> completeArrays;
completeArrays->SetInputConnection(threshold->GetOutputPort());
vtkNew<vtkXMLPUnstructuredGridWriter> writer;
writer->SetInputConnection(completeArrays->GetOutputPort());
std::ostringstream o;
o << dataDescription->GetTimeStep();
std::string name = this->FileName + o.str() + ".pvtu";
writer->SetFileName(name.c_str());
writer->Update();
return 1;
}
//----------------------------------------------------------------------------
void vtkCPVTKPipeline::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent << "OutputFrequency: " << this->OutputFrequency << "\n";
os << indent << "FileName: " << this->FileName << "\n";
}
|