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 13 of file NonLinearEnergy.cpp.

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

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