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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::EigenValuesAdvection:
Collaboration graph
[legend]

Public Member Functions

virtual ~EigenValuesAdvection ()
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 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 NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
boost::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 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 void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name. More...
 
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow. 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 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 WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $. More...
 
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form. 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 ScanForHistoryPoints ()
 Builds map of which element holds each history point. More...
 
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to 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 GetFinalTime ()
 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 GetNumElmVelocity ()
 
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, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void SetStepsToOne ()
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const string &s1, const string &s2)
 Perform a case-insensitive string comparison. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Static Public Member Functions

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

Static Public Attributes

static std::string className = GetEquationSystemFactory().RegisterCreatorFunction("EigenValuesAdvection", EigenValuesAdvection::create, "Eigenvalues of the weak advection operator.")
 Name of class. More...
 

Protected Member Functions

 EigenValuesAdvection (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual void v_InitObject ()
 Initialisation object for EquationSystem. More...
 
virtual void v_DoInitialise ()
 Virtual function for initialisation implementation. More...
 
virtual void v_DoSolve ()
 Virtual function for solve implementation. More...
 
virtual void v_GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
virtual void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const string &s1, const string &s2)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. 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)
 
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
 
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

Array< OneD, Array< OneD, NekDouble > > m_velocity
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
map< std::string, Array< OneD, Array< OneD, float > > > m_interpWeights
 Map of the interpolation weights for a specific filename. More...
 
map< std::string, Array< OneD, Array< OneD, unsigned int > > > m_interpInds
 Map of the interpolation indices for a specific filename. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_base
 Base fields. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_derivedfields
 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...
 
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...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. 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, Array< OneD, Array< OneD, NekDouble > > > m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_tanbasis
 2 x m_spacedim x nq 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...
 
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...
 
int m_NumMode
 Mode to use in case of single mode analysis. 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...
 

Detailed Description

Definition at line 45 of file EigenValuesAdvection.h.

Constructor & Destructor Documentation

Nektar::EigenValuesAdvection::~EigenValuesAdvection ( )
virtual

Definition at line 70 of file EigenValuesAdvection.cpp.

71  {
72 
73  }
Nektar::EigenValuesAdvection::EigenValuesAdvection ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 44 of file EigenValuesAdvection.cpp.

46  : EquationSystem(pSession)
47  {
48  }
SOLVER_UTILS_EXPORT EquationSystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises EquationSystem class members.

Member Function Documentation

static EquationSystemSharedPtr Nektar::EigenValuesAdvection::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 51 of file EigenValuesAdvection.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

