Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::CFLtester Class Reference

#include <CFLtester.h>

Inheritance diagram for Nektar::CFLtester:
Inheritance graph
[legend]
Collaboration diagram for Nektar::CFLtester:
Collaboration graph
[legend]

Public Member Functions

virtual ~CFLtester ()
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable.
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor.
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object.
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem.
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem.
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space.
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space.
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve.
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation.
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name.
template<class T >
boost::shared_ptr< T > as ()
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name.
SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr 
GetSession ()
 Get Session name.
SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
GetPressure ()
 Get pressure field if available.
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics.
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda.
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.
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.
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.
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.
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow.
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields.
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution.
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.
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields.
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf].
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $.
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) $.
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.
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.
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.
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields.
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.
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields.
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename.
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.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file.
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.
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.
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.
SOLVER_UTILS_EXPORT void ScanForHistoryPoints ()
 Builds map of which element holds each history point.
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to file.
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary.
SOLVER_UTILS_EXPORT Array
< OneD,
MultiRegions::ExpListSharedPtr > & 
UpdateFields ()
SOLVER_UTILS_EXPORT
LibUtilities::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated.
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time.
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.

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)

Static Public Attributes

static std::string className = GetEquationSystemFactory().RegisterCreatorFunction("CFLtester", CFLtester::create, "Testing CFL restriction")

Protected Member Functions

 CFLtester (const LibUtilities::SessionReaderSharedPtr &pSession)
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
virtual void v_InitObject ()
 Init object for UnsteadySystem class.
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)
virtual void v_GenerateSummary (SummaryList &s)
 Print a summary of time stepping parameters.
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members.
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control.
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file.
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Return the timestep to be used for the next step in the time-marching loop.
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members.
int nocase_cmp (const string &s1, const string &s2)
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time.
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.
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.
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space.
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space.
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::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status information.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme.
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use.
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
NekDouble m_epsilon
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used.
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used.
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used.
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state.
std::vector< int > m_intVariables
std::vector< FilterSharedPtrm_filters
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator.
LibUtilities::SessionReaderSharedPtr m_session
 The session reader.
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 Array holding all dependent variables.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object.
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh.
std::string m_sessionName
 Name of the session.
NekDouble m_time
 Current time of simulation.
NekDouble m_fintime
 Finish time of the simulation.
NekDouble m_timestep
 Time step size.
NekDouble m_lambda
 Lambda constant in real system if one required.
NekDouble m_checktime
 Time between checkpoints.
int m_steps
 Number of steps to take.
int m_checksteps
 Number of steps between checkpoints.
int m_spacedim
 Spatial dimension (>= expansion dim).
int m_expdim
 Expansion dimension.
bool m_SingleMode
 Flag to determine if single homogeneous mode is used.
bool m_HalfMode
 Flag to determine if half homogeneous mode is used.
bool m_MultipleModes
 Flag to determine if use multiple homogenenous modes are used.
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform.
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations.
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous.
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields.
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error.
enum HomogeneousType m_HomogeneousType
NekDouble m_LhomX
 physical length in X direction (if homogeneous)
NekDouble m_LhomY
 physical length in Y direction (if homogeneous)
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous)
int m_npointsX
 number of points in X direction (if homogeneous)
int m_npointsY
 number of points in Y direction (if homogeneous)
int m_npointsZ
 number of points in Z direction (if homogeneous)
int m_HomoDirec
 number of homogenous directions
int m_NumMode
 Mode to use in case of single mode analysis.

Private Member Functions

virtual NekDouble v_GetTimeStep (const Array< OneD, int > ExpOrder, const Array< OneD, NekDouble > CFL, NekDouble timeCFL)
virtual NekDouble v_GetTimeStep (int ExpOrder, NekDouble CFL, NekDouble TimeStability)
Array< OneD, NekDoubleGetStdVelocity (const Array< OneD, Array< OneD, NekDouble > > inarray)

Friends

class MemoryManager< CFLtester >

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1).
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous }
 Parameter for homogeneous expansions. More...

Detailed Description

Definition at line 81 of file CFLtester.h.

Constructor & Destructor Documentation

Nektar::CFLtester::~CFLtester ( )
virtual

Definition at line 78 of file CFLtester.cpp.

{
}
Nektar::CFLtester::CFLtester ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 48 of file CFLtester.cpp.

: UnsteadySystem(pSession)
{
}

