1 #include <vtkPolyDataReader.h> 
    2 #include <vtkPolyData.h> 
    3 #include <vtkPointData.h> 
    5 #include <vtkCellArray.h> 
    6 #include <vtkCellDataToPointData.h> 
    7 #include <vtkContourFilter.h> 
    9 #include <boost/unordered_set.hpp> 
   10 #include <boost/program_options.hpp> 
   11 namespace po = boost::program_options;
 
   34     Vertex(
double pX, 
double pY, 
double pZ, 
double pScalar, 
double factor)
 
   35         : x((int)floor(pX*factor)),
 
   36           y((int)floor(pY*factor)),
 
   37           z((int)floor(pZ*factor)),
 
   47         return (x == v.
x && y == v.
y && z == v.
z);
 
   55     return v1->operator==(*v2);
 
   62 struct VertexHash : std::unary_function<VertexSharedPtr, std::size_t>
 
   67         boost::hash_combine(seed, p -> x);
 
   68         boost::hash_combine(seed, p -> y);
 
   69         boost::hash_combine(seed, p -> z);
 
   75 typedef boost::unordered_set<VertexSharedPtr, VertexHash> 
VertexSet;
 
   83 int main(
int argc, 
char* argv[])
 
   86     po::options_description desc(
"Available options");
 
   88         (
"help,h",         
"Produce this help message.")
 
   89         (
"name,n", po::value<string>()->default_value(
"Intensity"),
 
   90                 "Name of field in VTK file to use for intensity.")
 
   91         (
"outname,m", po::value<string>()->default_value(
"intensity"),
 
   92                 "Name of field in output FLD file.")
 
   93         (
"precision,p",  po::value<double>()->default_value(1),
 
   94              "Precision of vertex matching.");
 
   96     po::options_description hidden(
"Hidden options");
 
   98         (
"file",   po::value<vector<string> >(), 
"Input filename");
 
  100     po::options_description cmdline_options;
 
  101     cmdline_options.add(desc).add(hidden);
 
  103     po::positional_options_description p;
 
  106     po::variables_map vm;
 
  111         po::store(po::command_line_parser(argc, argv).
 
  112                   options(cmdline_options).positional(p).run(), vm);
 
  115     catch (
const std::exception& e)
 
  117         cerr << e.what() << endl;
 
  122     if ( vm.count(
"help") || vm.count(
"file") == 0 ||
 
  123                              vm[
"file"].as<vector<string> >().size() != 3) {
 
  124         cerr << 
"Usage: VtkToFld session.xml intensity.vtk output.fld [options]" 
  131     std::vector<std::string> vFiles = vm[
"file"].as<vector<string> >();
 
  132     const string infile  = vFiles[1];
 
  133     const string outfile = vFiles[2];
 
  134     const double factor  = vm[
"precision"].as<
double>();
 
  135     const string name    = vm[
"name"].as<
string>();
 
  136     const string outname = vm[
"outname"].as<
string>();
 
  138     std::vector<std::string> vFilenames;
 
  143     vFilenames.push_back(vFiles[0]);
 
  144     vSession = LibUtilities::SessionReader::CreateInstance(2, argv, vFilenames);
 
  162         int coordim = Exp->GetCoordim(0);
 
  163         int nq      = Exp->GetNpoints();
 
  172             Exp->GetCoords(xc0,xc1);
 
  175             Exp->GetCoords(xc0,xc1,xc2);
 
  178             ASSERTL0(
false,
"Coordim not valid");
 
  183         vtkPolyDataReader *vtkMeshReader = vtkPolyDataReader::New();
 
  184         vtkMeshReader->SetFileName(infile.c_str());
 
  185         vtkMeshReader->Update();
 
  187         vtkPolyData *vtkMesh = vtkMeshReader->GetOutput();
 
  188         vtkCellDataToPointData* c2p = vtkCellDataToPointData::New();
 
  189 #if VTK_MAJOR_VERSION <= 5 
  190         c2p->SetInput(vtkMesh);
 
  192         c2p->SetInputData(vtkMesh);
 
  194         c2p->PassCellDataOn();
 
  196         vtkPolyData *vtkDataAtPoints = c2p->GetPolyDataOutput();
 
  198         vtkPoints *vtkPoints = vtkMesh->GetPoints();
 
  199         ASSERTL0(vtkPoints, 
"ERROR: cannot get points from mesh.");
 
  201         vtkCellArray *vtkPolys = vtkMesh->GetPolys();
 
  202         ASSERTL0(vtkPolys,  
"ERROR: cannot get polygons from mesh.");
 
  204         vtkPointData *vtkPData = vtkDataAtPoints->GetPointData();
 
  205         ASSERTL0(vtkPolys,  
"ERROR: cannot get point data from file.");
 
  215         if (!vtkDataAtPoints->GetPointData()->HasArray(name.c_str())) {
 
  216             n = vtkDataAtPoints->GetPointData()->GetNumberOfArrays();
 
  217             cerr << 
"Input file '" << infile
 
  218                  << 
"' does not have a field named '" 
  219                  << name << 
"'" << endl;
 
  220             cerr << 
"There are " << n << 
" arrays in this file." << endl;
 
  221             for (
int i = 0; i < n; ++i)
 
  224                      << vtkDataAtPoints->GetPointData()->GetArray(i)->GetName()
 
  233         n = vtkPoints->GetNumberOfPoints();
 
  234         for (i = 0; i < n; ++i)
 
  236             vtkPoints->GetPoint(i,p);
 
  237             val = vtkPData->GetScalars(name.c_str())->GetTuple1(i);
 
  238             boost::shared_ptr<Vertex> v(
new Vertex(p[0],p[1],p[2],val,factor));
 
  244         for (i = 0; i < Exp->GetNumElmts(); ++i)
 
  247             for (j = 0; j < e->GetNverts(); ++j)
 
  250                 coeff_idx = Exp->GetCoeff_Offset(i) + e->GetVertexMap(j);
 
  258                 boost::shared_ptr<Vertex> v(
new Vertex(x,y,z,0.0,factor));
 
  259                 vIter = points.find(v);
 
  264                 if (vIter == points.end())
 
  266                     cerr << 
"Vertex " << i << 
" not found. Looking for (" 
  267                             << x << 
", " << y << 
", " << z << 
")" << endl;
 
  271                     Exp->UpdateCoeffs()[coeff_idx] = (*vIter)->scalar;
 
  275         Exp->SetPhysState(
false);
 
  279         std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
 
  280                                                 = Exp->GetFieldDefinitions();
 
  281         std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
 
  283         for(i = 0; i < FieldDef.size(); ++i)
 
  285             FieldDef[i]->m_fields.push_back(outname);
 
  286             Exp->AppendFieldData(FieldDef[i], FieldData[i]);
 
  290         vFld.
Write(outfile, FieldDef, FieldData);
 
  294         cout << 
"An error occurred." << endl;
 
std::size_t operator()(VertexSharedPtr const &p) const 
#define ASSERTL0(condition, msg)
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
boost::unordered_set< VertexSharedPtr, VertexHash > VertexSet
Define a set of Vertices with associated hash function. 
int main(int argc, char *argv[])
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
void GetCoords(Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2=NullNekDouble1DArray, Array< OneD, NekDouble > &coords_3=NullNekDouble1DArray)
this function returns the physical coordinates of the quadrature points of the expansion ...
boost::shared_ptr< Vertex > VertexSharedPtr
Vertex(double pX, double pY, double pZ, double pScalar, double factor)
boost::shared_ptr< ExpList2D > ExpList2DSharedPtr
Shared pointer to an ExpList2D object. 
Represents a vertex in the mesh. 
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
Class for operating on FLD files. 
bool operator==(const Vertex &v)
void Write(const std::string &outFile, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata, const FieldMetaDataMap &fieldinfomap=NullFieldMetaDataMap)
Write data in FLD format. 
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
Hash function for the Vertex struct used for defining sets. 
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
bool operator==(const VertexSharedPtr &v1, const VertexSharedPtr &v2)
Define comparison operator for the vertex struct. 
boost::shared_ptr< PointGeom > PointGeomSharedPtr