Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
Nektar::EigenValuesAdvection Class Reference

#include <EigenValuesAdvection.h>

Inheritance diagram for Nektar::EigenValuesAdvection:
[legend]

Public Member Functions

 ~EigenValuesAdvection () override=default
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, const Array< OneD, const NekDouble > &input)
 
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > & UpdatePhysField (const int i)
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 

Protected Member Functions

 EigenValuesAdvection (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
Array< OneD, NekDouble > & GetNormalVelocity ()
 Get the normal velocity. More...
 
void v_InitObject (bool DeclareFields=true) override
 Initialisation object for EquationSystem. More...
 
void v_DoInitialise (bool dumpInitialConditions=false) override
 Virtual function for initialisation implementation. More...
 
void v_DoSolve () override
 Virtual function for solve implementation. More...
 
void GetFluxVector (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
Array< OneD, Array< OneD, NekDouble > > m_velocity
 
SolverUtils::AdvectionSharedPtr m_advObject
 
Array< OneD, NekDoublem_traceVn
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
Array< OneD, NekDoublem_movingFrameVelsxyz
 Moving frame of reference velocities (u, v, w, omega_x, omega_y, omega_z, a_x, a_y, a_z, domega_x, domega_y, domega_z) More...
 
Array< OneD, NekDoublem_movingFrameData
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 

Friends

class MemoryManager< EigenValuesAdvection >
 

Additional Inherited Members

- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

Definition at line 46 of file EigenValuesAdvection.h.

Constructor & Destructor Documentation

◆ ~EigenValuesAdvection()

Nektar::EigenValuesAdvection::~EigenValuesAdvection ( )
overridedefault

◆ EigenValuesAdvection()

Nektar::EigenValuesAdvection::EigenValuesAdvection ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
protected

Definition at line 48 of file EigenValuesAdvection.cpp.

51 : EquationSystem(pSession, pGraph)
52{
53}
SOLVER_UTILS_EXPORT EquationSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises EquationSystem class members.

Member Function Documentation

◆ create()

static EquationSystemSharedPtr Nektar::EigenValuesAdvection::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file EigenValuesAdvection.h.

55 {
58 pGraph);
59 p->InitObject();
60 return p;
61 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ GetFluxVector()

void Nektar::EigenValuesAdvection::GetFluxVector ( const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Definition at line 285 of file EigenValuesAdvection.cpp.

288{
289 ASSERTL1(flux[0].size() == m_velocity.size(),
290 "Dimension of flux array and velocity array do not match");
291
292 int nq = physfield[0].size();
293
294 for (int i = 0; i < flux.size(); ++i)
295 {
296 for (int j = 0; j < flux[0].size(); ++j)
297 {
298 Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
299 }
300 }
301}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Array< OneD, Array< OneD, NekDouble > > m_velocity
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72

References ASSERTL1, m_velocity, and Vmath::Vmul().

Referenced by v_InitObject().

◆ GetNormalVelocity()

Array< OneD, NekDouble > & Nektar::EigenValuesAdvection::GetNormalVelocity ( )
protected

Get the normal velocity.

Definition at line 123 of file EigenValuesAdvection.cpp.

124{
125 // Number of trace (interface) points
126 int nTracePts = GetTraceNpoints();
127
128 // Auxiliary variable to compute the normal velocity
129 Array<OneD, NekDouble> tmp(nTracePts);
130
131 // Reset the normal velocity
132 Vmath::Zero(nTracePts, m_traceVn, 1);
133
134 for (int i = 0; i < m_velocity.size(); ++i)
135 {
136 m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
137
138 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
139 m_traceVn, 1);
140 }
141
142 return m_traceVn;
143}
Array< OneD, NekDouble > m_traceVn
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, m_traceVn, m_velocity, Vmath::Vvtvp(), and Vmath::Zero().

Referenced by v_InitObject().

◆ v_DoInitialise()

void Nektar::EigenValuesAdvection::v_DoInitialise ( bool  dumpInitialConditions = false)
overrideprotectedvirtual

Virtual function for initialisation implementation.

By default, nothing needs initialising at the EquationSystem level.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 145 of file EigenValuesAdvection.cpp.

147{
148}

◆ v_DoSolve()

void Nektar::EigenValuesAdvection::v_DoSolve ( void  )
overrideprotectedvirtual

Virtual function for solve implementation.

Feeding the weak Advection oprator with a vector (inarray) Looping on inarray and changing the position of the only non-zero entry we simulate the multiplication by the identity matrix. The results stored in outarray is one of the columns of the weak advection oprators which are then stored in MATRIX for the futher eigenvalues calculation.

The result is stored in outarray (is the j-th columns of the weak advection operator). We now store it in MATRIX(j)

Set the j-th entry of inarray back to zero

Calulating the eigenvalues of the weak advection operator stored in (MATRIX) using Lapack routines

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 150 of file EigenValuesAdvection.cpp.

151{
152 int nvariables = 1;
153 int dofs = GetNcoeffs();
154 // bool UseContCoeffs = false;
155
156 Array<OneD, Array<OneD, NekDouble>> inarray(nvariables);
157 Array<OneD, Array<OneD, NekDouble>> tmp(nvariables);
158 Array<OneD, Array<OneD, NekDouble>> outarray(nvariables);
159 Array<OneD, Array<OneD, NekDouble>> WeakAdv(nvariables);
160
161 int npoints = GetNpoints();
162 int ncoeffs = GetNcoeffs();
163
164 switch (m_projectionType)
165 {
167 {
168 dofs = ncoeffs;
169 break;
170 }
173 {
174 // dofs = GetContNcoeffs();
175 // UseContCoeffs = true;
176 break;
177 }
178 }
179
180 cout << endl;
181 cout << "Num Phys Points = " << npoints << endl; // phisical points
182 cout << "Num Coeffs = " << ncoeffs << endl; //
183 cout << "Num Cont Coeffs = " << dofs << endl;
184
185 inarray[0] = Array<OneD, NekDouble>(npoints, 0.0);
186 outarray[0] = Array<OneD, NekDouble>(npoints, 0.0);
187 tmp[0] = Array<OneD, NekDouble>(npoints, 0.0);
188
189 WeakAdv[0] = Array<OneD, NekDouble>(ncoeffs, 0.0);
190 Array<OneD, NekDouble> MATRIX(npoints * npoints, 0.0);
191
192 for (int j = 0; j < npoints; j++)
193 {
194
195 inarray[0][j] = 1.0;
196
197 /// Feeding the weak Advection oprator with a vector (inarray)
198 /// Looping on inarray and changing the position of the only non-zero
199 /// entry we simulate the multiplication by the identity matrix. The
200 /// results stored in outarray is one of the columns of the weak
201 /// advection oprators which are then stored in MATRIX for the futher
202 /// eigenvalues calculation.
203 m_advObject->Advect(nvariables, m_fields, m_velocity, inarray, outarray,
204 0.0);
205 Vmath::Neg(npoints, outarray[0], 1);
206 switch (m_projectionType)
207 {
209 {
210 break;
211 }
214 {
215 for (int i = 0; i < nvariables; ++i)
216 {
217 // Projection
218 m_fields[i]->FwdTrans(outarray[i], WeakAdv[i]);
219
220 m_fields[i]->BwdTrans(WeakAdv[i], outarray[i]);
221 }
222 break;
223 }
224 }
225
226 /// The result is stored in outarray (is the j-th columns of the weak
227 /// advection operator). We now store it in MATRIX(j)
228 Vmath::Vcopy(npoints, &(outarray[0][0]), 1, &(MATRIX[j]), npoints);
229
230 /// Set the j-th entry of inarray back to zero
231 inarray[0][j] = 0.0;
232 }
233
234 ////////////////////////////////////////////////////////////////////////////////
235 /// Calulating the eigenvalues of the weak advection operator stored in
236 /// (MATRIX) using Lapack routines
237
238 char jobvl = 'N';
239 char jobvr = 'N';
240 int info = 0, lwork = 3 * npoints;
241 NekDouble dum;
242
243 Array<OneD, NekDouble> EIG_R(npoints);
244 Array<OneD, NekDouble> EIG_I(npoints);
245
246 Array<OneD, NekDouble> work(lwork);
247
248 Lapack::Dgeev(jobvl, jobvr, npoints, MATRIX.get(), npoints, EIG_R.get(),
249 EIG_I.get(), &dum, 1, &dum, 1, &work[0], lwork, info);
250
251 ////////////////////////////////////////////////////////
252 // Print Matrix
253 FILE *mFile;
254
255 mFile = fopen("WeakAdvMatrix.txt", "w");
256 for (int j = 0; j < npoints; j++)
257 {
258 for (int k = 0; k < npoints; k++)
259 {
260 fprintf(mFile, "%e ", MATRIX[j * npoints + k]);
261 }
262 fprintf(mFile, "\n");
263 }
264 fclose(mFile);
265
266 ////////////////////////////////////////////////////////
267 // Output of the EigenValues
268 FILE *pFile;
269
270 pFile = fopen("Eigenvalues.txt", "w");
271 for (int j = 0; j < npoints; j++)
272 {
273 fprintf(pFile, "%e %e\n", EIG_R[j], EIG_I[j]);
274 }
275 fclose(pFile);
276
277 cout << "\nEigenvalues : " << endl;
278 for (int j = 0; j < npoints; j++)
279 {
280 cout << EIG_R[j] << "\t" << EIG_I[j] << endl;
281 }
282 cout << endl;
283}
SolverUtils::AdvectionSharedPtr m_advObject
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
static void Dgeev(const char &uplo, const char &lrev, const int &n, const double *a, const int &lda, double *wr, double *wi, double *rev, const int &ldr, double *lev, const int &ldv, double *work, const int &lwork, int &info)
Solve general real matrix eigenproblem.
Definition: Lapack.hpp:329
double NekDouble
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Lapack::Dgeev(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetNpoints(), m_advObject, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, m_velocity, Vmath::Neg(), and Vmath::Vcopy().

◆ v_InitObject()

void Nektar::EigenValuesAdvection::v_InitObject ( bool  DeclareFeld = true)
overrideprotectedvirtual

Initialisation object for EquationSystem.

Continuous field

Setting up the normals

Setting up the normals

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 55 of file EigenValuesAdvection.cpp.

56{
57 EquationSystem::v_InitObject(DeclareFields);
58
59 // Define Velocity fields
60 m_velocity = Array<OneD, Array<OneD, NekDouble>>(m_spacedim);
61 std::vector<std::string> vel;
62 vel.push_back("Vx");
63 vel.push_back("Vy");
64 vel.push_back("Vz");
65 vel.resize(m_spacedim);
66
67 GetFunction("AdvectionVelocity")->Evaluate(vel, m_velocity);
68
69 // Type of advection class to be used
70 switch (m_projectionType)
71 {
72 // Discontinuous field
74 {
75 // Define the normal velocity fields
76 if (m_fields[0]->GetTrace())
77 {
78 m_traceVn = Array<OneD, NekDouble>(GetTraceNpoints());
79 }
80
81 string advName;
82 string riemName;
83 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
85 advName, advName);
87 this);
88 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
91 riemName, m_session);
92 m_riemannSolver->SetScalar(
94
95 m_advObject->SetRiemannSolver(m_riemannSolver);
96 m_advObject->InitObject(m_session, m_fields);
97 break;
98 }
99 // Continuous field
102 {
103 string advName;
104 m_session->LoadSolverInfo("AdvectionType", advName,
105 "NonConservative");
107 advName, advName);
109 this);
110 break;
111 }
112 default:
113 {
114 ASSERTL0(false, "Unsupported projection type.");
115 break;
116 }
117 }
118}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
int m_spacedim
Spatial dimension (>= expansion dim).
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareFeld=true)
Initialisation object for EquationSystem.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
RiemannSolverFactory & GetRiemannSolverFactory()

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::GetAdvectionFactory(), GetFluxVector(), Nektar::SolverUtils::EquationSystem::GetFunction(), GetNormalVelocity(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), m_advObject, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_traceVn, m_velocity, and Nektar::SolverUtils::EquationSystem::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< EigenValuesAdvection >

