Nektar++
Functions
CFLStep.cpp File Reference
#include <cstdio>
#include <cstdlib>
#include <LibUtilities/BasicUtils/SessionReader.h>
#include <SolverUtils/Driver.h>
#include <SpatialDomains/MeshGraphIO.h>
#include <IncNavierStokesSolver/EquationSystems/IncNavierStokes.h>

Go to the source code of this file.

Functions

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

Function Documentation

◆ main()

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

Definition at line 48 of file CFLStep.cpp.

49{
50 if (argc != 2)
51 {
52 fprintf(stderr, "Usage: ./CflStep file.xml \n");
53 fprintf(stderr, "\t Method will read initial conditions section of "
54 ".xml file for input \n");
55 exit(1);
56 }
57
60
61 string vDriverModule;
63 try
64 {
65 // Create session reader.
66 session = LibUtilities::SessionReader::CreateInstance(argc, argv);
67
68 // Create MeshGraph.
69 graph = SpatialDomains::MeshGraphIO::Read(session);
70
71 // Create driver
72 session->LoadSolverInfo("Driver", vDriverModule, "Standard");
73 drv = GetDriverFactory().CreateInstance(vDriverModule, session, graph);
74
75 EquationSystemSharedPtr EqSys = drv->GetEqu()[0];
76 IncNavierStokesSharedPtr IncNav = EqSys->as<IncNavierStokes>();
77
78 IncNav->SetInitialConditions(0.0, false);
79 Array<OneD, NekDouble> cfl = IncNav->GetElmtCFLVals();
80
81 // Reset Pressure field with CFL values
83 IncNav->UpdateFields();
84 int i, n, nquad, cnt;
85 int nfields = fields.size();
86 int nexp = fields[0]->GetExpSize();
87
88 int elmtid = Vmath::Imax(nexp, cfl, 1);
89
90 cout << "Max CFL: " << cfl[elmtid] << " In element " << elmtid << endl;
91
92 for (n = 0; n < nfields; ++n)
93 {
94 if (session->GetVariable(n) == "p")
95 {
96 break;
97 }
98 }
99
100 ASSERTL0(n != nfields, "Could not find field named p in m_fields");
101
102 Array<OneD, NekDouble> phys = fields[n]->UpdatePhys();
103
104 cnt = 0;
105 for (i = 0; i < fields[n]->GetExpSize(); ++i)
106 {
107 nquad = fields[n]->GetExp(i)->GetTotPoints();
108 Vmath::Fill(nquad, cfl[i], &phys[cnt], 1);
109 cnt += nquad;
110 }
111
112 fields[n]->FwdTransLocalElmt(fields[n]->GetPhys(),
113 fields[n]->UpdateCoeffs());
114
115 // Need to reset varibale name for output
116 session->SetVariable(n, "CFL");
117
118 // Reset session name for output file
119 std::string outname = IncNav->GetSessionName();
120
121 outname += "_CFLStep";
122 IncNav->ResetSessionName(outname);
123 IncNav->Output();
124 }
125 catch (const std::runtime_error &)
126 {
127 return 1;
128 }
129 catch (const std::string &eStr)
130 {
131 cout << "Error: " << eStr << endl;
132 }
133
134 return 0;
135}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
This class is the base class for Navier Stokes problems.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
SOLVER_UTILS_EXPORT void SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Initialise the data in the dependent fields.
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:66
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< IncNavierStokes > IncNavierStokesSharedPtr
int Imax(int n, const T *x, const int incx)
Return the index of the maximum element in x.
Definition: Vmath.hpp:623
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Vmath::Fill(), Nektar::SolverUtils::GetDriverFactory(), Vmath::Imax(), and Nektar::SolverUtils::EquationSystem::SetInitialConditions().