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

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 .xml file for input \n");
20  exit(1);
21  }
22 
25  string vDriverModule;
26  DriverSharedPtr drv;
27  try
28  {
29  // Create session reader.
30  session = LibUtilities::SessionReader::CreateInstance(argc, argv);
31 
32  // Create MeshGraph.
33  graph = SpatialDomains::MeshGraph::Read(session);
34 
35  // Create driver
36  session->LoadSolverInfo("Driver", vDriverModule, "Standard");
37  drv = GetDriverFactory().CreateInstance(vDriverModule, session, graph);
38 
39 
40  EquationSystemSharedPtr EqSys = drv->GetEqu()[0];
41  IncNavierStokesSharedPtr IncNav = EqSys->as<IncNavierStokes>();
42 
43  IncNav->SetInitialConditions(0.0,false);
44  Array<OneD, MultiRegions::ExpListSharedPtr> fields = IncNav->UpdateFields();
45 
46  int i;
47  int nConvectiveFields = IncNav->GetNConvectiveFields();
48  int nphys = fields[0]->GetTotPoints();
49  Array<OneD, Array<OneD, NekDouble> > VelFields(nConvectiveFields);
50  Array<OneD, Array<OneD, NekDouble> > NonLinear(nConvectiveFields);
51 
52  for(i = 0; i < nConvectiveFields; ++i)
53  {
54  VelFields[i] = fields[i]->UpdatePhys();
55  NonLinear[i] = Array<OneD, NekDouble> (nphys);
56  }
57 
58  std::shared_ptr<NavierStokesAdvection> A
59  = std::dynamic_pointer_cast<NavierStokesAdvection>(IncNav->GetAdvObject());
60 
61  if (!A)
62  {
63  cout << "Must use non-linear Navier-Stokes advection" << endl;
64  exit(-1);
65  }
66 
67  // calculate non-linear terms
68  A->Advect(nConvectiveFields, fields,
69  VelFields, VelFields,
70  NonLinear, 0.0);
71 
72 
73  // Evaulate Difference and put into fields;
74  for(i = 0; i < nConvectiveFields; ++i)
75  {
76  fields[i]->FwdTrans_IterPerExp(NonLinear[i],fields[i]->UpdateCoeffs());
77 
78  // subtract off all modes but top from orthogonal projection
79  int cnt = 0;
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(fields[i]->GetExp(n)->GetBasisNumModes(0)-1,fields[i]->GetCoeffs()+offset,coeffsred);
87 
88  Vmath::Vsub(ncoeffs,fields[i]->GetCoeffs()+offset,1,
89  coeffsred,1,tmp = fields[i]->UpdateCoeffs()+
90  offset,1);
91  cnt += ncoeffs;
92  }
93 
94  // Need to reset varibale name for output
95  string name = "NL_TopMode_"+session->GetVariable(i);
96  session->SetVariable(i,name.c_str());
97  }
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  }
108  catch (const std::runtime_error&)
109  {
110  return 1;
111  }
112  catch (const std::string& eStr)
113  {
114  cout << "Error: " << eStr << endl;
115  }
116 
117 
118  return 0;
119 }
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:145
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: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.cpp:372

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