Nektar++
Classes | Typedefs | Functions
VtkToFld.cpp File Reference
#include <LibUtilities/BasicUtils/VtkUtil.hpp>
#include <vtkPolyDataReader.h>
#include <vtkCellArray.h>
#include <vtkCellDataToPointData.h>
#include <vtkContourFilter.h>
#include <vtkPointData.h>
#include <vtkPoints.h>
#include <vtkPolyData.h>
#include <boost/program_options.hpp>
#include <LibUtilities/BasicUtils/FieldIO.h>
#include <LibUtilities/BasicUtils/HashUtils.hpp>
#include <LibUtilities/BasicUtils/SessionReader.h>
#include <LibUtilities/Communication/Comm.h>
#include <LocalRegions/Expansion2D.h>
#include <MultiRegions/ExpList.h>
#include <SpatialDomains/MeshGraphIO.h>

Go to the source code of this file.

Classes

struct  Vertex
 Represents a vertex in the mesh. More...
 
struct  VertexHash
 Hash function for the Vertex struct used for defining sets. More...
 

Typedefs

typedef std::shared_ptr< VertexVertexSharedPtr
 
typedef std::unordered_set< VertexSharedPtr, VertexHashVertexSet
 Define a set of Vertices with associated hash function. More...
 

Functions

bool operator== (const VertexSharedPtr &v1, const VertexSharedPtr &v2)
 Define comparison operator for the vertex struct. More...
 
int main (int argc, char *argv[])
 

Typedef Documentation

◆ VertexSet

typedef std::unordered_set<VertexSharedPtr, VertexHash> VertexSet

Define a set of Vertices with associated hash function.

Definition at line 112 of file VtkToFld.cpp.

◆ VertexSharedPtr

typedef std::shared_ptr<Vertex> VertexSharedPtr

Definition at line 92 of file VtkToFld.cpp.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Main function.

Usage: VtkToFld session.xml input.vtk output.fld [options]

Definition at line 119 of file VtkToFld.cpp.

