Nektar++
NonLinearEnergy.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: NonLinearEnergy.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description:
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <cstdio>
36#include <cstdlib>
37
39#include <SolverUtils/Driver.h>
41
44
45using namespace std;
46using namespace Nektar;
47using namespace Nektar::SolverUtils;
48
49int main(int argc, char *argv[])
50{
51 if (argc != 2)
52 {
53 fprintf(stderr, "Usage: ./NonLinearEnergy 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
88 for (i = 0; i < nConvectiveFields; ++i)
89 {
90 VelFields[i] = fields[i]->UpdatePhys();
91 NonLinear[i] = Array<OneD, NekDouble>(nphys);
92 }
93
94 std::shared_ptr<NavierStokesAdvection> A =
95 std::dynamic_pointer_cast<NavierStokesAdvection>(
96 IncNav->GetAdvObject());
97
98 if (!A)
99 {
100 cout << "Must use non-linear Navier-Stokes advection" << endl;
101 exit(-1);
102 }
103
104 // calculate non-linear terms
105 A->Advect(nConvectiveFields, fields, VelFields, VelFields, NonLinear,
106 0.0);
107
108 // Evaulate Difference and put into fields;
109 for (i = 0; i < nConvectiveFields; ++i)
110 {
111 fields[i]->FwdTransLocalElmt(NonLinear[i],
112 fields[i]->UpdateCoeffs());
113
114 // subtract off all modes but top from orthogonal projection
115 for (int n = 0; n < fields[i]->GetExpSize(); ++n)
116 {
117 int offset = fields[i]->GetCoeff_Offset(n);
118 int ncoeffs = fields[i]->GetExp(n)->GetNcoeffs();
119 Array<OneD, NekDouble> coeffs(ncoeffs), coeffsred(ncoeffs), tmp;
120
121 fields[i]->GetExp(n)->ReduceOrderCoeffs(
122 fields[i]->GetExp(n)->GetBasisNumModes(0) - 1,
123 fields[i]->GetCoeffs() + offset, coeffsred);
124
125 Vmath::Vsub(ncoeffs, fields[i]->GetCoeffs() + offset, 1,
126 coeffsred, 1,
127 tmp = fields[i]->UpdateCoeffs() + offset, 1);
128 }
129
130 // Need to reset varibale name for output
131 string name = "NL_TopMode_" + session->GetVariable(i);
132 session->SetVariable(i, name.c_str());
133 }
134
135 // Reset session name for output file
136 std::string outname = IncNav->GetSessionName();
137
138 outname += "_NonLinear_Energy";
139 IncNav->ResetSessionName(outname);
140 IncNav->Output();
141 }
142 catch (const std::runtime_error &)
143 {
144 return 1;
145 }
146 catch (const std::string &eStr)
147 {
148 cout << "Error: " << eStr << endl;
149 }
150
151 return 0;
152}
int main(int argc, char *argv[])
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
STL namespace.