52  {
54  p->InitObject();
55  return p;
56  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
void Nektar::EigenValuesAdvection::v_DoInitialise ( void  )
protectedvirtual

Virtual function for initialisation implementation.

By default, nothing needs initialising at the EquationSystem level.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 65 of file EigenValuesAdvection.cpp.

66  {
67 
68  }
void Nektar::EigenValuesAdvection::v_DoSolve ( void  )
protectedvirtual

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 75 of file EigenValuesAdvection.cpp.

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

76  {
77  int nvariables = 1;
78  int i,dofs = GetNcoeffs();
79  //bool UseContCoeffs = false;
80 
81  Array<OneD, Array<OneD, NekDouble> > inarray(nvariables);
82  Array<OneD, Array<OneD, NekDouble> > tmp(nvariables);
83  Array<OneD, Array<OneD, NekDouble> > outarray(nvariables);
84  Array<OneD, Array<OneD, NekDouble> > WeakAdv(nvariables);
85 
86  int npoints = GetNpoints();
87  int ncoeffs = GetNcoeffs();
88 
89  switch (m_projectionType)
90  {
92  {
93  dofs = ncoeffs;
94  break;
95  }
98  {
99  //dofs = GetContNcoeffs();
100  //UseContCoeffs = true;
101  break;
102  }
103  }
104 
105  cout << endl;
106  cout << "Num Phys Points = " << npoints << endl; // phisical points
107  cout << "Num Coeffs = " << ncoeffs << endl; //
108  cout << "Num Cont Coeffs = " << dofs << endl;
109 
110  inarray[0] = Array<OneD, NekDouble>(npoints,0.0);
111  outarray[0] = Array<OneD, NekDouble>(npoints,0.0);
112  tmp[0] = Array<OneD, NekDouble>(npoints,0.0);
113 
114  WeakAdv[0] = Array<OneD, NekDouble>(ncoeffs,0.0);
115  Array<OneD, NekDouble> MATRIX(npoints*npoints,0.0);
116 
117  for (int j = 0; j < npoints; j++)
118  {
119 
120  inarray[0][j] = 1.0;
121 
122  /// Feeding the weak Advection oprator with a vector (inarray)
123  /// Looping on inarray and changing the position of the only non-zero entry
124  /// we simulate the multiplication by the identity matrix.
125  /// The results stored in outarray is one of the columns of the weak advection oprators
126  /// which are then stored in MATRIX for the futher eigenvalues calculation.
127 
128  switch (m_projectionType)
129  {
131  {
132  WeakDGAdvection(inarray, WeakAdv,true,true,1);
133 
134  m_fields[0]->MultiplyByElmtInvMass(WeakAdv[0],WeakAdv[0]);
135 
136  m_fields[0]->BwdTrans(WeakAdv[0],outarray[0]);
137 
138  Vmath::Neg(npoints,outarray[0],1);
139  break;
140  }
143  {
144  // Calculate -V\cdot Grad(u);
145  for(i = 0; i < nvariables; ++i)
146  {
147  //Projection
148  m_fields[i]->FwdTrans(inarray[i],WeakAdv[i]);
149 
150  m_fields[i]->BwdTrans_IterPerExp(WeakAdv[i],tmp[i]);
151 
152  //Advection operator
153  AdvectionNonConservativeForm(m_velocity,tmp[i],outarray[i]);
154 
155  Vmath::Neg(npoints,outarray[i],1);
156 
157  //m_fields[i]->MultiplyByInvMassMatrix(WeakAdv[i],WeakAdv[i]);
158  //Projection
159  m_fields[i]->FwdTrans(outarray[i],WeakAdv[i]);
160 
161  m_fields[i]->BwdTrans_IterPerExp(WeakAdv[i],outarray[i]);
162  }
163  break;
164  }
165  }
166 
167  /// The result is stored in outarray (is the j-th columns of the weak advection operator).
168  /// We now store it in MATRIX(j)
169  Vmath::Vcopy(npoints,&(outarray[0][0]),1,&(MATRIX[j]),npoints);
170 
171  /// Set the j-th entry of inarray back to zero
172  inarray[0][j] = 0.0;
173  }
174 
175  ////////////////////////////////////////////////////////////////////////////////
176  /// Calulating the eigenvalues of the weak advection operator stored in (MATRIX)
177  /// using Lapack routines
178 
179  char jobvl = 'N';
180  char jobvr = 'N';
181  int info = 0, lwork = 3*npoints;
182  NekDouble dum;
183 
184  Array<OneD, NekDouble> EIG_R(npoints);
185  Array<OneD, NekDouble> EIG_I(npoints);
186 
187  Array<OneD, NekDouble> work(lwork);
188 
189  Lapack::Dgeev(jobvl,jobvr,npoints,MATRIX.get(),npoints,EIG_R.get(),EIG_I.get(),&dum,1,&dum,1,&work[0],lwork,info);
190 
191  ////////////////////////////////////////////////////////
192  //Print Matrix
193  FILE *mFile;
194 
195  mFile = fopen ("WeakAdvMatrix.txt","w");
196  for(int j = 0; j<npoints; j++)
197  {
198  for(int k = 0; k<npoints; k++)
199  {
200  fprintf(mFile,"%e ",MATRIX[j*npoints+k]);
201  }
202  fprintf(mFile,"\n");
203  }
204  fclose (mFile);
205 
206  ////////////////////////////////////////////////////////
207  //Output of the EigenValues
208  FILE *pFile;
209 
210  pFile = fopen ("Eigenvalues.txt","w");
211  for(int j = 0; j<npoints; j++)
212  {
213  fprintf(pFile,"%e %e\n",EIG_R[j],EIG_I[j]);
214  }
215  fclose (pFile);
216 
217  cout << "\nEigenvalues : " << endl;
218  for(int j = 0; j<npoints; j++)
219  {
220  cout << EIG_R[j] << "\t" << EIG_I[j] << endl;
221  }
222  cout << endl;
223  }
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm(const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
Compute the non-conservative advection.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
double NekDouble
Array< OneD, Array< OneD, NekDouble > > m_velocity
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNcoeffs()
SOLVER_UTILS_EXPORT void WeakDGAdvection(const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
Calculate the weak discontinuous Galerkin advection.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
void Nektar::EigenValuesAdvection::v_GetFluxVector ( const int  i,
Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 225 of file EigenValuesAdvection.cpp.

References ASSERTL1, Nektar::SolverUtils::EquationSystem::GetNpoints(), m_velocity, and Vmath::Vmul().

227  {
228  ASSERTL1(flux.num_elements() == m_velocity.num_elements(),"Dimension of flux array and velocity array do not match");
229 
230  for(int j = 0; j < flux.num_elements(); ++j)
231  {
232  Vmath::Vmul(GetNpoints(),physfield[i],1,
233  m_velocity[j],1,flux[j],1);
234  }
235  }
Array< OneD, Array< OneD, NekDouble > > m_velocity
SOLVER_UTILS_EXPORT int GetNpoints()
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
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.cpp:169
void Nektar::EigenValuesAdvection::v_InitObject ( )
protectedvirtual

Initialisation object for EquationSystem.

Continuous field

Setting up the normals

Setting up the normals

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 50 of file EigenValuesAdvection.cpp.

References Nektar::SolverUtils::EquationSystem::EvaluateFunction(), Nektar::SolverUtils::EquationSystem::m_spacedim, and m_velocity.

51  {
52  EquationSystem::v_InitObject();
53 
54  // Define Velocity fields
56  std::vector<std::string> vel;
57  vel.push_back("Vx");
58  vel.push_back("Vy");
59  vel.push_back("Vz");
60  vel.resize(m_spacedim);
61 
62  EvaluateFunction(vel, m_velocity, "AdvectionVelocity");
63  }
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT void EvaluateFunction(Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
Evaluates a function as specified in the session file.
Array< OneD, Array< OneD, NekDouble > > m_velocity
void Nektar::EigenValuesAdvection::v_NumericalFlux ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numflux 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 237 of file EigenValuesAdvection.cpp.

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

238  {
239  int i;
240 
241  int nTraceNumPoints = GetTraceNpoints();
242  int nvel = m_spacedim; //m_velocity.num_elements();
243 
244  Array<OneD, NekDouble > Fwd(nTraceNumPoints);
245  Array<OneD, NekDouble > Bwd(nTraceNumPoints);
246  Array<OneD, NekDouble > Vn (nTraceNumPoints,0.0);
247 
248  //Get Edge Velocity - Could be stored if time independent
249  for(i = 0; i < nvel; ++i)
250  {
251  m_fields[0]->ExtractTracePhys(m_velocity[i], Fwd);
252  Vmath::Vvtvp(nTraceNumPoints,m_traceNormals[i],1,Fwd,1,Vn,1,Vn,1);
253  }
254 
255  for(i = 0; i < numflux.num_elements(); ++i)
256  {
257  m_fields[i]->GetFwdBwdTracePhys(physfield[i],Fwd,Bwd);
258  m_fields[i]->GetTrace()->Upwind(Vn,Fwd,Bwd,numflux[i]);
259  // calculate m_fields[i]*Vn
260  Vmath::Vmul(nTraceNumPoints,numflux[i],1,Vn,1,numflux[i],1);
261  }
262  }
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.cpp:428
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
int m_spacedim
Spatial dimension (>= expansion dim).
Array< OneD, Array< OneD, NekDouble > > m_velocity
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
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.cpp:169

Friends And Related Function Documentation

friend class MemoryManager< EigenValuesAdvection >
friend

Definition at line 48 of file EigenValuesAdvection.h.

Member Data Documentation

string Nektar::EigenValuesAdvection::className = GetEquationSystemFactory().RegisterCreatorFunction("EigenValuesAdvection", EigenValuesAdvection::create, "Eigenvalues of the weak advection operator.")
static

Name of class.

Definition at line 58 of file EigenValuesAdvection.h.

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

Definition at line 63 of file EigenValuesAdvection.h.

Referenced by v_DoSolve(), v_GetFluxVector(), v_InitObject(), and v_NumericalFlux().