Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 Nektar;
10 using namespace Nektar::SolverUtils;
11 
13  "CharateristicVariables","c","Output characteristic variables");
14 
15 static std::string SetToOneD = LibUtilities::SessionReader::RegisterCmdLineArgument("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 
25 
27  string vDriverModule;
28  DriverSharedPtr drv;
29 
30 
31  try
32  {
33 
34  // Define new input with extra argument to intialisae -OneD=false
35  int newargc = argc+1;
36  char **newargv = new char*[newargc];
37 
38  newargv[0] = argv[0];
39  newargv[1] = new char[30];
40  strcpy(newargv[1], "--SetToOneSpaceDimension=false");
41 
42  for(int i = 1; i < argc; ++i)
43  {
44  newargv[i+1] = argv[i];
45  }
46 
47  // Create session reader.
48  session = LibUtilities::SessionReader::CreateInstance(newargc, newargv);
49  delete[] newargv;
50 
51  bool CalcCharacteristicVariables = false;
52 
53  if(session->DefinesCmdLineArgument(cvar))
54  {
55  CalcCharacteristicVariables = true;
56  }
57 
58 
59  // Create driver
60  session->LoadSolverInfo("Driver", vDriverModule, "Standard");
61  drv = GetDriverFactory().CreateInstance(vDriverModule, session);
62 
63  EquationSystemSharedPtr EqSys = drv->GetEqu()[0];
64 
65  PulseWaveSystemSharedPtr PulseWave;
66  if(!(PulseWave = boost::dynamic_pointer_cast
67  <PulseWaveSystem>(EqSys)))
68  {
69  ASSERTL0(false,"Failed to dynamically cast to PulseWaveSystemOutput");
70  }
71 
72  std::string fname(argv[argc-1]);
73  Array<OneD, MultiRegions::ExpListSharedPtr> Vessels;
74 
75  int ndomains = PulseWave->GetNdomains();
76 
77  PulseWave->ImportFldToMultiDomains(fname,Vessels = PulseWave->UpdateVessels(),
78  ndomains);
79  int fdot = fname.find_last_of('.');
80 
81  if (fdot != std::string::npos)
82  {
83  string ending = fname.substr(fdot);
84 
85  // If .chk or .fld we exchange the extension in the output file.
86  // For all other files (e.g. .bse) we append the extension to avoid
87  // conflicts.
88  if (ending == ".chk" || ending == ".fld")
89  {
90  fname = fname.substr(0,fdot);
91  }
92  }
93 
94  fname = fname + ".dat";
95 
96  ofstream outfile(fname.c_str());
97  int nvariables = session->GetVariables().size();
98  std::string var = "";
99  int j;
100  for(j = 0; j < nvariables-1; ++j)
101  {
102  var += session->GetVariable(j) + ", ";
103  }
104  var += session->GetVariable(j);
105 
106  if(CalcCharacteristicVariables)
107  {
108  var += ", Char1, Char2";
109  }
110 
111  Vessels[0]->WriteTecplotHeader(outfile,var);
112 
113  for(int n = 0; n < ndomains; ++n)
114  {
115  Vessels[n*nvariables]->WriteTecplotZone(outfile);
116  for(int j = 0; j < nvariables; ++j)
117  {
118  Vessels[n*nvariables+j]->WriteTecplotField(outfile);
119  }
120 
121  if(CalcCharacteristicVariables)
122  {
123  PulseWave->CalcCharacteristicVariables(n*nvariables);
124 
125  for(int j = 0; j < nvariables; ++j)
126  {
127  Vessels[n*nvariables+j]->WriteTecplotField(outfile);
128  }
129  }
130  Vessels[n*nvariables]->WriteTecplotConnectivity(outfile);
131  }
132  }
133 
134  catch (const std::runtime_error&)
135  {
136  return 1;
137  }
138  catch (const std::string& eStr)
139  {
140  cout << "Error: " << eStr << endl;
141  }
142 
143  return 0;
144 }