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

#include <NonlinearPeregrine.h>

Inheritance diagram for Nektar::NonlinearPeregrine:
[legend]

Public Member Functions

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

Static Public Member Functions

static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::NonlinearSWE
static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::ShallowWaterSystem
static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::NonlinearSWE
static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::ShallowWaterSystem
static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 

Protected Member Functions

 NonlinearPeregrine (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
void v_InitObject (bool DeclareFields=true) override
 Init object for UnsteadySystem class. More...
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0) override
 Set the initial conditions. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
void LaitoneSolitaryWave (NekDouble amp, NekDouble d, NekDouble time, NekDouble x_offset)
 
void WCESolve (Array< OneD, NekDouble > &fce, NekDouble lambda)
 
void NumericalFluxForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &numfluxX, Array< OneD, NekDouble > &numfluxY)
 
void SetBoundaryConditionsForcing (Array< OneD, Array< OneD, NekDouble > > &inarray, NekDouble time)
 
void WallBoundaryForcing (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &inarray)
 
void SetBoundaryConditionsContVariables (Array< OneD, NekDouble > &inarray, NekDouble time)
 
void WallBoundaryContVariables (int bcRegion, int cnt, Array< OneD, NekDouble > &inarray)
 
void NumericalFluxConsVariables (Array< OneD, NekDouble > &physfield, Array< OneD, NekDouble > &outX, Array< OneD, NekDouble > &outY)
 
