Nektar++
Functions | Variables
Fld2Tecplot.cpp File Reference
#include <cstdio>
#include <cstdlib>
#include <LibUtilities/BasicUtils/SessionReader.h>
#include <SolverUtils/Driver.h>
#include <PulseWaveSolver/EquationSystems/PulseWaveSystem.h>

Go to the source code of this file.

Functions

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

Variables

static std::string SetToOneD
 

Function Documentation

◆ main()

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

Definition at line 51 of file Fld2Tecplot.cpp.

52{
53 if ((argc < 3) || (argc > 4))
54 {
55 fprintf(stderr, "Usage: ./Fld2Tecplot [-c] file.xml file.fld\n");
56 exit(1);
57 }
58
61 string vDriverModule;
63
64 try
65 {
66
67 // Define new input with extra argument to intialisae -OneD=false
68 int newargc = argc + 1;
69 char **newargv = new char *[newargc];
70
71 newargv[0] = argv[0];
72 newargv[1] = new char[31];
73 strcpy(newargv[1], "--SetToOneSpaceDimension=false");
74
75 for (int i = 1; i < argc; ++i)
76 {
77 newargv[i + 1] = argv[i];
78 }
79
80 // Create session reader and MeshGraph.
81 session = LibUtilities::SessionReader::CreateInstance(newargc, newargv);
82 graph = SpatialDomains::MeshGraph::Read(session);
83 delete[] newargv;
84
85 // Create driver
86 session->LoadSolverInfo("Driver", vDriverModule, "Standard");
87 drv = GetDriverFactory().CreateInstance(vDriverModule, session, graph);
88
89 EquationSystemSharedPtr EqSys = drv->GetEqu()[0];
90
92 if (!(PulseWave = std::dynamic_pointer_cast<PulseWaveSystem>(EqSys)))
93 {
94 ASSERTL0(false,
95 "Failed to dynamically cast to PulseWaveSystemOutput");
96 }
97
98 std::string fname(argv[argc - 1]);
100
101 int ndomains = PulseWave->GetNdomains();
102
103 PulseWave->ImportFldToMultiDomains(
104 fname, Vessels = PulseWave->UpdateVessels(), ndomains);
105 int fdot = fname.find_last_of('.');
106
107 if (fdot != std::string::npos)
108 {
109 string ending = fname.substr(fdot);
110
111 // If .chk or .fld we exchange the extension in the output file.
112 // For all other files (e.g. .bse) we append the extension to avoid
113 // conflicts.
114 if (ending == ".chk" || ending == ".fld")
115 {
116 fname = fname.substr(0, fdot);
117 }
118 }
119
120 fname += ".dat";
121
122 ofstream outfile(fname.c_str());
123 int nvariables = session->GetVariables().size();
124 std::string var = "";
125 int j;
126 for (j = 0; j < nvariables - 1; ++j)
127 {
128 var += session->GetVariable(j) + ", ";
129 }
130 var += session->GetVariable(j);
131
132 Vessels[0]->WriteTecplotHeader(outfile, var);
133
134 for (int n = 0; n < ndomains; ++n)
135 {
136 Vessels[n * nvariables]->WriteTecplotZone(outfile);
137 for (int j = 0; j < nvariables; ++j)
138 {
139 Vessels[n * nvariables + j]->WriteTecplotField(outfile);
140 }
141
142 Vessels[n * nvariables]->WriteTecplotConnectivity(outfile);
143 }
144 }
145
146 catch (const std::runtime_error &)
147 {
148 return 1;
149 }
150 catch (const std::string &eStr)
151 {
152 cout << "Error: " << eStr << endl;
153 }
154
155 return 0;
156}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Driver > DriverSharedPtr
A shared pointer to a Driver object.
Definition: Driver.h:52
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
std::shared_ptr< PulseWaveSystem > PulseWaveSystemSharedPtr

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), and Nektar::SolverUtils::GetDriverFactory().

Variable Documentation

◆ SetToOneD

std::string SetToOneD
static
Initial value:
=
LibUtilities::SessionReader::RegisterCmdLineArgument(
"SetToOneSpaceDimension", "1", "Redefine mesh to be aligned to x-axis")

Definition at line 47 of file Fld2Tecplot.cpp.