Member Function Documentation

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

Definition at line 86 of file CFLtester.h.

{
p->InitObject();
return p;}
void Nektar::CFLtester::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Definition at line 131 of file CFLtester.cpp.

References ASSERTL0, 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, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Vcopy().

Referenced by v_InitObject().

{
int j;
int nvariables = inarray.num_elements();
{
{
// Just copy over array
int npoints = GetNpoints();
for(j = 0; j < nvariables; ++j)
{
Vmath::Vcopy(npoints,inarray[j],1,outarray[j],1);
}
}
break;
{
Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs());
for(j = 0; j < nvariables; ++j)
{
m_fields[j]->FwdTrans(inarray[j],coeffs);
m_fields[j]->BwdTrans_IterPerExp(coeffs,outarray[j]);
}
break;
}
default:
ASSERTL0(false,"Unknown projection scheme");
break;
}
}
void Nektar::CFLtester::DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Definition at line 82 of file CFLtester.cpp.

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

Referenced by v_InitObject().

{
int j;
int nvariables = inarray.num_elements();
int npoints = GetNpoints();
{
{
int ncoeffs = inarray[0].num_elements();
Array<OneD, Array<OneD, NekDouble> > WeakAdv(nvariables);
WeakAdv[0] = Array<OneD, NekDouble>(ncoeffs*nvariables);
for(j = 1; j < nvariables; ++j)
{
WeakAdv[j] = WeakAdv[j-1] + ncoeffs;
}
WeakDGAdvection(inarray, WeakAdv, true, true);
for(j = 0; j < nvariables; ++j)
{
m_fields[j]->MultiplyByElmtInvMass(WeakAdv[j], WeakAdv[j]);
m_fields[j]->BwdTrans(WeakAdv[j],outarray[j]);
Vmath::Neg(npoints,outarray[j],1);
}
break;
}
{
// Calculate -V\cdot Grad(u);
for(j = 0; j < nvariables; ++j)
{
inarray[j],
outarray[j]);
Vmath::Neg(npoints, outarray[j], 1);
}
break;
}
}
}
Array< OneD, NekDouble > Nektar::CFLtester::GetStdVelocity ( const Array< OneD, Array< OneD, NekDouble > >  inarray)
private

Definition at line 329 of file CFLtester.cpp.

References ASSERTL0, Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_expdim, and Nektar::SolverUtils::EquationSystem::m_fields.

Referenced by v_GetTimeStep().

