Nektar++
Functions
Aliasing.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 Aliasing.cpp.

49 {
50  if (argc != 2)
51  {
52  fprintf(stderr, "Usage: ./Aliasing 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;
61  DriverSharedPtr drv;
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  Array<OneD, Array<OneD, NekDouble>> NonLinearDealiased(
87  nConvectiveFields);
88 
89  for (i = 0; i < nConvectiveFields; ++i)
90  {
91  VelFields[i] = fields[i]->UpdatePhys();
92  NonLinear[i] = Array<OneD, NekDouble>(nphys);
93  NonLinearDealiased[i] = Array<OneD, NekDouble>(nphys);
94  }
95 
96  std::shared_ptr<NavierStokesAdvection> A =
97  std::dynamic_pointer_cast<NavierStokesAdvection>(
98  IncNav->GetAdvObject());
99 
100  if (!A)
101  {
102  cout << "Must use non-linear Navier-Stokes advection" << endl;
103  exit(-1);
104  }
105 
106  // calculate non-linear terms without dealiasing
107  A->SetSpecHPDealiasing(false);
108  A->Advect(nConvectiveFields, fields, VelFields, VelFields, NonLinear,
109  0.0);
110 
111  // calculate non-linear terms with dealiasing
112  A->SetSpecHPDealiasing(true);
113  A->Advect(nConvectiveFields, fields, VelFields, VelFields,
114  NonLinearDealiased, 0.0);
115 
116  // Evaulate Difference and put into fields;
117  for (i = 0; i < nConvectiveFields; ++i)
118  {
119  Vmath::Vsub(nphys, NonLinearDealiased[i], 1, NonLinear[i], 1,
120  NonLinear[i], 1);
121  fields[i]->FwdTransLocalElmt(NonLinear[i],
122  fields[i]->UpdateCoeffs());
123  // Need to reset varibale name for output
124  string name = "NL_Aliasing_" + session->GetVariable(i);
125  session->SetVariable(i, name.c_str());
126  }
127 
128  // Reset session name for output file
129  std::string outname = IncNav->GetSessionName();
130 
131  outname += "_NonLinear_Aliasing";
132  IncNav->ResetSessionName(outname);
133  IncNav->Output();
134  }
135  catch (const std::runtime_error &)
136  {
137  return 1;
138  }
139  catch (const std::string &eStr)
140  {
141  cout << "Error: " << eStr << endl;
142  }
143 
144  return 0;
145 }
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
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

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