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