Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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>
Include dependency graph for NonLinearEnergy.cpp:

Go to the source code of this file.

Functions

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

Function Documentation

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

Definition at line 14 of file NonLinearEnergy.cpp.

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

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 
24  string vDriverModule;
25  DriverSharedPtr drv;
26  try
27  {
28  // Create session reader.
29  session = LibUtilities::SessionReader::CreateInstance(argc, argv);
30 
31  // Create driver
32  session->LoadSolverInfo("Driver", vDriverModule, "Standard");
33  drv = GetDriverFactory().CreateInstance(vDriverModule, session);
34 
35 
36  EquationSystemSharedPtr EqSys = drv->GetEqu()[0];
37  IncNavierStokesSharedPtr IncNav = EqSys->as<IncNavierStokes>();
38 
39  IncNav->SetInitialConditions(0.0,false);
40  Array<OneD, MultiRegions::ExpListSharedPtr> fields = IncNav->UpdateFields();
41 
42  int i;
43  int nConvectiveFields = IncNav->GetNConvectiveFields();
44  int nphys = fields[0]->GetTotPoints();
45  Array<OneD, Array<OneD, NekDouble> > VelFields(nConvectiveFields);
46  Array<OneD, Array<OneD, NekDouble> > NonLinear(nConvectiveFields);
47 
48  for(i = 0; i < nConvectiveFields; ++i)
49  {
50  VelFields[i] = fields[i]->UpdatePhys();
51  NonLinear[i] = Array<OneD, NekDouble> (nphys);
52  }
53 
54  boost::shared_ptr<NavierStokesAdvection> A
55  = boost::dynamic_pointer_cast<NavierStokesAdvection>(IncNav->GetAdvObject());
56 
57  if (!A)
58  {
59  cout << "Must use non-linear Navier-Stokes advection" << endl;
60  exit(-1);
61  }
62 
63  // calculate non-linear terms
64  A->Advect(nConvectiveFields, fields,
65  VelFields, VelFields,
66  NonLinear, 0.0);
67 
68 
69  // Evaulate Difference and put into fields;
70  for(i = 0; i < nConvectiveFields; ++i)
71  {
72  fields[i]->FwdTrans_IterPerExp(NonLinear[i],fields[i]->UpdateCoeffs());
73 
74  // subtract off all modes but top from orthogonal projection
75  int cnt = 0;
76  for(int n = 0; n < fields[i]->GetExpSize(); ++n)
77  {
78  int offset = fields[i]->GetCoeff_Offset(n);
79  int ncoeffs = fields[i]->GetExp(n)->GetNcoeffs();
80  Array<OneD, NekDouble> coeffs(ncoeffs),coeffsred(ncoeffs),tmp;
81 
82  fields[i]->GetExp(n)->ReduceOrderCoeffs(fields[i]->GetExp(n)->GetBasisNumModes(0)-1,fields[i]->GetCoeffs()+offset,coeffsred);
83 
84  Vmath::Vsub(ncoeffs,fields[i]->GetCoeffs()+offset,1,
85  coeffsred,1,tmp = fields[i]->UpdateCoeffs()+
86  offset,1);
87  cnt += ncoeffs;
88  }
89 
90  // Need to reset varibale name for output
91  string name = "NL_TopMode_"+session->GetVariable(i);
92  session->SetVariable(i,name.c_str());
93  }
94 
95 
96  // Reset session name for output file
97  std::string outname = IncNav->GetSessionName();
98 
99  outname += "_NonLinear_Energy";
100  IncNav->ResetSessionName(outname);
101  IncNav->Output();
102 
103  }
104  catch (const std::runtime_error&)
105  {
106  return 1;
107  }
108  catch (const std::string& eStr)
109  {
110  cout << "Error: " << eStr << endl;
111  }
112 
113 
114  return 0;
115 }
boost::shared_ptr< Driver > DriverSharedPtr
A shared pointer to a Driver object.
Definition: Driver.h:52
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
boost::shared_ptr< IncNavierStokes > IncNavierStokesSharedPtr
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
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:329
This class is the base class for Navier Stokes problems.
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:66
SOLVER_UTILS_EXPORT void SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Initialise the data in the dependent fields.