{
// Checking if the problem is 2D
ASSERTL0(m_expdim >= 2, "Method not implemented for 1D");
int nTotQuadPoints = GetTotPoints();
int n_element = m_fields[0]->GetExpSize();
int nvel = inarray.num_elements();
NekDouble pntVelocity;
// Getting the standard velocity vector on the 2D normal space
Array<OneD, Array<OneD, NekDouble> > stdVelocity(nvel);
Array<OneD, NekDouble> stdV(n_element, 0.0);
for (int i = 0; i < nvel; ++i)
{
stdVelocity[i] = Array<OneD, NekDouble>(nTotQuadPoints);
}
if (nvel == 2)
{
for (int el = 0; el < n_element; ++el)
{
m_fields[0]->GetExp(el)->as<LocalRegions::Expansion2D>();
int n_points = el2D->GetTotPoints();
Array<OneD, const NekDouble> jac =
el2D->GetGeom2D()->GetMetricInfo()->GetJac();
Array<TwoD, const NekDouble> gmat =
el2D->GetGeom2D()->GetMetricInfo()->GetDerivFactors();
if (el2D->GetGeom2D()->GetMetricInfo()->GetGtype()
{
for (int i = 0; i < n_points; i++)
{
stdVelocity[0][i] = gmat[0][i]*inarray[0][i]
+ gmat[2][i]*inarray[1][i];
stdVelocity[1][i] = gmat[1][i]*inarray[0][i]
+ gmat[3][i]*inarray[1][i];
}
}
else
{
for (int i = 0; i < n_points; i++)
{
stdVelocity[0][i] = gmat[0][0]*inarray[0][i]
+ gmat[2][0]*inarray[1][i];
stdVelocity[1][i] = gmat[1][0]*inarray[0][i]
+ gmat[3][0]*inarray[1][i];
}
}
for (int i = 0; i < n_points; i++)
{
pntVelocity = sqrt(stdVelocity[0][i]*stdVelocity[0][i]
+ stdVelocity[1][i]*stdVelocity[1][i]);
if (pntVelocity>stdV[el])
{
stdV[el] = pntVelocity;
}
}
}
}
else
{
for (int el = 0; el < n_element; ++el)
{
m_fields[0]->GetExp(el)->as<LocalRegions::Expansion3D>();
int n_points = el3D->GetTotPoints();
Array<OneD, const NekDouble> jac =
el3D->GetGeom3D()->GetMetricInfo()->GetJac();
Array<TwoD, const NekDouble> gmat =
el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
if (el3D->GetGeom3D()->GetMetricInfo()->GetGtype()
{
for (int i = 0; i < n_points; i++)
{
stdVelocity[0][i] = gmat[0][i]*inarray[0][i]
+ gmat[3][i]*inarray[1][i]
+ gmat[6][i]*inarray[2][i];
stdVelocity[1][i] = gmat[1][i]*inarray[0][i]
+ gmat[4][i]*inarray[1][i]
+ gmat[7][i]*inarray[2][i];
stdVelocity[2][i] = gmat[2][i]*inarray[0][i]
+ gmat[5][i]*inarray[1][i]
+ gmat[8][i]*inarray[2][i];
}
}
else
{
Array<OneD, const NekDouble> jac =
el3D->GetGeom3D()->GetMetricInfo()->GetJac();
Array<TwoD, const NekDouble> gmat =
el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
for (int i = 0; i < n_points; i++)
{
stdVelocity[0][i] = gmat[0][0]*inarray[0][i]
+ gmat[3][0]*inarray[1][i]
+ gmat[6][0]*inarray[2][i];
stdVelocity[1][i] = gmat[1][0]*inarray[0][i]
+ gmat[4][0]*inarray[1][i]
+ gmat[7][0]*inarray[2][i];
stdVelocity[2][i] = gmat[2][0]*inarray[0][i]
+ gmat[5][0]*inarray[1][i]
+ gmat[8][0]*inarray[2][i];
}
}
for (int i = 0; i < n_points; i++)
{
pntVelocity = sqrt(stdVelocity[0][i]*stdVelocity[0][i]
+ stdVelocity[1][i]*stdVelocity[1][i]
+ stdVelocity[2][i]*stdVelocity[2][i]);
if (pntVelocity > stdV[el])
{
stdV[el] = pntVelocity;
}
}
}
}
return stdV;
}
void Nektar::CFLtester::v_GenerateSummary ( SummaryList s)
protectedvirtual

Print a summary of time stepping parameters.

Prints a summary with some information regards the time-stepping.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 222 of file CFLtester.cpp.

void Nektar::CFLtester::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 172 of file CFLtester.cpp.

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

{
ASSERTL1(flux.num_elements() == m_velocity.num_elements(),
"Dimension of flux array and velocity array do not match");
for (int j = 0; j < flux.num_elements(); ++j)
{
physfield[i], 1,
m_velocity[j], 1,
flux[j], 1);
}
}
NekDouble Nektar::CFLtester::v_GetTimeStep ( const Array< OneD, int >  ExpOrder,
const Array< OneD, NekDouble CFL,
NekDouble  timeCFL 
)
privatevirtual

Definition at line 229 of file CFLtester.cpp.

References GetStdVelocity(), Nektar::SolverUtils::EquationSystem::m_fields, m_velocity, and Vmath::Vmin().

{
int n_element = m_fields[0]->GetExpSize();
//const NekDouble minLengthStdTri = 1.414213;
//const NekDouble minLengthStdQuad = 2.0;
//const NekDouble cLambda = 0.2; // Spencer book pag. 317
Array<OneD, NekDouble> tstep (n_element, 0.0);
Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
stdVelocity = GetStdVelocity(m_velocity);
for (int el = 0; el < n_element; ++el)
{
int npoints = m_fields[0]->GetExp(el)->GetTotPoints();
Array<OneD, NekDouble> one2D(npoints, 1.0);
//NekDouble Area = m_fields[0]->GetExp(el)->Integral(one2D);
if(boost::dynamic_pointer_cast<LocalRegions::TriExp>(m_fields[0]->GetExp(el)))
{
//tstep[el] = timeCFL/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
//tstep[el] = timeCFL*minLengthStdTri/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
tstep[el] = CFL[el]/(stdVelocity[el]);
}
else if(boost::dynamic_pointer_cast<LocalRegions::QuadExp>(m_fields[0]->GetExp(el)))
{
//tstep[el] = timeCFL/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
//tstep[el] = timeCFL*minLengthStdQuad/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
tstep[el] = CFL[el]/(stdVelocity[el]);
}
}
NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
return TimeStep;
}
NekDouble Nektar::CFLtester::v_GetTimeStep ( int  ExpOrder,
NekDouble  CFL,
NekDouble  TimeStability 
)
privatevirtual

Definition at line 271 of file CFLtester.cpp.

References ASSERTL0, Nektar::EigenvaluesAnaMeshesAB2, Nektar::EigenvaluesAnaMeshesRK2, Nektar::EigenvaluesAnaMeshesRK4, and Nektar::SolverUtils::EquationSystem::m_fields.

{
//================================================================
// This function has been created just to test specific problems, hence is not general
// and it has been implemented in a rude fashion, as the full CFLtester class.
// For real CFL calculations refer to the general implementation above. (A.Bolis)
//================================================================
NekDouble TimeStep;
NekDouble SpatialStability;
int n_elements = m_fields[0]->GetExpSize();
//solve ambiguity in windows
NekDouble n_elem = n_elements;
NekDouble DH = sqrt(n_elem);
int H = (int)DH;
int P = ExpOrder-1;
//================================================================
// Regular meshes
//SpatialStability = EigenvaluesRegMeshes[H-1][P-1];
//================================================================
// Anisotropic meshes
if (TimeStability == 1.0)
{
SpatialStability = EigenvaluesAnaMeshesAB2[H/2][P-1];
}
else if (TimeStability == 2.0)
{
SpatialStability = EigenvaluesAnaMeshesRK2[H/2][P-1];
}
else if (TimeStability == 2.784)
{
SpatialStability = EigenvaluesAnaMeshesRK4[H/2][P-1];
}
else
{
ASSERTL0(false,"Dominant eigenvalues database not present for this time-stepping scheme")
}
//================================================================
TimeStep = (TimeStability/SpatialStability)*CFL;
//================================================================
return TimeStep;
}
void Nektar::CFLtester::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 54 of file CFLtester.cpp.

References ASSERTL0, Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoOdeProjection(), DoOdeRhs(), Nektar::SolverUtils::EquationSystem::EvaluateFunction(), Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_spacedim, and m_velocity.

{
m_velocity = Array<OneD, Array<OneD, NekDouble> >(m_spacedim);
std::vector<std::string> vel;
vel.push_back("Vx");
vel.push_back("Vy");
vel.push_back("Vz");
vel.resize(m_spacedim);
EvaluateFunction(vel, m_velocity, "AdvectionVelocity");
{
}
else
{
ASSERTL0(false, "Implicit unsteady Advection not set up.");
}
}
void Nektar::CFLtester::v_NumericalFlux ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numflux 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 189 of file CFLtester.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().

{
int i;
int nTraceNumPoints = GetTraceNpoints();
int nvel = m_spacedim;
Array<OneD, NekDouble > Fwd(nTraceNumPoints);
Array<OneD, NekDouble > Bwd(nTraceNumPoints);
Array<OneD, NekDouble > Vn (nTraceNumPoints,0.0);
// Get Edge Velocity - Could be stored if time independent
for (i = 0; i < nvel; ++i)
{
m_fields[0]->ExtractTracePhys(m_velocity[i], Fwd);
Vmath::Vvtvp(nTraceNumPoints,
Fwd, 1, Vn, 1, Vn, 1);
}
for (i = 0; i < numflux.num_elements(); ++i)
{
m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
//evaulate upwinded m_fields[i]
m_fields[i]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux[i]);
// calculate m_fields[i]*Vn
Vmath::Vmul(nTraceNumPoints, numflux[i], 1, Vn, 1, numflux[i], 1);
}
}

Friends And Related Function Documentation

friend class MemoryManager< CFLtester >
friend

Definition at line 84 of file CFLtester.h.

Member Data Documentation

string Nektar::CFLtester::className = GetEquationSystemFactory().RegisterCreatorFunction("CFLtester", CFLtester::create, "Testing CFL restriction")
static

Definition at line 91 of file CFLtester.h.

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

Definition at line 97 of file CFLtester.h.

Referenced by DoOdeRhs(), v_GetFluxVector(), v_GetTimeStep(), v_InitObject(), and v_NumericalFlux().