friend class MemoryManager< EigenValuesAdvection >
friend

Definition at line 1 of file EigenValuesAdvection.h.

Member Data Documentation

◆ className

string Nektar::EigenValuesAdvection::className
static
Initial value:
=
"EigenValuesAdvection", EigenValuesAdvection::create,
"Eigenvalues of the weak advection operator.")
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 63 of file EigenValuesAdvection.h.

◆ m_advObject

SolverUtils::AdvectionSharedPtr Nektar::EigenValuesAdvection::m_advObject
protected

Definition at line 70 of file EigenValuesAdvection.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_riemannSolver

SolverUtils::RiemannSolverSharedPtr Nektar::EigenValuesAdvection::m_riemannSolver
protected

Definition at line 68 of file EigenValuesAdvection.h.

Referenced by v_InitObject().

◆ m_traceVn

Array<OneD, NekDouble> Nektar::EigenValuesAdvection::m_traceVn
protected

Definition at line 71 of file EigenValuesAdvection.h.

Referenced by GetNormalVelocity(), and v_InitObject().

◆ m_velocity

Array<OneD, Array<OneD, NekDouble> > Nektar::EigenValuesAdvection::m_velocity
protected

Definition at line 69 of file EigenValuesAdvection.h.

Referenced by GetFluxVector(), GetNormalVelocity(), v_DoSolve(), and v_InitObject().