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

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 >
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 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 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 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 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 void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. 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, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 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, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
Array< OneD, NekDouble > & GetNormalVelocity ()
 Get the normal velocity. More...
 
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...
 
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 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 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...
 
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...
 
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_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...
 
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 []
 

Detailed Description

Definition at line 46 of file EigenValuesAdvection.h.

Constructor & Destructor Documentation

◆ ~EigenValuesAdvection()

Nektar::EigenValuesAdvection::~EigenValuesAdvection ( )
virtual

Definition at line 152 of file EigenValuesAdvection.cpp.

153  {
154 
155  }

◆ EigenValuesAdvection()

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

Definition at line 45 of file EigenValuesAdvection.cpp.

48  : EquationSystem(pSession, pGraph)
49  {
50  }
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.

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

55  {
57  ::AllocateSharedPtr(pSession, pGraph);
58  p->InitObject();
59  return p;
60  }
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.

◆ GetFluxVector()

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

Definition at line 291 of file EigenValuesAdvection.cpp.

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

Referenced by v_InitObject().

294  {
295  ASSERTL1(flux[0].num_elements() == m_velocity.num_elements(),
296  "Dimension of flux array and velocity array do not match");
297 
298  int nq = physfield[0].num_elements();
299 
300  for (int i = 0; i < flux.num_elements(); ++i)
301  {
302  for (int j = 0; j < flux[0].num_elements(); ++j)
303  {
304  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
305  flux[i][j], 1);
306  }
307  }
308  }
Array< OneD, Array< OneD, NekDouble > > m_velocity
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
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:186

◆ GetNormalVelocity()

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

Get the normal velocity.

Definition at line 122 of file EigenValuesAdvection.cpp.

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().

123  {
124  // Number of trace (interface) points
125  int nTracePts = GetTraceNpoints();
126 
127  // Auxiliary variable to compute the normal velocity
128  Array<OneD, NekDouble> tmp(nTracePts);
129 
130  // Reset the normal velocity
131  Vmath::Zero(nTracePts, m_traceVn, 1);
132 
133  for (int i = 0; i < m_velocity.num_elements(); ++i)
134  {
135  m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
136 
137  Vmath::Vvtvp(nTracePts,
138  m_traceNormals[i], 1,
139  tmp, 1,
140  m_traceVn, 1,
141  m_traceVn, 1);
142  }
143 
144  return m_traceVn;
145  }
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:445
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, Array< OneD, NekDouble > > m_velocity
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
Array< OneD, NekDouble > m_traceVn

◆ v_DoInitialise()

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

148  {
149 
150  }

◆ v_DoSolve()

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

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().

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

◆ v_InitObject()

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

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().

53  {
55 
56  // Define Velocity fields
57  m_velocity = Array<OneD, Array<OneD, NekDouble> >(m_spacedim);
58  std::vector<std::string> vel;
59  vel.push_back("Vx");
60  vel.push_back("Vy");
61  vel.push_back("Vz");
62  vel.resize(m_spacedim);
63 
64  GetFunction( "AdvectionVelocity")->Evaluate(vel, m_velocity);
65 
66  // Type of advection class to be used
67  switch(m_projectionType)
68  {
69  // Discontinuous field
71  {
72  // Define the normal velocity fields
73  if (m_fields[0]->GetTrace())
74  {
75  m_traceVn = Array<OneD, NekDouble>(GetTraceNpoints());
76  }
77 
78  string advName;
79  string riemName;
80  m_session->LoadSolverInfo(
81  "AdvectionType", advName, "WeakDG");
83  GetAdvectionFactory().CreateInstance(advName, advName);
84  m_advObject->SetFluxVector(
86  m_session->LoadSolverInfo(
87  "UpwindType", riemName, "Upwind");
90  riemName, m_session);
91  m_riemannSolver->SetScalar(
93 
94  m_advObject->SetRiemannSolver(m_riemannSolver);
95  m_advObject->InitObject(m_session, m_fields);
96  break;
97  }
98  // Continuous field
101  {
102  string advName;
103  m_session->LoadSolverInfo(
104  "AdvectionType", advName, "NonConservative");
106  GetAdvectionFactory().CreateInstance(advName, advName);
107  m_advObject->SetFluxVector(
109  break;
110  }
111  default:
112  {
113  ASSERTL0(false, "Unsupported projection type.");
114  break;
115  }
116  }
117  }
SolverUtils::AdvectionSharedPtr m_advObject
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Initialisation object for EquationSystem.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
RiemannSolverFactory & GetRiemannSolverFactory()
int m_spacedim
Spatial dimension (>= expansion dim).
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
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.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
Array< OneD, Array< OneD, NekDouble > > m_velocity
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, NekDouble > m_traceVn

Friends And Related Function Documentation

◆ MemoryManager< EigenValuesAdvection >

friend class MemoryManager< EigenValuesAdvection >
friend

Definition at line 49 of file EigenValuesAdvection.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 62 of file EigenValuesAdvection.h.

◆ m_advObject

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

Definition at line 69 of file EigenValuesAdvection.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_riemannSolver

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

Definition at line 67 of file EigenValuesAdvection.h.

Referenced by v_InitObject().

◆ m_traceVn

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

Definition at line 70 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 68 of file EigenValuesAdvection.h.

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