120{
121 // Set up available options
122 po::options_description desc("Available options");
123
124 // clang-format off
125 desc.add_options()
126 ("help,h", "Produce this help message.")
127 ("name,n", po::value<string>()->default_value("Intensity"),
128 "Name of field in VTK file to use for intensity.")
129 ("outname,m", po::value<string>()->default_value("intensity"),
130 "Name of field in output FLD file.")
131 ("precision,p", po::value<double>()->default_value(1),
132 "Precision of vertex matching.");
133 // clang-format on
134
135 po::options_description hidden("Hidden options");
136
137 // clang-format off
138 hidden.add_options()
139 ("file", po::value<vector<string> >(), "Input filename");
140 // clang-format on
141
142 po::options_description cmdline_options;
143 cmdline_options.add(desc).add(hidden);
144
145 po::positional_options_description p;
146 p.add("file", -1);
147
148 po::variables_map vm;
149
150 // Parse command-line options
151 try
152 {
153 po::store(po::command_line_parser(argc, argv)
154 .options(cmdline_options)
155 .positional(p)
156 .run(),
157 vm);
158 po::notify(vm);
159 }
160 catch (const std::exception &e)
161 {
162 cerr << e.what() << endl;
163 cerr << desc;
164 return 1;
165 }
166
167 if (vm.count("help") || vm.count("file") == 0 ||
168 vm["file"].as<vector<string>>().size() != 3)
169 {
170 cerr << "Usage: VtkToFld session.xml intensity.vtk output.fld [options]"
171 << endl;
172 cerr << desc;
173 return 1;
174 }
175
176 // Extract command-line argument values
177 std::vector<std::string> vFiles = vm["file"].as<vector<string>>();
178 const string infile = vFiles[1];
179 const string outfile = vFiles[2];
180 const double factor = vm["precision"].as<double>();
181 const string name = vm["name"].as<string>();
182 const string outname = vm["outname"].as<string>();
183
184 std::vector<std::string> vFilenames;
188
189 vFilenames.push_back(vFiles[0]);
190 vSession = LibUtilities::SessionReader::CreateInstance(2, argv, vFilenames);
191
192 try
193 {
194 //----------------------------------------------
195 // Read in mesh from input file
196 graph2D = SpatialDomains::MeshGraphIO::Read(vSession);
197 //----------------------------------------------
198
199 //----------------------------------------------
200 // Define Expansion
202 graph2D);
203 //----------------------------------------------
204
205 //----------------------------------------------
206 // Set up coordinates of mesh
207 int coordim = Exp->GetCoordim(0);
208 int nq = Exp->GetNpoints();
209
210 Array<OneD, NekDouble> xc0(nq, 0.0);
211 Array<OneD, NekDouble> xc1(nq, 0.0);
212 Array<OneD, NekDouble> xc2(nq, 0.0);
213
214 switch (coordim)
215 {
216 case 2:
217 Exp->GetCoords(xc0, xc1);
218 break;
219 case 3:
220 Exp->GetCoords(xc0, xc1, xc2);
221 break;
222 default:
223 ASSERTL0(false, "Coordim not valid");
224 break;
225 }
226 //----------------------------------------------
227
228 vtkPolyDataReader *vtkMeshReader = vtkPolyDataReader::New();
229 vtkMeshReader->SetFileName(infile.c_str());
230 vtkMeshReader->Update();
231
232 vtkPolyData *vtkMesh = vtkMeshReader->GetOutput();
233 vtkCellDataToPointData *c2p = vtkCellDataToPointData::New();
234#if VTK_MAJOR_VERSION <= 5
235 c2p->SetInput(vtkMesh);
236#else
237 c2p->SetInputData(vtkMesh);
238#endif
239 c2p->PassCellDataOn();
240 c2p->Update();
241 vtkPolyData *vtkDataAtPoints = c2p->GetPolyDataOutput();
242
243 vtkPoints *vtkPoints = vtkMesh->GetPoints();
244 ASSERTL0(vtkPoints, "ERROR: cannot get points from mesh.");
245
246 vtkCellArray *vtkPolys = vtkMesh->GetPolys();
247 ASSERTL0(vtkPolys, "ERROR: cannot get polygons from mesh.");
248
249 vtkPointData *vtkPData = vtkDataAtPoints->GetPointData();
250 ASSERTL0(vtkPolys, "ERROR: cannot get point data from file.");
251
252 VertexSet points;
253 VertexSet::iterator vIter;
254 double p[3];
255 double val;
256 double x, y, z;
257 int coeff_idx;
258 int i, j, n;
259
260 if (!vtkDataAtPoints->GetPointData()->HasArray(name.c_str()))
261 {
262 n = vtkDataAtPoints->GetPointData()->GetNumberOfArrays();
263 cerr << "Input file '" << infile
264 << "' does not have a field named '" << name << "'" << endl;
265 cerr << "There are " << n << " arrays in this file." << endl;
266 for (int i = 0; i < n; ++i)
267 {
268 cerr << " "
269 << vtkDataAtPoints->GetPointData()->GetArray(i)->GetName()
270 << endl;
271 }
272 return 1;
273 }
274
275 // Build up an unordered set of vertices from the VTK file. For each
276 // vertex a hashed value of the coordinates is generated to within a
277 // given tolerance.
278 n = vtkPoints->GetNumberOfPoints();
279 for (i = 0; i < n; ++i)
280 {
281 vtkPoints->GetPoint(i, p);
282 val = vtkPData->GetScalars(name.c_str())->GetTuple1(i);
283 std::shared_ptr<Vertex> v(
284 new Vertex(p[0], p[1], p[2], val, factor));
285 points.insert(v);
286 }
287
288 // Now process each vertex of each element in the mesh
290 for (i = 0; i < Exp->GetNumElmts(); ++i)
291 {
292 StdRegions::StdExpansionSharedPtr e = Exp->GetExp(i);
293 for (j = 0; j < e->GetNverts(); ++j)
294 {
295 // Get the index of the coefficient corresponding to this vertex
296 coeff_idx = Exp->GetCoeff_Offset(i) + e->GetVertexMap(j);
297
298 // Get the coordinates of the vertex
299 vert =
300 e->as<LocalRegions::Expansion2D>()->GetGeom2D()->GetVertex(
301 j);
302 vert->GetCoords(x, y, z);
303
304 // Look up the vertex in the VertexSet
305 std::shared_ptr<Vertex> v(new Vertex(x, y, z, 0.0, factor));
306 vIter = points.find(v);
307
308 // If not found, maybe the tolerance should be reduced?
309 // If found, record the scalar value from the VTK file in the
310 // corresponding coefficient.
311 if (vIter == points.end())
312 {
313 cerr << "Vertex " << i << " not found. Looking for (" << x
314 << ", " << y << ", " << z << ")" << endl;
315 }
316 else
317 {
318 Exp->UpdateCoeffs()[coeff_idx] = (*vIter)->scalar;
319 }
320 }
321 }
322 Exp->SetPhysState(false);
323
324 //-----------------------------------------------
325 // Write solution to file
326 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
327 Exp->GetFieldDefinitions();
328 std::vector<std::vector<NekDouble>> FieldData(FieldDef.size());
329
330 for (i = 0; i < FieldDef.size(); ++i)
331 {
332 FieldDef[i]->m_fields.push_back(outname);
333 Exp->AppendFieldData(FieldDef[i], FieldData[i]);
334 }
335
337 LibUtilities::FieldIO::CreateDefault(vSession);
338 vFld->Write(outfile, FieldDef, FieldData);
339 //-----------------------------------------------
340 }
341 catch (...)
342 {
343 cout << "An error occurred." << endl;
344 }
345}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
std::unordered_set< VertexSharedPtr, VertexHash > VertexSet
Define a set of Vertices with associated hash function.
Definition: VtkToFld.cpp:112
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
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
Definition: StdExpansion.h:579
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::vector< double > z(NPUPPER)
Represents a vertex in the mesh.
Definition: VtkToFld.cpp:75

References ASSERTL0, Nektar::StdRegions::StdExpansion::GetCoords(), CellMLToNektar.pycml::name, CellMLToNektar.cellml_metadata::p, CellMLToNektar.translators::run(), and Nektar::UnitTests::z().

◆ operator==()

bool operator== ( const VertexSharedPtr v1,
const VertexSharedPtr v2 
)

Define comparison operator for the vertex struct.

Definition at line 95 of file VtkToFld.cpp.

96{
97 return v1->operator==(*v2);
98}