Nektar++
NonLinearEnergy.cpp
Go to the documentation of this file.
1 #include <cstdio>
2 #include <cstdlib>
3 
5 #include <SolverUtils/Driver.h>
6 
9 
10 using namespace std;
11 using namespace Nektar;
12 using namespace Nektar::SolverUtils;
13 
14 int main(int argc, char *argv[])
15 {
16  if (argc != 2)
17  {
18  fprintf(stderr, "Usage: ./NonLinearEnergy file.xml \n");
19  fprintf(stderr, "\t Method will read intiial conditions section of "
20  ".xml file for input \n");
21  exit(1);
22  }
23 
26  string vDriverModule;
27  DriverSharedPtr drv;
28  try
29  {
30  // Create session reader.
31  session = LibUtilities::SessionReader::CreateInstance(argc, argv);
32 
33  // Create MeshGraph.
34  graph = SpatialDomains::MeshGraph::Read(session);
35 
36  // Create driver
37  session->LoadSolverInfo("Driver", vDriverModule, "Standard");
38  drv = GetDriverFactory().CreateInstance(vDriverModule, session, graph);
39 
40  EquationSystemSharedPtr EqSys = drv->GetEqu()[0];
41  IncNavierStokesSharedPtr IncNav = EqSys->as<IncNavierStokes>();
42 
43  IncNav->SetInitialConditions(0.0, false);
45  IncNav->UpdateFields();
46 
47  int i;
48  int nConvectiveFields = IncNav->GetNConvectiveFields();
49  int nphys = fields[0]->GetTotPoints();
50  Array<OneD, Array<OneD, NekDouble>> VelFields(nConvectiveFields);
51  Array<OneD, Array<OneD, NekDouble>> NonLinear(nConvectiveFields);
52 
53  for (i = 0; i < nConvectiveFields; ++i)
54  {
55  VelFields[i] = fields[i]->UpdatePhys();
56  NonLinear[i] = Array<OneD, NekDouble>(nphys);
57  }
58 
59  std::shared_ptr<NavierStokesAdvection> A =
60  std::dynamic_pointer_cast<NavierStokesAdvection>(
61  IncNav->GetAdvObject());
62 
63  if (!A)
64  {
65  cout << "Must use non-linear Navier-Stokes advection" << endl;
66  exit(-1);
67  }
68 
69  // calculate non-linear terms
70  A->Advect(nConvectiveFields, fields, VelFields, VelFields, NonLinear,
71  0.0);
72 
73  // Evaulate Difference and put into fields;
74  for (i = 0; i < nConvectiveFields; ++i)
75  {
76  fields[i]->FwdTransLocalElmt(NonLinear[i],
77  fields[i]->UpdateCoeffs());
78 
79  // subtract off all modes but top from orthogonal projection
80  for (int n = 0; n < fields[i]->GetExpSize(); ++n)
81  {
82  int offset = fields[i]->GetCoeff_Offset(n);
83  int ncoeffs = fields[i]->GetExp(n)->GetNcoeffs();
84  Array<OneD, NekDouble> coeffs(ncoeffs), coeffsred(ncoeffs), tmp;
85 
86  fields[i]->GetExp(n)->ReduceOrderCoeffs(
87  fields[i]->GetExp(n)->GetBasisNumModes(0) - 1,
88  fields[i]->GetCoeffs() + offset, coeffsred);
89 
90  Vmath::Vsub(ncoeffs, fields[i]->GetCoeffs() + offset, 1,
91  coeffsred, 1,
92  tmp = fields[i]->UpdateCoeffs() + offset, 1);
93  }
94 
95  // Need to reset varibale name for output
96  string name = "NL_TopMode_" + session->GetVariable(i);
97  session->SetVariable(i, name.c_str());
98  }
99 
100  // Reset session name for output file
101  std::string outname = IncNav->GetSessionName();
102 
103  outname += "_NonLinear_Energy";
104  IncNav->ResetSessionName(outname);
105  IncNav->Output();
106  }
107  catch (const std::runtime_error &)
108  {
109  return 1;
110  }
111  catch (const std::string &eStr)
112  {
113  cout << "Error: " << eStr << endl;
114  }
115 
116  return 0;
117 }
int main(int argc, char *argv[])
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:144
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:51
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:64
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
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.cpp:419