- Protected Member Functions inherited from Nektar::NonlinearSWE
 NonlinearSWE (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
void v_InitObject (bool DeclareFields=true) override
 Init object for UnsteadySystem class. More...
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
void GetFluxVector (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 
void GetVelocityVector (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
 Compute the velocity field \( \mathbf{v} \) given the momentum \( h\mathbf{v} \). More...
 
void AddVariableDepth (const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
- Protected Member Functions inherited from Nektar::ShallowWaterSystem
 ShallowWaterSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
void v_InitObject (bool DeclareFields=true) override
 Init object for UnsteadySystem class. More...
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
void SetBoundaryConditions (Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
 
void WallBoundary2D (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd)
 
void WallBoundary (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
 
void AddCoriolis (const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void PrimitiveToConservative ()
 
void ConservativeToPrimitive ()
 
NekDouble GetGravity ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals ()
 
const Array< OneD, NekDouble > & GetDepth ()
 
bool IsConstantDepth ()
 
void CopyBoundaryTrace (const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics. More...
 
SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
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. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble > > vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform dummy projection. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Private Attributes

StdRegions::ConstFactorMap m_factors
 problem type selector More...
 
ProblemType m_problemType
 
NekDouble m_const_depth
 

Friends

class MemoryManager< NonlinearPeregrine >
 

Additional Inherited Members

- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Protected Attributes inherited from Nektar::ShallowWaterSystem
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
SolverUtils::AdvectionSharedPtr m_advection
 
SolverUtils::DiffusionSharedPtr m_diffusion
 
bool m_primitive
 Indicates if variables are primitive or conservative. More...
 
bool m_constantDepth
 Indicates if constant depth case. More...
 
NekDouble m_g
 Acceleration of gravity. More...
 
Array< OneD, NekDoublem_depth
 Still water depth. More...
 
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
 
Array< OneD, NekDoublem_coriolis
 Coriolis force. More...
 
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
Array< OneD, NekDoublem_movingFrameVelsxyz
 Moving frame of reference velocities (u, v, w, omega_x, omega_y, omega_z, a_x, a_y, a_z, domega_x, domega_y, domega_z) More...
 
Array< OneD, NekDoublem_movingFrameData
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 
- Protected Attributes inherited from Nektar::SolverUtils::ALEHelper
Array< OneD, MultiRegions::ExpListSharedPtrm_fieldsALE
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 
std::vector< ALEBaseShPtrm_ALEs
 
bool m_ALESolver = false
 
bool m_ImplicitALESolver = false
 
NekDouble m_prevStageTime = 0.0
 
int m_spaceDim
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

Definition at line 53 of file NonlinearPeregrine.h.

Constructor & Destructor Documentation

◆ ~NonlinearPeregrine()

Nektar::NonlinearPeregrine::~NonlinearPeregrine ( )
overridedefault

◆ NonlinearPeregrine()

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

Definition at line 50 of file NonlinearPeregrine.cpp.

53 : NonlinearSWE(pSession, pGraph), m_factors()
54{
57 // note: eFactorTau = 1.0 becomes unstable...
58 // we need to investigate the behaviuor w.r.t. tau
59}
StdRegions::ConstFactorMap m_factors
problem type selector
NonlinearSWE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)

References Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorTau, and m_factors.

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 59 of file NonlinearPeregrine.h.

62 {
65 pGraph);
66 p->InitObject();
67 return p;
68 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

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

◆ DoOdeRhs()

void Nektar::NonlinearPeregrine::DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Definition at line 211 of file NonlinearPeregrine.cpp.

214{
215 int i;
216 int nvariables = inarray.size();
217 int ncoeffs = GetNcoeffs();
218 int nq = GetTotPoints();
219
220 switch (m_projectionType)
221 {
223 {
224
225 //-------------------------------------------------------
226 // inarray in physical space
227 Array<OneD, Array<OneD, NekDouble>> modarray(2);
228 for (i = 0; i < 2; ++i)
229 {
230 modarray[i] = Array<OneD, NekDouble>(ncoeffs, 0.0);
231 }
232 //-------------------------------------------------------
233
234 //-------------------------------------------------------
235 // Compute the DG advection including the numerical flux
236 // by using SolverUtils/Advection
237 // Input and output in physical space
238 m_advection->Advect(nvariables - 1, m_fields,
239 NullNekDoubleArrayOfArray, inarray, outarray,
240 time);
241 //-------------------------------------------------------
242
243 //-------------------------------------------------------
244 // negate the outarray since moving terms to the rhs
245 for (i = 0; i < nvariables - 1; ++i)
246 {
247 Vmath::Neg(nq, outarray[i], 1);
248 }
249 //-------------------------------------------------------
250
251 //-------------------------------------------------
252 // Add "source terms"
253 // Input and output in physical space
254
255 // Coriolis forcing
256 if (m_coriolis.size() != 0)
257 {
258 AddCoriolis(inarray, outarray);
259 }
260
261 // Variable Depth
262 if (m_constantDepth != true)
263 {
264 ASSERTL0(false,
265 "Variable depth not supported for the Peregrine "
266 "equations");
267 }
268
269 //-------------------------------------------------
270
271 //---------------------------------------
272 // As no more terms is required for the
273 // continuity equation and we have aleady evaluated
274 // the values for h_t we are done for h
275 //---------------------------------------
276
277 //-------------------------------------------------
278 // go to modal space
279 m_fields[0]->IProductWRTBase(outarray[1], modarray[0]);
280 m_fields[0]->IProductWRTBase(outarray[2], modarray[1]);
281
282 // store f1 and f2 for later use (modal space)
283 Array<OneD, NekDouble> f1(ncoeffs);
284 Array<OneD, NekDouble> f2(ncoeffs);
285
286 Vmath::Vcopy(ncoeffs, modarray[0], 1, f1, 1); // f1
287 Vmath::Vcopy(ncoeffs, modarray[1], 1, f2, 1); // f2
288
289 // Solve the remaining block-diagonal systems
290 m_fields[0]->MultiplyByElmtInvMass(modarray[0], modarray[0]);
291 m_fields[0]->MultiplyByElmtInvMass(modarray[1], modarray[1]);
292 //---------------------------------------------
293
294 //---------------------------------------------
295
296 //-------------------------------------------------
297 // create tmp fields to be used during
298 // the dispersive section
299
300 Array<OneD, Array<OneD, NekDouble>> coeffsfield(2);
301 Array<OneD, Array<OneD, NekDouble>> physfield(2);
302
303 for (i = 0; i < 2; ++i)
304 {
305 coeffsfield[i] = Array<OneD, NekDouble>(ncoeffs);
306 physfield[i] = Array<OneD, NekDouble>(nq);
307 }
308 //---------------------------------------------
309
310 //---------------------------------------------
311 // Go from modal to physical space
312 Vmath::Vcopy(nq, outarray[1], 1, physfield[0], 1);
313 Vmath::Vcopy(nq, outarray[2], 1, physfield[1], 1);
314 //---------------------------------------
315
316 //---------------------------------------
317 // Start for solve of mixed dispersive terms
318 // using the 'WCE method'
319 // (Eskilsson & Sherwin, JCP 2006)
320
321 // constant depth case
322 // \nabla \cdot (\nabla z) - invgamma z
323 // = - invgamma (\nabla \cdot {\bf f}_(2,3)
324
325 NekDouble gamma = (m_const_depth * m_const_depth) * (1.0 / 3.0);
326 NekDouble invgamma = 1.0 / gamma;
327
328 int nTraceNumPoints = GetTraceTotPoints();
329 Array<OneD, NekDouble> upwindX(nTraceNumPoints);
330 Array<OneD, NekDouble> upwindY(nTraceNumPoints);
331 //--------------------------------------------
332
333 //--------------------------------------------
334 // Compute the forcing function for the
335 // wave continuity equation
336
337 // Set boundary condidtions for z
338 SetBoundaryConditionsForcing(physfield, time);
339
340 // \nabla \phi \cdot f_{2,3}
341 m_fields[0]->IProductWRTDerivBase(0, physfield[0], coeffsfield[0]);
342 m_fields[0]->IProductWRTDerivBase(1, physfield[1], coeffsfield[1]);
343 Vmath::Vadd(ncoeffs, coeffsfield[0], 1, coeffsfield[1], 1,
344 coeffsfield[0], 1);
345 Vmath::Neg(ncoeffs, coeffsfield[0], 1);
346
347 // Evaluate upwind numerical flux (physical space)
348 NumericalFluxForcing(physfield, upwindX, upwindY);
349 Array<OneD, NekDouble> normflux(nTraceNumPoints);
350 Vmath::Vvtvvtp(nTraceNumPoints, upwindX, 1, m_traceNormals[0], 1,
351 upwindY, 1, m_traceNormals[1], 1, normflux, 1);
352 m_fields[0]->AddTraceIntegral(normflux, coeffsfield[0]);
353 m_fields[0]->MultiplyByElmtInvMass(coeffsfield[0], coeffsfield[0]);
354 m_fields[0]->BwdTrans(coeffsfield[0], physfield[0]);
355
356 Vmath::Smul(nq, -invgamma, physfield[0], 1, physfield[0], 1);
357
358 // ok: forcing function for HelmSolve... done!
359 //--------------------------------------
360
361 //--------------------------------------
362 // Solve the Helmhotz-type equation
363 // for the wave continuity equation
364 // (missing slope terms...)
365
366 // note: this is just valid for the constant depth case:
367
368 // as of now we need not to specify any
369 // BC routine for the WCE: periodic
370 // and zero Neumann (for walls)
371
372 WCESolve(physfield[0], invgamma);
373
374 Vmath::Vcopy(nq, physfield[0], 1, outarray[3], 1); // store z
375
376 // ok: Wave Continuity Equation... done!
377 //------------------------------------
378
379 //------------------------------------
380 // Return to the primary variables
381
382 // (h {\bf u})_t = gamma \nabla z + {\bf f}_{2,3}
383
384 Vmath::Smul(nq, gamma, physfield[0], 1, physfield[0], 1);
385
386 // Set boundary conditions
387 SetBoundaryConditionsContVariables(physfield[0], time);
388
389 m_fields[0]->IProductWRTDerivBase(0, physfield[0], coeffsfield[0]);
390 m_fields[0]->IProductWRTDerivBase(1, physfield[0], coeffsfield[1]);
391
392 Vmath::Neg(ncoeffs, coeffsfield[0], 1);
393 Vmath::Neg(ncoeffs, coeffsfield[1], 1);
394
395 // Evaluate upwind numerical flux (physical space)
396 NumericalFluxConsVariables(physfield[0], upwindX, upwindY);
397
398 {
399 Array<OneD, NekDouble> uptemp(nTraceNumPoints, 0.0);
400 Vmath::Vvtvvtp(nTraceNumPoints, upwindX, 1, m_traceNormals[0],
401 1, uptemp, 1, m_traceNormals[1], 1, normflux, 1);
402 m_fields[0]->AddTraceIntegral(normflux, coeffsfield[0]);
403 Vmath::Vvtvvtp(nTraceNumPoints, uptemp, 1, m_traceNormals[0], 1,
404 upwindY, 1, m_traceNormals[1], 1, normflux, 1);
405 m_fields[0]->AddTraceIntegral(normflux, coeffsfield[1]);
406 }
407
408 Vmath::Vadd(ncoeffs, f1, 1, coeffsfield[0], 1, modarray[0], 1);
409 Vmath::Vadd(ncoeffs, f2, 1, coeffsfield[1], 1, modarray[1], 1);
410
411 m_fields[1]->MultiplyByElmtInvMass(modarray[0], modarray[0]);
412 m_fields[2]->MultiplyByElmtInvMass(modarray[1], modarray[1]);
413
414 m_fields[1]->BwdTrans(modarray[0], outarray[1]);
415 m_fields[2]->BwdTrans(modarray[1], outarray[2]);
416
417 // ok: returned to conservative variables... done!
418 //---------------------
419
420 break;
421 }
424 ASSERTL0(false, "Unknown projection scheme for the Peregrine");
425 break;
426 default:
427 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
428 break;
429 }
430}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void SetBoundaryConditionsForcing(Array< OneD, Array< OneD, NekDouble > > &inarray, NekDouble time)
void NumericalFluxForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &numfluxX, Array< OneD, NekDouble > &numfluxY)
void SetBoundaryConditionsContVariables(Array< OneD, NekDouble > &inarray, NekDouble time)
void NumericalFluxConsVariables(Array< OneD, NekDouble > &physfield, Array< OneD, NekDouble > &outX, Array< OneD, NekDouble > &outY)
void WCESolve(Array< OneD, NekDouble > &fce, NekDouble lambda)
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
SolverUtils::AdvectionSharedPtr m_advection
bool m_constantDepth
Indicates if constant depth case.
Array< OneD, NekDouble > m_coriolis
Coriolis force.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT int GetTotPoints()
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
double NekDouble
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
Definition: Vmath.hpp:439
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Nektar::ShallowWaterSystem::AddCoriolis(), ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::ShallowWaterSystem::m_advection, m_const_depth, Nektar::ShallowWaterSystem::m_constantDepth, Nektar::ShallowWaterSystem::m_coriolis, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_traceNormals, Vmath::Neg(), Nektar::NullNekDoubleArrayOfArray, NumericalFluxConsVariables(), NumericalFluxForcing(), SetBoundaryConditionsContVariables(), SetBoundaryConditionsForcing(), Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vvtvvtp(), and WCESolve().

Referenced by v_InitObject().

◆ LaitoneSolitaryWave()

void Nektar::NonlinearPeregrine::LaitoneSolitaryWave ( NekDouble  amp,
NekDouble  d,
NekDouble  time,
NekDouble  x_offset 
)
protected

Definition at line 433 of file NonlinearPeregrine.cpp.

435{
436 int nq = GetTotPoints();
437
438 NekDouble A = 1.0;
439 NekDouble C = sqrt(m_g * d) * (1.0 + 0.5 * (amp / d));
440
441 Array<OneD, NekDouble> x0(nq);
442 Array<OneD, NekDouble> x1(nq);
443 Array<OneD, NekDouble> zeros(nq, 0.0);
444
445 // get the coordinates (assuming all fields have the same
446 // discretisation)
447 m_fields[0]->GetCoords(x0, x1);
448
449 for (int i = 0; i < nq; i++)
450 {
451 (m_fields[0]->UpdatePhys())[i] =
452 amp * pow((1.0 / cosh(sqrt(0.75 * (amp / (d * d * d))) *
453 (A * (x0[i] + x_offset) - C * time))),
454 2.0);
455 (m_fields[1]->UpdatePhys())[i] =
456 (amp / d) *
457 pow((1.0 / cosh(sqrt(0.75 * (amp / (d * d * d))) *
458 (A * (x0[i] + x_offset) - C * time))),
459 2.0) *
460 sqrt(m_g * d);
461 }
462
463 Vmath::Sadd(nq, d, m_fields[0]->GetPhys(), 1, m_fields[0]->UpdatePhys(), 1);
464 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
465 m_fields[1]->UpdatePhys(), 1);
466 Vmath::Vcopy(nq, zeros, 1, m_fields[2]->UpdatePhys(), 1);
467 Vmath::Vcopy(nq, zeros, 1, m_fields[3]->UpdatePhys(), 1);
468
469 // Forward transform to fill the coefficient space
470 for (int i = 0; i < 4; ++i)
471 {
472 m_fields[i]->SetPhysState(true);
473 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
474 m_fields[i]->UpdateCoeffs());
475 }
476}
NekDouble m_g
Acceleration of gravity.
std::vector< double > d(NPUPPER *NPUPPER)
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.hpp:194
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::UnitTests::d(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::ShallowWaterSystem::m_g, Vmath::Sadd(), tinysimd::sqrt(), Vmath::Vcopy(), and Vmath::Vmul().

Referenced by v_SetInitialConditions().

◆ NumericalFluxConsVariables()

void Nektar::NonlinearPeregrine::NumericalFluxConsVariables ( Array< OneD, NekDouble > &  physfield,
Array< OneD, NekDouble > &  outX,
Array< OneD, NekDouble > &  outY 
)
protected

Definition at line 698 of file NonlinearPeregrine.cpp.

701{
702 int i;
703 int nTraceNumPoints = GetTraceTotPoints();
704
705 //-----------------------------------------------------
706 // get temporary arrays
707 Array<OneD, NekDouble> Fwd(nTraceNumPoints);
708 Array<OneD, NekDouble> Bwd(nTraceNumPoints);
709 //-----------------------------------------------------
710
711 //-----------------------------------------------------
712 // get the physical values at the trace
713 // (we have put any time-dependent BC in field[1])
714
715 m_fields[1]->GetFwdBwdTracePhys(physfield, Fwd, Bwd);
716 //-----------------------------------------------------
717
718 //-----------------------------------------------------
719 // use centred fluxes for the numerical flux
720 for (i = 0; i < nTraceNumPoints; ++i)
721 {
722 outX[i] = 0.5 * (Fwd[i] + Bwd[i]);
723 outY[i] = 0.5 * (Fwd[i] + Bwd[i]);
724 }
725 //-----------------------------------------------------
726}

References Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), and Nektar::SolverUtils::EquationSystem::m_fields.

Referenced by DoOdeRhs().

◆ NumericalFluxForcing()

void Nektar::NonlinearPeregrine::NumericalFluxForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, NekDouble > &  numfluxX,
Array< OneD, NekDouble > &  numfluxY 
)
protected

Definition at line 501 of file NonlinearPeregrine.cpp.

504{
505 int i;
506 int nTraceNumPoints = GetTraceTotPoints();
507
508 //-----------------------------------------------------
509 // get temporary arrays
510 Array<OneD, Array<OneD, NekDouble>> Fwd(2);
511 Array<OneD, Array<OneD, NekDouble>> Bwd(2);
512
513 for (i = 0; i < 2; ++i)
514 {
515 Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
516 Bwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
517 }
518 //-----------------------------------------------------
519
520 //-----------------------------------------------------
521 // get the physical values at the trace
522 // (any time-dependent BC previously put in fields[1] and [2]
523
524 m_fields[1]->GetFwdBwdTracePhys(inarray[0], Fwd[0], Bwd[0]);
525 m_fields[2]->GetFwdBwdTracePhys(inarray[1], Fwd[1], Bwd[1]);
526 //-----------------------------------------------------
527
528 //-----------------------------------------------------
529 // use centred fluxes for the numerical flux
530 for (i = 0; i < nTraceNumPoints; ++i)
531 {
532 numfluxX[i] = 0.5 * (Fwd[0][i] + Bwd[0][i]);
533 numfluxY[i] = 0.5 * (Fwd[1][i] + Bwd[1][i]);
534 }
535 //-----------------------------------------------------
536}

References Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), and Nektar::SolverUtils::EquationSystem::m_fields.

Referenced by DoOdeRhs().

◆ SetBoundaryConditionsContVariables()

void Nektar::NonlinearPeregrine::SetBoundaryConditionsContVariables ( Array< OneD, NekDouble > &  inarray,
NekDouble  time 
)
protected

Definition at line 644 of file NonlinearPeregrine.cpp.

646{
647 int cnt = 0;
648
649 // loop over Boundary Regions
650 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
651 {
652 // Use wall for all
653 // Wall Boundary Condition
654 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
655 "Wall"))
656 {
657 WallBoundaryContVariables(n, cnt, inarray);
658 }
659
660 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
661 {
662 WallBoundaryContVariables(n, cnt, inarray);
663 }
664
665 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize() - 1;
666 }
667}
void WallBoundaryContVariables(int bcRegion, int cnt, Array< OneD, NekDouble > &inarray)

