Nektar++
Functions
Aliasing.cpp File Reference
#include <cstdio>
#include <cstdlib>
#include <LibUtilities/BasicUtils/SessionReader.h>
#include <SolverUtils/Driver.h>
#include <SpatialDomains/MeshGraphIO.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 49 of file Aliasing.cpp.

50{
51 if (argc != 2)
52 {
53 fprintf(stderr, "Usage: ./Aliasing file.xml \n");
54 fprintf(stderr, "\t Method will read intiial conditions section of "
55 ".xml file for input \n");
56 exit(1);
57 }
58
61 string vDriverModule;
63 try
64 {
65 // Create session reader.
66 session = LibUtilities::SessionReader::CreateInstance(argc, argv);
67
68 // Create MeshGraph.
69 graph = SpatialDomains::MeshGraphIO::Read(session);
70
71 // Create driver
72 session->LoadSolverInfo("Driver", vDriverModule, "Standard");
73 drv = GetDriverFactory().CreateInstance(vDriverModule, session, graph);
74
75 EquationSystemSharedPtr EqSys = drv->GetEqu()[0];
76 IncNavierStokesSharedPtr IncNav = EqSys->as<IncNavierStokes>();
77
78 IncNav->SetInitialConditions(0.0, false);
80 IncNav->UpdateFields();
81
82 int i;
83 int nConvectiveFields = IncNav->GetNConvectiveFields();
84 int nphys = fields[0]->GetTotPoints();
85 Array<OneD, Array<OneD, NekDouble>> VelFields(nConvectiveFields);
86 Array<OneD, Array<OneD, NekDouble>> NonLinear(nConvectiveFields);
87 Array<OneD, Array<OneD, NekDouble>> NonLinearDealiased(
88 nConvectiveFields);
89
90 for (i = 0; i < nConvectiveFields; ++i)
91 {
92 VelFields[i] = fields[i]->UpdatePhys();
93 NonLinear[i] = Array<OneD, NekDouble>(nphys);
94 NonLinearDealiased[i] = Array<OneD, NekDouble>(nphys);
95 }
96
97 std::shared_ptr<NavierStokesAdvection> A =
98 std::dynamic_pointer_cast<NavierStokesAdvection>(
99 IncNav->GetAdvObject());
100
101 if (!A)
102 {
103 cout << "Must use non-linear Navier-Stokes advection" << endl;
104 exit(-1);
105 }
106
107 // calculate non-linear terms without dealiasing
108 A->SetSpecHPDealiasing(false);
109 A->Advect(nConvectiveFields, fields, VelFields, VelFields, NonLinear,
110 0.0);
111
112 // calculate non-linear terms with dealiasing
113 A->SetSpecHPDealiasing(true);
114 A->Advect(nConvectiveFields, fields, VelFields, VelFields,
115 NonLinearDealiased, 0.0);
116
117 // Evaulate Difference and put into fields;
118 for (i = 0; i < nConvectiveFields; ++i)
119 {
120 Vmath::Vsub(nphys, NonLinearDealiased[i], 1, NonLinear[i], 1,
121 NonLinear[i], 1);
122 fields[i]->FwdTransLocalElmt(NonLinear[i],
123 fields[i]->UpdateCoeffs());
124 // Need to reset varibale name for output
125 string name = "NL_Aliasing_" + session->GetVariable(i);
126 session->SetVariable(i, name.c_str());
127 }
128
129 // Reset session name for output file
130 std::string outname = IncNav->GetSessionName();
131
132 outname += "_NonLinear_Aliasing";
133 IncNav->ResetSessionName(outname);
134 IncNav->Output();
135 }
136 catch (const std::runtime_error &)
137 {
138 return 1;
139 }
140 catch (const std::string &eStr)
141 {
142 cout << "Error: " << eStr << endl;
143 }
144
145 return 0;
146}
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.
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:52
DriverFactory & GetDriverFactory()
Definition: Driver.cpp:66
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.hpp:220

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