Nektar++
Functions
NonLinearEnergy.cpp File Reference
#include <cstdio>
#include <cstdlib>
#include <LibUtilities/BasicUtils/SessionReader.h>
#include <SolverUtils/Driver.h>
#include <IncNavierStokesSolver/AdvectionTerms/NavierStokesAdvection.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 NonLinearEnergy.cpp.

49{
50 if (argc != 2)
51 {
52 fprintf(stderr, "Usage: ./NonLinearEnergy file.xml \n");
53 fprintf(stderr, "\t Method will read intiial conditions section of "
54 ".xml file for input \n");
55 exit(1);
56 }
57
60 string vDriverModule;
62 try
63 {
64 // Create session reader.
65 session = LibUtilities::SessionReader::CreateInstance(argc, argv);
66
67 // Create MeshGraph.
68 graph = SpatialDomains::MeshGraph::Read(session);
69
70 // Create driver
71 session->LoadSolverInfo("Driver", vDriverModule, "Standard");
72 drv = GetDriverFactory().CreateInstance(vDriverModule, session, graph);
73
74 EquationSystemSharedPtr EqSys = drv->GetEqu()[0];
75 IncNavierStokesSharedPtr IncNav = EqSys->as<IncNavierStokes>();
76
77 IncNav->SetInitialConditions(0.0, false);
79 IncNav->UpdateFields();
80
81 int i;
82 int nConvectiveFields = IncNav->GetNConvectiveFields();
83 int nphys = fields[0]->GetTotPoints();
84 Array<OneD, Array<OneD, NekDouble>> VelFields(nConvectiveFields);
85 Array<OneD, Array<OneD, NekDouble>> NonLinear(nConvectiveFields);
86
87 for (i = 0; i < nConvectiveFields; ++i)
88 {
89 VelFields[i] = fields[i]->UpdatePhys();
90 NonLinear[i] = Array<OneD, NekDouble>(nphys);
91 }
92
93 std::shared_ptr<NavierStokesAdvection> A =
94 std::dynamic_pointer_cast<NavierStokesAdvection>(
95 IncNav->GetAdvObject());
96
97 if (!A)
98 {
99 cout << "Must use non-linear Navier-Stokes advection" << endl;
100 exit(-1);
101 }
102
103 // calculate non-linear terms
104 A->Advect(nConvectiveFields, fields, VelFields, VelFields, NonLinear,
105 0.0);
106
107 // Evaulate Difference and put into fields;
108 for (i = 0; i < nConvectiveFields; ++i)
109 {
110 fields[i]->FwdTransLocalElmt(NonLinear[i],
111 fields[i]->UpdateCoeffs());
112
113 // subtract off all modes but top from orthogonal projection
114 for (int n = 0; n < fields[i]->GetExpSize(); ++n)
115 {
116 int offset = fields[i]->GetCoeff_Offset(n);
117 int ncoeffs = fields[i]->GetExp(n)->GetNcoeffs();
118 Array<OneD, NekDouble> coeffs(ncoeffs), coeffsred(ncoeffs), tmp;
119
120 fields[i]->GetExp(n)->ReduceOrderCoeffs(
121 fields[i]->GetExp(n)->GetBasisNumModes(0) - 1,
122 fields[i]->GetCoeffs() + offset, coeffsred);
123
124 Vmath::Vsub(ncoeffs, fields[i]->GetCoeffs() + offset, 1,
125 coeffsred, 1,
126 tmp = fields[i]->UpdateCoeffs() + offset, 1);
127 }
128
129 // Need to reset varibale name for output
130 string name = "NL_TopMode_" + session->GetVariable(i);
131 session->SetVariable(i, name.c_str());
132 }
133
134 // Reset session name for output file
135 std::string outname = IncNav->GetSessionName();
136
137 outname += "_NonLinear_Energy";
138 IncNav->ResetSessionName(outname);
139 IncNav->Output();
140 }
141 catch (const std::runtime_error &)
142 {
143 return 1;
144 }
145 catch (const std::string &eStr)
146 {
147 cout << "Error: " << eStr << endl;
148 }
149
150 return 0;
151}
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.
Definition: NekFactory.hpp:143
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: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< IncNavierStokes > IncNavierStokesSharedPtr
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::SolverUtils::GetDriverFactory(), CellMLToNektar.pycml::name, Nektar::SolverUtils::EquationSystem::SetInitialConditions(), and Vmath::Vsub().