References Nektar::SolverUtils::EquationSystem::m_fields, and WallBoundaryContVariables().

Referenced by DoOdeRhs().

◆ SetBoundaryConditionsForcing()

void Nektar::NonlinearPeregrine::SetBoundaryConditionsForcing ( Array< OneD, Array< OneD, NekDouble > > &  inarray,
NekDouble  time 
)
protected

Definition at line 538 of file NonlinearPeregrine.cpp.

541{
542 int cnt = 0;
543
544 // loop over Boundary Regions
545 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
546 {
547 // Use wall for all BC...
548 // Wall Boundary Condition
549 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
550 "Wall"))
551 {
552 WallBoundaryForcing(n, cnt, inarray);
553 }
554
555 // Timedependent Boundary Condition
556 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
557 {
558 ASSERTL0(false, "time-dependent BC not implemented for Boussinesq");
559 }
560 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
561 }
562}
void WallBoundaryForcing(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &inarray)

References ASSERTL0, Nektar::SolverUtils::EquationSystem::m_fields, and WallBoundaryForcing().

Referenced by DoOdeRhs().

◆ v_GenerateSummary()

void Nektar::NonlinearPeregrine::v_GenerateSummary ( SolverUtils::SummaryList s)
overrideprotectedvirtual

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 163 of file NonlinearPeregrine.cpp.

