Nektar++
Fld2Tecplot.cpp
Go to the documentation of this file.
1 #include <cstdio>
2 #include <cstdlib>
3 
5 #include <SolverUtils/Driver.h>
6 
8 
9 using namespace std;
10 using namespace Nektar;
11 using namespace Nektar::SolverUtils;
12 
13 static std::string SetToOneD =
14  LibUtilities::SessionReader::RegisterCmdLineArgument(
15  "SetToOneSpaceDimension", "1", "Redefine mesh to be aligned to x-axis");
16 
17 int main(int argc, char *argv[])
18 {
19  if ((argc < 3) || (argc > 4))
20  {
21  fprintf(stderr, "Usage: ./Fld2Tecplot [-c] file.xml file.fld\n");
22  exit(1);
23  }
24 
27  string vDriverModule;
28  DriverSharedPtr drv;
29 
30  try
31  {
32 
33  // Define new input with extra argument to intialisae -OneD=false
34  int newargc = argc + 1;
35  char **newargv = new char *[newargc];
36 
37  newargv[0] = argv[0];
38  newargv[1] = new char[31];
39  strcpy(newargv[1], "--SetToOneSpaceDimension=false");
40 
41  for (int i = 1; i < argc; ++i)
42  {
43  newargv[i + 1] = argv[i];
44  }
45 
46  // Create session reader and MeshGraph.
47  session = LibUtilities::SessionReader::CreateInstance(newargc, newargv);
48  graph = SpatialDomains::MeshGraph::Read(session);
49  delete[] newargv;
50 
51  // Create driver
52  session->LoadSolverInfo("Driver", vDriverModule, "Standard");
53  drv = GetDriverFactory().CreateInstance(vDriverModule, session, graph);
54 
55  EquationSystemSharedPtr EqSys = drv->GetEqu()[0];
56 
57  PulseWaveSystemSharedPtr PulseWave;
58  if (!(PulseWave = std::dynamic_pointer_cast<PulseWaveSystem>(EqSys)))
59  {
60  ASSERTL0(false,
61  "Failed to dynamically cast to PulseWaveSystemOutput");
62  }
63 
64  std::string fname(argv[argc - 1]);
66 
67  int ndomains = PulseWave->GetNdomains();
68 
69  PulseWave->ImportFldToMultiDomains(
70  fname, Vessels = PulseWave->UpdateVessels(), ndomains);
71  int fdot = fname.find_last_of('.');
72 
73  if (fdot != std::string::npos)
74  {
75  string ending = fname.substr(fdot);
76 
77  // If .chk or .fld we exchange the extension in the output file.
78  // For all other files (e.g. .bse) we append the extension to avoid
79  // conflicts.
80  if (ending == ".chk" || ending == ".fld")
81  {
82  fname = fname.substr(0, fdot);
83  }
84  }
85 
86  fname += ".dat";
87 
88  ofstream outfile(fname.c_str());
89  int nvariables = session->GetVariables().size();
90  std::string var = "";
91  int j;
92  for (j = 0; j < nvariables - 1; ++j)
93  {
94  var += session->GetVariable(j) + ", ";
95  }
96  var += session->GetVariable(j);
97 
98  Vessels[0]->WriteTecplotHeader(outfile, var);
99 
100  for (int n = 0; n < ndomains; ++n)
101  {
102  Vessels[n * nvariables]->WriteTecplotZone(outfile);
103  for (int j = 0; j < nvariables; ++j)
104  {
105  Vessels[n * nvariables + j]->WriteTecplotField(outfile);
106  }
107 
108  Vessels[n * nvariables]->WriteTecplotConnectivity(outfile);
109  }
110  }
111 
112  catch (const std::runtime_error &)
113  {
114  return 1;
115  }
116  catch (const std::string &eStr)
117  {
118  cout << "Error: " << eStr << endl;
119  }
120 
121  return 0;
122 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
int main(int argc, char *argv[])
Definition: Fld2Tecplot.cpp:17
static std::string SetToOneD
Definition: Fld2Tecplot.cpp:13
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Driver > DriverSharedPtr
A shared pointer to a Driver object.
Definition: Driver.h:51
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:64
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< PulseWaveSystem > PulseWaveSystemSharedPtr