164{
166 SolverUtils::AddSummaryItem(s, "", "z should be in field[3]");
167}
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47

References Nektar::SolverUtils::AddSummaryItem(), and Nektar::NonlinearSWE::v_GenerateSummary().

◆ v_InitObject()

void Nektar::NonlinearPeregrine::v_InitObject ( bool  DeclareField = true)
overrideprotectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 61 of file NonlinearPeregrine.cpp.

62{
64
65 if (m_session->DefinesSolverInfo("PROBLEMTYPE"))
66 {
67 int i;
68 std::string ProblemTypeStr = m_session->GetSolverInfo("PROBLEMTYPE");
69 for (i = 0; i < (int)SIZE_ProblemType; ++i)
70 {
71 if (boost::iequals(ProblemTypeMap[i], ProblemTypeStr))
72 {
74 break;
75 }
76 }
77 }
78 else
79 {
81 }
82
84 {
87 }
88 else
89 {
90 ASSERTL0(false, "Implicit Peregrine not set up.");
91 }
92
93 // NB! At the moment only the constant depth case is
94 // supported for the Peregrine eq.
95 if (m_session->DefinesParameter("ConstDepth"))
96 {
97 m_const_depth = m_session->GetParameter("ConstDepth");
98 }
99 else
100 {
101 ASSERTL0(false, "Constant Depth not specified");
102 }
103
104 // Type of advection class to be used
105 switch (m_projectionType)
106 {
107 // Continuous field
109 {
110 ASSERTL0(false,
111 "Continuous projection type not supported for Peregrine.");
112 break;
113 }
114 // Discontinuous field
116 {
117 std::string advName;
118 std::string diffName;
119 std::string riemName;
120
121 //---------------------------------------------------------------
122 // Setting up advection and diffusion operators
123 // NB: diffusion not set up for SWE at the moment
124 // but kept here for future use ...
125 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
126 // m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
128 advName, advName);
129
131 this);
132
133 // Setting up Riemann solver for advection operator
134 m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
135
138 riemName, m_session);
139
140 // Setting up parameters for advection operator Riemann solver
141 m_riemannSolver->SetParam("gravity",
143 m_riemannSolver->SetAuxVec("vecLocs",
146 this);
148 this);
149
150 // Concluding initialisation of advection / diffusion operators
151 m_advection->SetRiemannSolver(m_riemannSolver);
152 m_advection->InitObject(m_session, m_fields);
153 break;
154 }
155 default:
156 {
157 ASSERTL0(false, "Unsupported projection type.");
158 break;
159 }
160 }
161}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
const Array< OneD, NekDouble > & GetDepth()
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
RiemannSolverFactory & GetRiemannSolverFactory()
const char *const ProblemTypeMap[]
@ SIZE_ProblemType
Length of enum list.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), Nektar::ShallowWaterSystem::DoOdeProjection(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::SolverUtils::GetAdvectionFactory(), Nektar::ShallowWaterSystem::GetDepth(), Nektar::NonlinearSWE::GetFluxVector(), Nektar::ShallowWaterSystem::GetGravity(), Nektar::ShallowWaterSystem::GetNormals(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::ShallowWaterSystem::GetVecLocs(), Nektar::ShallowWaterSystem::m_advection, m_const_depth, Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_ode, m_problemType, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::ShallowWaterSystem::m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, Nektar::ProblemTypeMap, Nektar::SIZE_ProblemType, and Nektar::ShallowWaterSystem::v_InitObject().

◆ v_SetInitialConditions()

void Nektar::NonlinearPeregrine::v_SetInitialConditions ( NekDouble  initialtime = 0.0,
bool  dumpInitialConditions = true,
const int  domain = 0 
)
overrideprotectedvirtual

Set the initial conditions.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 172 of file NonlinearPeregrine.cpp.

175{
176 switch (m_problemType)
177 {
178 case eSolitaryWave:
179 {
180 LaitoneSolitaryWave(0.1, m_const_depth, 0.0, 0.0);
181 break;
182 }
183 default:
184 {
185 EquationSystem::v_SetInitialConditions(initialtime, false);
186 m_nchk--; // Note: m_nchk has been incremented in EquationSystem.
187 break;
188 }
189 }
190
191 if (dumpInitialConditions && m_checksteps && m_nchk == 0 &&
192 !m_comm->IsParallelInTime())
193 {
195 }
196 else if (dumpInitialConditions && m_nchk == 0 && m_comm->IsParallelInTime())
197 {
198 std::string newdir = m_sessionName + ".pit";
199 if (!fs::is_directory(newdir))
200 {
201 fs::create_directory(newdir);
202 }
203 if (m_comm->GetTimeComm()->GetRank() == 0)
204 {
205 WriteFld(newdir + "/" + m_sessionName + "_0.fld");
206 }
207 }
208 m_nchk++;
209}
void LaitoneSolitaryWave(NekDouble amp, NekDouble d, NekDouble time, NekDouble x_offset)
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
LibUtilities::CommSharedPtr m_comm
Communicator.
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
std::string m_sessionName
Name of the session.
int m_nchk
Number of checkpoints written so far.
int m_checksteps
Number of steps between checkpoints.
@ eSolitaryWave
First order Laitone solitary wave.

References Nektar::SolverUtils::EquationSystem::Checkpoint_Output(), Nektar::eSolitaryWave, LaitoneSolitaryWave(), Nektar::SolverUtils::EquationSystem::m_checksteps, Nektar::SolverUtils::EquationSystem::m_comm, m_const_depth, Nektar::SolverUtils::EquationSystem::m_nchk, m_problemType, Nektar::SolverUtils::EquationSystem::m_sessionName, Nektar::SolverUtils::EquationSystem::v_SetInitialConditions(), and Nektar::SolverUtils::EquationSystem::WriteFld().

◆ WallBoundaryContVariables()

void Nektar::NonlinearPeregrine::WallBoundaryContVariables ( int  bcRegion,
int  cnt,
Array< OneD, NekDouble > &  inarray 
)
protected

Definition at line 669 of file NonlinearPeregrine.cpp.

671{
672 int nTraceNumPoints = GetTraceTotPoints();
673
674 // get physical values of z for the forward trace
675 Array<OneD, NekDouble> z(nTraceNumPoints);
676 m_fields[0]->ExtractTracePhys(inarray, z);
677
678 // Adjust the physical values of the trace to take
679 // user defined boundaries into account
680 int e, id1, id2, npts;
682 m_fields[0]->GetBndCondExpansions()[bcRegion];
683
684 for (e = 0; e < bcexp->GetExpSize(); ++e)
685 {
686 npts = bcexp->GetExp(e)->GetTotPoints();
687 id1 = bcexp->GetPhys_Offset(e);
688 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
689 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
690
691 // copy boundary adjusted values into the boundary expansion
692 // field[1] and field[2]
693 bcexp = m_fields[1]->GetBndCondExpansions()[bcRegion];
694 Vmath::Vcopy(npts, &z[id2], 1, &(bcexp->UpdatePhys())[id1], 1);
695 }
696}
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::vector< double > z(NPUPPER)

References Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, Vmath::Vcopy(), and Nektar::UnitTests::z().

Referenced by SetBoundaryConditionsContVariables().

◆ WallBoundaryForcing()

void Nektar::NonlinearPeregrine::WallBoundaryForcing ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  inarray 
)
protected

Definition at line 565 of file NonlinearPeregrine.cpp.

567{
568 int nTraceNumPoints = GetTraceTotPoints();
569 int nvariables = 2;
570
571 // get physical values of f1 and f2 for the forward trace
572 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
573 for (int i = 0; i < nvariables; ++i)
574 {
575 Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
576 m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
577 }
578
579 // Adjust the physical values of the trace to take
580 // user defined boundaries into account
581 int e, id1, id2, npts;
583 m_fields[0]->GetBndCondExpansions()[bcRegion];
584 for (e = 0; e < bcexp->GetExpSize(); ++e)
585 {
586 npts = bcexp->GetExp(e)->GetTotPoints();
587 id1 = bcexp->GetPhys_Offset(e);
588 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
589 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
590
591 switch (m_expdim)
592 {
593 case 1:
594 {
595 ASSERTL0(false, "1D not yet implemented for Boussinesq");
596 break;
597 }
598 case 2:
599 {
600 Array<OneD, NekDouble> tmp_n(npts);
601 Array<OneD, NekDouble> tmp_t(npts);
602
603 Vmath::Vmul(npts, &Fwd[0][id2], 1, &m_traceNormals[0][id2], 1,
604 &tmp_n[0], 1);
605 Vmath::Vvtvp(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
606 &tmp_n[0], 1, &tmp_n[0], 1);
607
608 Vmath::Vmul(npts, &Fwd[0][id2], 1, &m_traceNormals[1][id2], 1,
609 &tmp_t[0], 1);
610 Vmath::Vvtvm(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
611 &tmp_t[0], 1, &tmp_t[0], 1);
612
613 // negate the normal flux
614 Vmath::Neg(npts, tmp_n, 1);
615
616 // rotate back to Cartesian
617 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
618 &Fwd[0][id2], 1);
619 Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
620 &Fwd[0][id2], 1, &Fwd[0][id2], 1);
621
622 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
623 &Fwd[1][id2], 1);
624 Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
625 &Fwd[1][id2], 1, &Fwd[1][id2], 1);
626 break;
627 }
628 case 3:
629 ASSERTL0(false, "3D not implemented for Boussinesq equations");
630 break;
631 default:
632 ASSERTL0(false, "Illegal expansion dimension");
633 }
634
635 // copy boundary adjusted values into the boundary expansion
636 bcexp = m_fields[1]->GetBndCondExpansions()[bcRegion];
637 Vmath::Vcopy(npts, &Fwd[0][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
638
639 bcexp = m_fields[2]->GetBndCondExpansions()[bcRegion];
640 Vmath::Vcopy(npts, &Fwd[1][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
641 }
642}
int m_expdim
Expansion dimension.
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Vvtvm(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)
vvtvm (vector times vector minus vector): z = w*x - y
Definition: Vmath.hpp:381

References ASSERTL0, Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::SolverUtils::EquationSystem::m_expdim, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, Vmath::Neg(), Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvm(), and Vmath::Vvtvp().

Referenced by SetBoundaryConditionsForcing().

◆ WCESolve()

void Nektar::NonlinearPeregrine::WCESolve ( Array< OneD, NekDouble > &  fce,
NekDouble  lambda 
)
protected

Definition at line 478 of file NonlinearPeregrine.cpp.

479{
480 int nq = GetTotPoints();
481
483
484 for (int j = 0; j < nq; j++)
485 {
486 (m_fields[3]->UpdatePhys())[j] = fce[j];
487 }
488
489 m_fields[3]->SetPhysState(true);
490
491 m_fields[3]->HelmSolve(m_fields[3]->GetPhys(), m_fields[3]->UpdateCoeffs(),
492 m_factors);
493
494 m_fields[3]->BwdTrans(m_fields[3]->GetCoeffs(), m_fields[3]->UpdatePhys());
495
496 m_fields[3]->SetPhysState(true);
497
498 Vmath::Vcopy(nq, m_fields[3]->GetPhys(), 1, fce, 1);
499}

References Nektar::StdRegions::eFactorLambda, Nektar::SolverUtils::EquationSystem::GetTotPoints(), m_factors, Nektar::SolverUtils::EquationSystem::m_fields, and Vmath::Vcopy().

Referenced by DoOdeRhs().

Friends And Related Function Documentation

◆ MemoryManager< NonlinearPeregrine >

friend class MemoryManager< NonlinearPeregrine >
friend

Definition at line 51 of file NonlinearPeregrine.h.

Member Data Documentation

◆ className

std::string Nektar::NonlinearPeregrine::className
static
Initial value:
=
"NonlinearPeregrine", NonlinearPeregrine::create,
"Nonlinear Peregrine equations in conservative variables.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 70 of file NonlinearPeregrine.h.

◆ m_const_depth

NekDouble Nektar::NonlinearPeregrine::m_const_depth
private

Definition at line 122 of file NonlinearPeregrine.h.

Referenced by DoOdeRhs(), v_InitObject(), and v_SetInitialConditions().

◆ m_factors

StdRegions::ConstFactorMap Nektar::NonlinearPeregrine::m_factors
private

problem type selector

Definition at line 117 of file NonlinearPeregrine.h.

Referenced by NonlinearPeregrine(), and WCESolve().

◆ m_problemType

ProblemType Nektar::NonlinearPeregrine::m_problemType
private

Definition at line 120 of file NonlinearPeregrine.h.

Referenced by v_InitObject(), and v_SetInitialConditions().