Nektar++
Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Protected Attributes | Private 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

virtual ~NonlinearPeregrine ()
 problem type selector More...
 
- Public Member Functions inherited from Nektar::ShallowWaterSystem
virtual ~ShallowWaterSystem ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 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 void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble >> &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, const Array< OneD, const NekDouble > &input)
 
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 int GetPararealIterationNumber ()
 
SOLVER_UTILS_EXPORT void SetPararealIterationNumber (int num)
 
SOLVER_UTILS_EXPORT bool GetUseInitialCondition ()
 
SOLVER_UTILS_EXPORT void SetUseInitialCondition (bool 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...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 
SOLVER_UTILS_EXPORT bool ParallelInTime ()
 Check if solver use Parallel-in-Time. More...
 

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::ShallowWaterSystem
static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Public Attributes

ProblemType m_problemType
 
- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth More...
 

Static Public Attributes

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)
 
virtual void v_InitObject (bool DeclareFields=true) override
 Init object for UnsteadySystem class. More...
 
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)
 
void GetFluxVector (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux)
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
virtual void v_PrimitiveToConservative () override
 
virtual void v_ConservativeToPrimitive () override
 
virtual void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0) override
 Set the initial conditions. More...
 
const Array< OneD, NekDouble > & GetDepthFwd ()
 
const Array< OneD, NekDouble > & GetDepthBwd ()
 
- Protected Member Functions inherited from Nektar::ShallowWaterSystem
 ShallowWaterSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
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 NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise () override
 Sets up initial conditions. 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)
 
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...
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
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 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_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

StdRegions::ConstFactorMap m_factors
 
Array< OneD, NekDoublem_dFwd
 Still water depth traces. More...
 
Array< OneD, NekDoublem_dBwd
 
- Protected Attributes inherited from Nektar::ShallowWaterSystem
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
SolverUtils::RiemannSolverSharedPtr m_riemannSolverLDG
 
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
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
NekDouble m_epsilon
 
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...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
NekDouble m_TimeIntegLambda = 0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm = -1.0
 
int m_TotLinItePerStep = 0
 
int m_StagesPerStep = 1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts = 0
 
int m_TotLinIts = 0
 
int m_TotImpStages = 0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity 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_timestepMax = -1.0
 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_pararealIter
 Number of parareal 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_useInitialCondition
 Flag to determine if IC are used. 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. More...
 
Array< OneD, NekDoublem_movingFrameTheta
 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...
 

Private Member Functions

void NumericalFlux1D (Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &numfluxX)
 
void NumericalFlux2D (Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &numfluxX, Array< OneD, Array< OneD, NekDouble >> &numfluxY)
 
void LaitoneSolitaryWave (NekDouble amp, NekDouble d, NekDouble time, NekDouble x_offset)
 
void SetBoundaryConditions (Array< OneD, Array< OneD, NekDouble >> &physarray, NekDouble time)
 
void WallBoundary2D (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble >> &Fwd, Array< OneD, Array< OneD, NekDouble >> &physarray)
 
void WallBoundary (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble >> &Fwd, Array< OneD, Array< OneD, NekDouble >> &physarray)
 Wall boundary condition. More...
 
void AddCoriolis (const Array< OneD, const Array< OneD, NekDouble >> &physarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 
void AddVariableDepth (const Array< OneD, const Array< OneD, NekDouble >> &physarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
 
void ConservativeToPrimitive (const Array< OneD, const Array< OneD, NekDouble >> &physin, Array< OneD, Array< OneD, NekDouble >> &physout)
 
void PrimitiveToConservative (const Array< OneD, const Array< OneD, NekDouble >> &physin, Array< OneD, Array< OneD, NekDouble >> &physout)
 
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 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)
 

Private Attributes

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...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

Detailed Description

Definition at line 58 of file NonlinearPeregrine.h.

Constructor & Destructor Documentation

◆ ~NonlinearPeregrine()

Nektar::NonlinearPeregrine::~NonlinearPeregrine ( )
virtual

problem type selector

Definition at line 168 of file NonlinearPeregrine.cpp.

169 {
170 }

◆ NonlinearPeregrine()

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

Definition at line 55 of file NonlinearPeregrine.cpp.

58  : ShallowWaterSystem(pSession, pGraph), m_factors()
59 {
61  m_factors[StdRegions::eFactorTau] = 1000000.0;
62  // note: eFactorTau = 1.0 becomes unstable...
63  // we need to investigate the behaviuor w.r.t. tau
64 }
StdRegions::ConstFactorMap m_factors
ShallowWaterSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

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

Member Function Documentation

◆ AddCoriolis()

void Nektar::NonlinearPeregrine::AddCoriolis ( const Array< OneD, const Array< OneD, NekDouble >> &  physarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
private

Definition at line 173 of file NonlinearPeregrine.cpp.

176 {
177 
178  int ncoeffs = GetNcoeffs();
179  int nq = GetTotPoints();
180 
181  Array<OneD, NekDouble> tmp(nq);
182  Array<OneD, NekDouble> mod(ncoeffs);
183 
184  switch (m_projectionType)
185  {
187  {
188  // add to hu equation
189  Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
190  m_fields[0]->IProductWRTBase(tmp, mod);
191  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
192  m_fields[0]->BwdTrans(mod, tmp);
193  Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
194 
195  // add to hv equation
196  Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
197  Vmath::Neg(nq, tmp, 1);
198  m_fields[0]->IProductWRTBase(tmp, mod);
199  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
200  m_fields[0]->BwdTrans(mod, tmp);
201  Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
202  break;
203  }
206  {
207  // add to hu equation
208  Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
209  Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
210 
211  // add to hv equation
212  Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
213  Vmath::Neg(nq, tmp, 1);
214  Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
215  break;
216  }
217  default:
218  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
219  break;
220  }
221 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
Array< OneD, NekDouble > m_coriolis
Coriolis force.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT int GetTotPoints()
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
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.cpp:359

References ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::ShallowWaterSystem::m_coriolis, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, Vmath::Neg(), Vmath::Vadd(), and Vmath::Vmul().

Referenced by DoOdeRhs().

◆ AddVariableDepth()

void Nektar::NonlinearPeregrine::AddVariableDepth ( const Array< OneD, const Array< OneD, NekDouble >> &  physarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray 
)
private

Definition at line 224 of file NonlinearPeregrine.cpp.

227 {
228 
229  int ncoeffs = GetNcoeffs();
230  int nq = GetTotPoints();
231 
232  Array<OneD, NekDouble> tmp(nq);
233  Array<OneD, NekDouble> mod(ncoeffs);
234 
235  switch (m_projectionType)
236  {
238  {
239  for (int i = 0; i < m_spacedim; ++i)
240  {
241  Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
242  Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
243  m_fields[0]->IProductWRTBase(tmp, mod);
244  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
245  m_fields[0]->BwdTrans(mod, tmp);
246  Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
247  }
248  break;
249  }
252  {
253  for (int i = 0; i < m_spacedim; ++i)
254  {
255  Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
256  Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
257  Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
258  }
259  break;
260  }
261  default:
262  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
263  break;
264  }
265 }
NekDouble m_g
Acceleration of gravity.
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
int m_spacedim
Spatial dimension (>= expansion dim).
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.cpp:248

References ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::ShallowWaterSystem::m_bottomSlope, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::ShallowWaterSystem::m_g, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Smul(), Vmath::Vadd(), and Vmath::Vmul().

◆ ConservativeToPrimitive()

void Nektar::NonlinearPeregrine::ConservativeToPrimitive ( const Array< OneD, const Array< OneD, NekDouble >> &  physin,
Array< OneD, Array< OneD, NekDouble >> &  physout 
)
private

Definition at line 742 of file NonlinearPeregrine.cpp.

745 {
746  int nq = GetTotPoints();
747 
748  if (physin.get() == physout.get())
749  {
750  // copy indata and work with tmp array
751  Array<OneD, Array<OneD, NekDouble>> tmp(3);
752  for (int i = 0; i < 3; ++i)
753  {
754  // deep copy
755  tmp[i] = Array<OneD, NekDouble>(nq);
756  Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
757  }
758 
759  // \eta = h - d
760  Vmath::Vsub(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
761 
762  // u = hu/h
763  Vmath::Vdiv(nq, tmp[1], 1, tmp[0], 1, physout[1], 1);
764 
765  // v = hv/ v
766  Vmath::Vdiv(nq, tmp[2], 1, tmp[0], 1, physout[2], 1);
767  }
768  else
769  {
770  // \eta = h - d
771  Vmath::Vsub(nq, physin[0], 1, m_depth, 1, physout[0], 1);
772 
773  // u = hu/h
774  Vmath::Vdiv(nq, physin[1], 1, physin[0], 1, physout[1], 1);
775 
776  // v = hv/ v
777  Vmath::Vdiv(nq, physin[2], 1, physin[0], 1, physout[2], 1);
778  }
779 }
Array< OneD, NekDouble > m_depth
Still water depth.
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.cpp:284
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:419

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::ShallowWaterSystem::m_depth, Vmath::Vcopy(), Vmath::Vdiv(), and Vmath::Vsub().

◆ 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 64 of file NonlinearPeregrine.h.

67  {
70  pGraph);
71  p->InitObject();
72  return p;
73  }
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.

◆ DoOdeProjection()

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

Definition at line 490 of file NonlinearPeregrine.cpp.

493 {
494  int i;
495  int nvariables = inarray.size();
496 
497  switch (m_projectionType)
498  {
500  {
501 
502  // Just copy over array
503  int npoints = GetNpoints();
504 
505  for (i = 0; i < nvariables; ++i)
506  {
507  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
508  }
509 
510  SetBoundaryConditions(outarray, time);
511  break;
512  }
515  {
516 
518  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
519 
520  for (i = 0; i < nvariables; ++i)
521  {
522  m_fields[i]->FwdTrans(inarray[i], coeffs);
523  m_fields[i]->BwdTrans(coeffs, outarray[i]);
524  }
525  break;
526  }
527  default:
528  ASSERTL0(false, "Unknown projection scheme");
529  break;
530  }
531 }
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble >> &physarray, NekDouble time)
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.

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

Referenced by v_InitObject().

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

270 {
271  int i;
272  int nvariables = inarray.size();
273  int ncoeffs = GetNcoeffs();
274  int nq = GetTotPoints();
275 
276  switch (m_projectionType)
277  {
279  {
280 
281  //-------------------------------------------------------
282  // inarray in physical space
283 
284  Array<OneD, Array<OneD, NekDouble>> modarray(nvariables);
285  for (i = 0; i < nvariables; ++i)
286  {
287  modarray[i] = Array<OneD, NekDouble>(ncoeffs, 0.0);
288  }
289  //-------------------------------------------------------
290 
291  //-------------------------------------------------------
292  // Compute the DG advection including the numerical flux
293  // by using SolverUtils/Advection
294  // Input and output in physical space
295  Array<OneD, Array<OneD, NekDouble>> advVel;
296 
297  m_advection->Advect(nvariables - 1, m_fields, advVel, inarray,
298  outarray, time);
299  //-------------------------------------------------------
300 
301  //-------------------------------------------------------
302  // negate the outarray since moving terms to the rhs
303  for (i = 0; i < nvariables - 1; ++i)
304  {
305  Vmath::Neg(nq, outarray[i], 1);
306  }
307  //-------------------------------------------------------
308 
309  //-------------------------------------------------
310  // Add "source terms"
311  // Input and output in physical space
312 
313  // Coriolis forcing
314  if (m_coriolis.size() != 0)
315  {
316  AddCoriolis(inarray, outarray);
317  }
318 
319  // Variable Depth
320  if (m_constantDepth != true)
321  {
322  ASSERTL0(false,
323  "Variable depth not supported for the Peregrine "
324  "equations");
325  }
326 
327  //-------------------------------------------------
328 
329  //---------------------------------------
330  // As no more terms is required for the
331  // continuity equation and we have aleady evaluated
332  // the values for h_t we are done for h
333  //---------------------------------------
334 
335  //-------------------------------------------------
336  // go to modal space
337  m_fields[0]->IProductWRTBase(outarray[1], modarray[1]);
338  m_fields[0]->IProductWRTBase(outarray[2], modarray[2]);
339 
340  // store f1 and f2 for later use (modal space)
341  Array<OneD, NekDouble> f1(ncoeffs);
342  Array<OneD, NekDouble> f2(ncoeffs);
343 
344  Vmath::Vcopy(ncoeffs, modarray[1], 1, f1, 1); // f1
345  Vmath::Vcopy(ncoeffs, modarray[2], 1, f2, 1); // f2
346 
347  // Solve the remaining block-diagonal systems
348  m_fields[0]->MultiplyByElmtInvMass(modarray[1], modarray[1]);
349  m_fields[0]->MultiplyByElmtInvMass(modarray[2], modarray[2]);
350  //---------------------------------------------
351 
352  //---------------------------------------------
353 
354  //-------------------------------------------------
355  // create tmp fields to be used during
356  // the dispersive section
357 
358  Array<OneD, Array<OneD, NekDouble>> coeffsfield(2);
359  Array<OneD, Array<OneD, NekDouble>> physfield(2);
360 
361  for (i = 0; i < 2; ++i)
362  {
363  coeffsfield[i] = Array<OneD, NekDouble>(ncoeffs);
364  physfield[i] = Array<OneD, NekDouble>(nq);
365  }
366  //---------------------------------------------
367 
368  //---------------------------------------------
369  // Go from modal to physical space
370  Vmath::Vcopy(nq, outarray[1], 1, physfield[0], 1);
371  Vmath::Vcopy(nq, outarray[2], 1, physfield[1], 1);
372  //---------------------------------------
373 
374  //---------------------------------------
375  // Start for solve of mixed dispersive terms
376  // using the 'WCE method'
377  // (Eskilsson & Sherwin, JCP 2006)
378 
379  // constant depth case
380  // \nabla \cdot (\nabla z) - invgamma z
381  // = - invgamma (\nabla \cdot {\bf f}_(2,3)
382 
383  NekDouble gamma = (m_const_depth * m_const_depth) * (1.0 / 3.0);
384  NekDouble invgamma = 1.0 / gamma;
385 
386  int nTraceNumPoints = GetTraceTotPoints();
387  Array<OneD, Array<OneD, NekDouble>> upwindX(1);
388  Array<OneD, Array<OneD, NekDouble>> upwindY(1);
389  upwindX[0] = Array<OneD, NekDouble>(nTraceNumPoints);
390  upwindY[0] = Array<OneD, NekDouble>(nTraceNumPoints);
391  //--------------------------------------------
392 
393  //--------------------------------------------
394  // Compute the forcing function for the
395  // wave continuity equation
396 
397  // Set boundary condidtions for z
398  SetBoundaryConditionsForcing(physfield, time);
399 
400  // \nabla \phi \cdot f_{2,3}
401  m_fields[0]->IProductWRTDerivBase(0, physfield[0], coeffsfield[0]);
402  m_fields[0]->IProductWRTDerivBase(1, physfield[1], coeffsfield[1]);
403  Vmath::Vadd(ncoeffs, coeffsfield[0], 1, coeffsfield[1], 1,
404  coeffsfield[0], 1);
405  Vmath::Neg(ncoeffs, coeffsfield[0], 1);
406 
407  // Evaluate upwind numerical flux (physical space)
408  NumericalFluxForcing(physfield, upwindX[0], upwindY[0]);
409 
410  m_fields[0]->AddTraceIntegral(upwindX[0], upwindY[0],
411  coeffsfield[0]);
412  m_fields[0]->MultiplyByElmtInvMass(coeffsfield[0], coeffsfield[0]);
413  m_fields[0]->BwdTrans(coeffsfield[0], physfield[0]);
414 
415  Vmath::Smul(nq, -invgamma, physfield[0], 1, physfield[0], 1);
416 
417  // ok: forcing function for HelmSolve... done!
418  //--------------------------------------
419 
420  //--------------------------------------
421  // Solve the Helmhotz-type equation
422  // for the wave continuity equation
423  // (missing slope terms...)
424 
425  // note: this is just valid for the constant depth case:
426 
427  // as of now we need not to specify any
428  // BC routine for the WCE: periodic
429  // and zero Neumann (for walls)
430 
431  WCESolve(physfield[0], invgamma);
432 
433  Vmath::Vcopy(nq, physfield[0], 1, outarray[3], 1); // store z
434 
435  // ok: Wave Continuity Equation... done!
436  //------------------------------------
437 
438  //------------------------------------
439  // Return to the primary variables
440 
441  // (h {\bf u})_t = gamma \nabla z + {\bf f}_{2,3}
442 
443  Vmath::Smul(nq, gamma, physfield[0], 1, physfield[0], 1);
444 
445  // Set boundary conditions
446  SetBoundaryConditionsContVariables(physfield[0], time);
447 
448  m_fields[0]->IProductWRTDerivBase(0, physfield[0], coeffsfield[0]);
449  m_fields[1]->IProductWRTDerivBase(1, physfield[0], coeffsfield[1]);
450 
451  Vmath::Neg(ncoeffs, coeffsfield[0], 1);
452  Vmath::Neg(ncoeffs, coeffsfield[1], 1);
453 
454  // Evaluate upwind numerical flux (physical space)
455  NumericalFluxConsVariables(physfield[0], upwindX[0], upwindY[0]);
456 
457  {
458  Array<OneD, NekDouble> uptemp(nTraceNumPoints, 0.0);
459 
460  m_fields[0]->AddTraceIntegral(upwindX[0], uptemp,
461  coeffsfield[0]);
462  m_fields[0]->AddTraceIntegral(uptemp, upwindY[0],
463  coeffsfield[1]);
464  }
465 
466  Vmath::Vadd(ncoeffs, f1, 1, coeffsfield[0], 1, modarray[1], 1);
467  Vmath::Vadd(ncoeffs, f2, 1, coeffsfield[1], 1, modarray[2], 1);
468 
469  m_fields[1]->MultiplyByElmtInvMass(modarray[1], modarray[1]);
470  m_fields[2]->MultiplyByElmtInvMass(modarray[2], modarray[2]);
471 
472  m_fields[1]->BwdTrans(modarray[1], outarray[1]);
473  m_fields[2]->BwdTrans(modarray[2], outarray[2]);
474 
475  // ok: returned to conservative variables... done!
476  //---------------------
477 
478  break;
479  }
482  ASSERTL0(false, "Unknown projection scheme for the Peregrine");
483  break;
484  default:
485  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
486  break;
487  }
488 }
void SetBoundaryConditionsContVariables(Array< OneD, NekDouble > &inarray, NekDouble time)
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble >> &physarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
void NumericalFluxConsVariables(Array< OneD, NekDouble > &physfield, Array< OneD, NekDouble > &outX, Array< OneD, NekDouble > &outY)
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 WCESolve(Array< OneD, NekDouble > &fce, NekDouble lambda)
SolverUtils::AdvectionSharedPtr m_advection
bool m_constantDepth
Indicates if constant depth case.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
double NekDouble

References 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, Vmath::Neg(), NumericalFluxConsVariables(), NumericalFluxForcing(), SetBoundaryConditionsContVariables(), SetBoundaryConditionsForcing(), Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), and WCESolve().

Referenced by v_InitObject().

◆ GetDepthBwd()

const Array<OneD, NekDouble>& Nektar::NonlinearPeregrine::GetDepthBwd ( )
inlineprotected

Definition at line 120 of file NonlinearPeregrine.h.

121  {
122  return m_dBwd;
123  }
Array< OneD, NekDouble > m_dBwd

References m_dBwd.

◆ GetDepthFwd()

const Array<OneD, NekDouble>& Nektar::NonlinearPeregrine::GetDepthFwd ( )
inlineprotected

Definition at line 116 of file NonlinearPeregrine.h.

117  {
118  return m_dFwd;
119  }
Array< OneD, NekDouble > m_dFwd
Still water depth traces.

References m_dFwd.

◆ GetFluxVector()

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

Definition at line 704 of file NonlinearPeregrine.cpp.

707 {
708  int i, j;
709  int nq = m_fields[0]->GetTotPoints();
710 
711  NekDouble g = m_g;
712  Array<OneD, Array<OneD, NekDouble>> velocity(m_spacedim);
713 
714  // Flux vector for the mass equation
715  for (i = 0; i < m_spacedim; ++i)
716  {
717  velocity[i] = Array<OneD, NekDouble>(nq);
718  Vmath::Vcopy(nq, physfield[i + 1], 1, flux[0][i], 1);
719  }
720 
721  GetVelocityVector(physfield, velocity);
722 
723  // Put (0.5 g h h) in tmp
724  Array<OneD, NekDouble> tmp(nq);
725  Vmath::Vmul(nq, physfield[0], 1, physfield[0], 1, tmp, 1);
726  Vmath::Smul(nq, 0.5 * g, tmp, 1, tmp, 1);
727 
728  // Flux vector for the momentum equations
729  for (i = 0; i < m_spacedim; ++i)
730  {
731  for (j = 0; j < m_spacedim; ++j)
732  {
733  Vmath::Vmul(nq, velocity[j], 1, physfield[i + 1], 1, flux[i + 1][j],
734  1);
735  }
736 
737  // Add (0.5 g h h) to appropriate field
738  Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
739  }
740 }
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity)
Compute the velocity field given the momentum .

References GetVelocityVector(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::ShallowWaterSystem::m_g, Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vmul().

Referenced by v_InitObject().

◆ GetVelocityVector()

void Nektar::NonlinearPeregrine::GetVelocityVector ( const Array< OneD, Array< OneD, NekDouble >> &  physfield,
Array< OneD, Array< OneD, NekDouble >> &  velocity 
)
private

Compute the velocity field \( \mathbf{v} \) given the momentum \( h\mathbf{v} \).

Parameters
physfieldMomentum field.
velocityVelocity field.

Definition at line 862 of file NonlinearPeregrine.cpp.

865 {
866  const int npts = physfield[0].size();
867 
868  for (int i = 0; i < m_spacedim; ++i)
869  {
870  Vmath::Vdiv(npts, physfield[1 + i], 1, physfield[0], 1, velocity[i], 1);
871  }
872 }

References Nektar::SolverUtils::EquationSystem::m_spacedim, and Vmath::Vdiv().

Referenced by GetFluxVector().

◆ LaitoneSolitaryWave()

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

Definition at line 1143 of file NonlinearPeregrine.cpp.

1145 {
1146  int nq = GetTotPoints();
1147 
1148  NekDouble A = 1.0;
1149  NekDouble C = sqrt(m_g * d) * (1.0 + 0.5 * (amp / d));
1150 
1151  Array<OneD, NekDouble> x0(nq);
1152  Array<OneD, NekDouble> x1(nq);
1153  Array<OneD, NekDouble> zeros(nq, 0.0);
1154 
1155  // get the coordinates (assuming all fields have the same
1156  // discretisation)
1157  m_fields[0]->GetCoords(x0, x1);
1158 
1159  for (int i = 0; i < nq; i++)
1160  {
1161  (m_fields[0]->UpdatePhys())[i] =
1162  amp * pow((1.0 / cosh(sqrt(0.75 * (amp / (d * d * d))) *
1163  (A * (x0[i] + x_offset) - C * time))),
1164  2.0);
1165  (m_fields[1]->UpdatePhys())[i] =
1166  (amp / d) *
1167  pow((1.0 / cosh(sqrt(0.75 * (amp / (d * d * d))) *
1168  (A * (x0[i] + x_offset) - C * time))),
1169  2.0) *
1170  sqrt(m_g * d);
1171  }
1172 
1173  Vmath::Sadd(nq, d, m_fields[0]->GetPhys(), 1, m_fields[0]->UpdatePhys(), 1);
1174  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
1175  m_fields[1]->UpdatePhys(), 1);
1176  Vmath::Vcopy(nq, zeros, 1, m_fields[2]->UpdatePhys(), 1);
1177  Vmath::Vcopy(nq, zeros, 1, m_fields[3]->UpdatePhys(), 1);
1178 
1179  // Forward transform to fill the coefficient space
1180  for (int i = 0; i < 4; ++i)
1181  {
1182  m_fields[i]->SetPhysState(true);
1183  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1184  m_fields[i]->UpdateCoeffs());
1185  }
1186 }
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:384
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

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

◆ NumericalFlux1D()

void Nektar::NonlinearPeregrine::NumericalFlux1D ( Array< OneD, Array< OneD, NekDouble >> &  physfield,
Array< OneD, Array< OneD, NekDouble >> &  numfluxX 
)
private

◆ NumericalFlux2D()

void Nektar::NonlinearPeregrine::NumericalFlux2D ( Array< OneD, Array< OneD, NekDouble >> &  physfield,
Array< OneD, Array< OneD, NekDouble >> &  numfluxX,
Array< OneD, Array< OneD, NekDouble >> &  numfluxY 
)
private

◆ NumericalFluxConsVariables()

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

Definition at line 1109 of file NonlinearPeregrine.cpp.

1112 {
1113  int i;
1114  int nTraceNumPoints = GetTraceTotPoints();
1115 
1116  //-----------------------------------------------------
1117  // get temporary arrays
1118  Array<OneD, Array<OneD, NekDouble>> Fwd(1);
1119  Array<OneD, Array<OneD, NekDouble>> Bwd(1);
1120 
1121  Fwd[0] = Array<OneD, NekDouble>(nTraceNumPoints);
1122  Bwd[0] = Array<OneD, NekDouble>(nTraceNumPoints);
1123  //-----------------------------------------------------
1124 
1125  //-----------------------------------------------------
1126  // get the physical values at the trace
1127  // (we have put any time-dependent BC in field[1])
1128 
1129  m_fields[1]->GetFwdBwdTracePhys(physfield, Fwd[0], Bwd[0]);
1130  //-----------------------------------------------------
1131 
1132  //-----------------------------------------------------
1133  // use centred fluxes for the numerical flux
1134  for (i = 0; i < nTraceNumPoints; ++i)
1135  {
1136  outX[i] = 0.5 * (Fwd[0][i] + Bwd[0][i]);
1137  outY[i] = 0.5 * (Fwd[0][i] + Bwd[0][i]);
1138  }
1139  //-----------------------------------------------------
1140 }

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 
)
private

Definition at line 906 of file NonlinearPeregrine.cpp.

909 {
910  int i;
911  int nTraceNumPoints = GetTraceTotPoints();
912 
913  //-----------------------------------------------------
914  // get temporary arrays
915  Array<OneD, Array<OneD, NekDouble>> Fwd(2);
916  Array<OneD, Array<OneD, NekDouble>> Bwd(2);
917 
918  for (i = 0; i < 2; ++i)
919  {
920  Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
921  Bwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
922  }
923  //-----------------------------------------------------
924 
925  //-----------------------------------------------------
926  // get the physical values at the trace
927  // (any time-dependent BC previuosly put in fields[1] and [2]
928 
929  m_fields[1]->GetFwdBwdTracePhys(inarray[0], Fwd[0], Bwd[0]);
930  m_fields[2]->GetFwdBwdTracePhys(inarray[1], Fwd[1], Bwd[1]);
931  //-----------------------------------------------------
932 
933  //-----------------------------------------------------
934  // use centred fluxes for the numerical flux
935  for (i = 0; i < nTraceNumPoints; ++i)
936  {
937  numfluxX[i] = 0.5 * (Fwd[0][i] + Bwd[0][i]);
938  numfluxY[i] = 0.5 * (Fwd[1][i] + Bwd[1][i]);
939  }
940  //-----------------------------------------------------
941 }

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

Referenced by DoOdeRhs().

◆ PrimitiveToConservative()

void Nektar::NonlinearPeregrine::PrimitiveToConservative ( const Array< OneD, const Array< OneD, NekDouble >> &  physin,
Array< OneD, Array< OneD, NekDouble >> &  physout 
)
private

Definition at line 798 of file NonlinearPeregrine.cpp.

801 {
802 
803  int nq = GetTotPoints();
804 
805  if (physin.get() == physout.get())
806  {
807  // copy indata and work with tmp array
808  Array<OneD, Array<OneD, NekDouble>> tmp(3);
809  for (int i = 0; i < 3; ++i)
810  {
811  // deep copy
812  tmp[i] = Array<OneD, NekDouble>(nq);
813  Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
814  }
815 
816  // h = \eta + d
817  Vmath::Vadd(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
818 
819  // hu = h * u
820  Vmath::Vmul(nq, physout[0], 1, tmp[1], 1, physout[1], 1);
821 
822  // hv = h * v
823  Vmath::Vmul(nq, physout[0], 1, tmp[2], 1, physout[2], 1);
824  }
825  else
826  {
827  // h = \eta + d
828  Vmath::Vadd(nq, physin[0], 1, m_depth, 1, physout[0], 1);
829 
830  // hu = h * u
831  Vmath::Vmul(nq, physout[0], 1, physin[1], 1, physout[1], 1);
832 
833  // hv = h * v
834  Vmath::Vmul(nq, physout[0], 1, physin[2], 1, physout[2], 1);
835  }
836 }

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::ShallowWaterSystem::m_depth, Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Vmul().

◆ SetBoundaryConditions()

void Nektar::NonlinearPeregrine::SetBoundaryConditions ( Array< OneD, Array< OneD, NekDouble >> &  physarray,
NekDouble  time 
)
private

Definition at line 534 of file NonlinearPeregrine.cpp.

536 {
537 
538  int nvariables = m_fields.size();
539  int cnt = 0;
540  int nTracePts = GetTraceTotPoints();
541 
542  // Extract trace for boundaries. Needs to be done on all processors to avoid
543  // deadlock.
544  Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
545  for (int i = 0; i < nvariables; ++i)
546  {
547  Fwd[i] = Array<OneD, NekDouble>(nTracePts);
548  m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
549  }
550 
551  // loop over Boundary Regions
552  for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
553  {
554 
555  // Wall Boundary Condition
556  if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
557  "Wall"))
558  {
559  WallBoundary2D(n, cnt, Fwd, inarray);
560  }
561 
562  // Time Dependent Boundary Condition (specified in meshfile)
563  if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
564  {
565  for (int i = 0; i < nvariables; ++i)
566  {
567  m_fields[i]->EvaluateBoundaryConditions(time);
568  }
569  }
570  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
571  }
572 }
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble >> &Fwd, Array< OneD, Array< OneD, NekDouble >> &physarray)

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

Referenced by DoOdeProjection().

◆ SetBoundaryConditionsContVariables()

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

Definition at line 1053 of file NonlinearPeregrine.cpp.

1055 {
1056  boost::ignore_unused(time);
1057 
1058  int cnt = 0;
1059 
1060  // loop over Boundary Regions
1061  for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
1062  {
1063  // Use wall for all
1064  // Wall Boundary Condition
1065  if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
1066  "Wall"))
1067  {
1068  WallBoundaryContVariables(n, cnt, inarray);
1069  }
1070 
1071  if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
1072  {
1073  WallBoundaryContVariables(n, cnt, inarray);
1074  }
1075 
1076  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize() - 1;
1077  }
1078 }
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 
)
private

Definition at line 943 of file NonlinearPeregrine.cpp.

945 {
946  boost::ignore_unused(time);
947 
948  int cnt = 0;
949 
950  // loop over Boundary Regions
951  for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
952  {
953  // Use wall for all BC...
954  // Wall Boundary Condition
955  if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
956  "Wall"))
957  {
958  WallBoundaryForcing(n, cnt, inarray);
959  }
960 
961  // Timedependent Boundary Condition
962  if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
963  {
964  ASSERTL0(false, "time-dependent BC not implemented for Boussinesq");
965  }
966  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
967  }
968 }
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_ConservativeToPrimitive()

void Nektar::NonlinearPeregrine::v_ConservativeToPrimitive ( )
overrideprotectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 781 of file NonlinearPeregrine.cpp.

782 {
783  int nq = GetTotPoints();
784 
785  // u = hu/h
786  Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
787  m_fields[1]->UpdatePhys(), 1);
788 
789  // v = hv/ v
790  Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
791  m_fields[2]->UpdatePhys(), 1);
792 
793  // \eta = h - d
794  Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
795  m_fields[0]->UpdatePhys(), 1);
796 }

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::ShallowWaterSystem::m_depth, Nektar::SolverUtils::EquationSystem::m_fields, Vmath::Vdiv(), and Vmath::Vsub().

◆ v_GenerateSummary()

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

Print a summary of time stepping parameters.

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 874 of file NonlinearPeregrine.cpp.

875 {
877  SolverUtils::AddSummaryItem(s, "Variables", "h should be in field[0]");
878  SolverUtils::AddSummaryItem(s, "", "hu should be in field[1]");
879  SolverUtils::AddSummaryItem(s, "", "hv should be in field[2]");
880  SolverUtils::AddSummaryItem(s, "", "z should be in field[3]");
881 }
virtual 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:49

References Nektar::SolverUtils::AddSummaryItem(), and Nektar::ShallowWaterSystem::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::ShallowWaterSystem.

Definition at line 66 of file NonlinearPeregrine.cpp.

67 {
68  ShallowWaterSystem::v_InitObject(DeclareFields);
69 
70  if (m_session->DefinesSolverInfo("PROBLEMTYPE"))
71  {
72  int i;
73  std::string ProblemTypeStr = m_session->GetSolverInfo("PROBLEMTYPE");
74  for (i = 0; i < (int)SIZE_ProblemType; ++i)
75  {
76  if (boost::iequals(ProblemTypeMap[i], ProblemTypeStr))
77  {
79  break;
80  }
81  }
82  }
83  else
84  {
86  }
87 
89  {
92  }
93  else
94  {
95  ASSERTL0(false, "Implicit Peregrine not set up.");
96  }
97 
98  // NB! At the moment only the constant depth case is
99  // supported for the Peregrine eq.
100  if (m_session->DefinesParameter("ConstDepth"))
101  {
102  m_const_depth = m_session->GetParameter("ConstDepth");
103  }
104  else
105  {
106  ASSERTL0(false, "Constant Depth not specified");
107  }
108 
109  // Type of advection class to be used
110  switch (m_projectionType)
111  {
112  // Continuous field
114  {
115  ASSERTL0(false,
116  "Continuous projection type not supported for Peregrine.");
117  break;
118  }
119  // Discontinuous field
121  {
122  string advName;
123  string diffName;
124  string riemName;
125 
126  //---------------------------------------------------------------
127  // Setting up advection and diffusion operators
128  // NB: diffusion not set up for SWE at the moment
129  // but kept here for future use ...
130  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
131  // m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
133  advName, advName);
134 
136  this);
137 
138  // Setting up Riemann solver for advection operator
139  m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
140 
143  riemName, m_session);
144 
145  // Setting up parameters for advection operator Riemann solver
146  m_riemannSolver->SetParam("gravity",
148  m_riemannSolver->SetAuxVec("vecLocs",
151  this);
152  m_riemannSolver->SetScalar("depth", &NonlinearPeregrine::GetDepth,
153  this);
154 
155  // Concluding initialisation of advection / diffusion operators
156  m_advection->SetRiemannSolver(m_riemannSolver);
157  m_advection->InitObject(m_session, m_fields);
158  break;
159  }
160  default:
161  {
162  ASSERTL0(false, "Unsupported projection type.");
163  break;
164  }
165  }
166 }
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
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()
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
virtual void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
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:47
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(), DoOdeProjection(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::SolverUtils::GetAdvectionFactory(), Nektar::ShallowWaterSystem::GetDepth(), 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_PrimitiveToConservative()

void Nektar::NonlinearPeregrine::v_PrimitiveToConservative ( )
overrideprotectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 838 of file NonlinearPeregrine.cpp.

839 {
840  int nq = GetTotPoints();
841 
842  // h = \eta + d
843  Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
844  m_fields[0]->UpdatePhys(), 1);
845 
846  // hu = h * u
847  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
848  m_fields[1]->UpdatePhys(), 1);
849 
850  // hv = h * v
851  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
852  m_fields[2]->UpdatePhys(), 1);
853 }

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::ShallowWaterSystem::m_depth, Nektar::SolverUtils::EquationSystem::m_fields, Vmath::Vadd(), and Vmath::Vmul().

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

1194 {
1195  boost::ignore_unused(domain);
1196 
1197  switch (m_problemType)
1198  {
1199  case eSolitaryWave:
1200  {
1201  LaitoneSolitaryWave(0.1, m_const_depth, 0.0, 0.0);
1202  break;
1203  }
1204  default:
1205  {
1206  EquationSystem::v_SetInitialConditions(initialtime, false);
1207  break;
1208  }
1209  }
1210 
1211  if (dumpInitialConditions)
1212  {
1213  // Dump initial conditions to file
1214  Checkpoint_Output(0);
1215  }
1216 }
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)
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
@ eSolitaryWave
First order Laitone solitary wave.

References Nektar::SolverUtils::EquationSystem::Checkpoint_Output(), Nektar::eSolitaryWave, LaitoneSolitaryWave(), m_const_depth, m_problemType, and Nektar::SolverUtils::EquationSystem::v_SetInitialConditions().

◆ WallBoundary()

void Nektar::NonlinearPeregrine::WallBoundary ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble >> &  Fwd,
Array< OneD, Array< OneD, NekDouble >> &  physarray 
)
private

Wall boundary condition.

Definition at line 578 of file NonlinearPeregrine.cpp.

581 {
582  int i;
583  int nvariables = physarray.size();
584 
585  // Adjust the physical values of the trace to take
586  // user defined boundaries into account
587  int e, id1, id2, npts;
589  m_fields[0]->GetBndCondExpansions()[bcRegion];
590  for (e = 0; e < bcexp->GetExpSize(); ++e)
591  {
592  npts = bcexp->GetExp(e)->GetTotPoints();
593  id1 = bcexp->GetPhys_Offset(e);
594  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
595  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
596 
597  // For 2D/3D, define: v* = v - 2(v.n)n
598  Array<OneD, NekDouble> tmp(npts, 0.0);
599 
600  // Calculate (v.n)
601  for (i = 0; i < m_spacedim; ++i)
602  {
603  Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
604  &tmp[0], 1, &tmp[0], 1);
605  }
606 
607  // Calculate 2.0(v.n)
608  Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
609 
610  // Calculate v* = v - 2.0(v.n)n
611  for (i = 0; i < m_spacedim; ++i)
612  {
613  Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
614  &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
615  }
616 
617  // copy boundary adjusted values into the boundary expansion
618  for (i = 0; i < nvariables; ++i)
619  {
620  bcexp = m_fields[i]->GetBndCondExpansions()[bcRegion];
621  Vmath::Vcopy(npts, &Fwd[i][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
622  }
623  }
624 }
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::EquationSystem::m_traceNormals, Vmath::Smul(), Vmath::Vcopy(), and Vmath::Vvtvp().

◆ WallBoundary2D()

void Nektar::NonlinearPeregrine::WallBoundary2D ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble >> &  Fwd,
Array< OneD, Array< OneD, NekDouble >> &  physarray 
)
private

Definition at line 626 of file NonlinearPeregrine.cpp.

629 {
630  boost::ignore_unused(physarray);
631 
632  int i;
633  int nvariables = 3;
634 
635  // Adjust the physical values of the trace to take
636  // user defined boundaries into account
637  int e, id1, id2, npts;
639  m_fields[0]->GetBndCondExpansions()[bcRegion];
640 
641  for (e = 0; e < bcexp->GetExpSize(); ++e)
642  {
643  npts = bcexp->GetExp(e)->GetNumPoints(0);
644  id1 = bcexp->GetPhys_Offset(e);
645  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
646  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
647 
648  switch (m_expdim)
649  {
650  case 1:
651  {
652  // negate the forward flux
653  Vmath::Neg(npts, &Fwd[1][id2], 1);
654  break;
655  }
656  case 2:
657  {
658  Array<OneD, NekDouble> tmp_n(npts);
659  Array<OneD, NekDouble> tmp_t(npts);
660 
661  Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
662  &tmp_n[0], 1);
663  Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
664  &tmp_n[0], 1, &tmp_n[0], 1);
665 
666  Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
667  &tmp_t[0], 1);
668  Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
669  &tmp_t[0], 1, &tmp_t[0], 1);
670 
671  // negate the normal flux
672  Vmath::Neg(npts, tmp_n, 1);
673 
674  // rotate back to Cartesian
675  Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
676  &Fwd[1][id2], 1);
677  Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
678  &Fwd[1][id2], 1, &Fwd[1][id2], 1);
679 
680  Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
681  &Fwd[2][id2], 1);
682  Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
683  &Fwd[2][id2], 1, &Fwd[2][id2], 1);
684  break;
685  }
686  case 3:
687  ASSERTL0(false,
688  "3D not implemented for Shallow Water Equations");
689  break;
690  default:
691  ASSERTL0(false, "Illegal expansion dimension");
692  }
693 
694  // copy boundary adjusted values into the boundary expansion
695  for (i = 0; i < nvariables; ++i)
696  {
697  bcexp = m_fields[i]->GetBndCondExpansions()[bcRegion];
698  Vmath::Vcopy(npts, &Fwd[i][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
699  }
700  }
701 }
int m_expdim
Expansion dimension.
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.cpp:598

References ASSERTL0, 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 SetBoundaryConditions().

◆ WallBoundaryContVariables()

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

Definition at line 1080 of file NonlinearPeregrine.cpp.

1082 {
1083  int nTraceNumPoints = GetTraceTotPoints();
1084 
1085  // get physical values of z for the forward trace
1086  Array<OneD, NekDouble> z(nTraceNumPoints);
1087  m_fields[0]->ExtractTracePhys(inarray, z);
1088 
1089  // Adjust the physical values of the trace to take
1090  // user defined boundaries into account
1091  int e, id1, id2, npts;
1093  m_fields[0]->GetBndCondExpansions()[bcRegion];
1094 
1095  for (e = 0; e < bcexp->GetExpSize(); ++e)
1096  {
1097  npts = bcexp->GetExp(e)->GetTotPoints();
1098  id1 = bcexp->GetPhys_Offset(e);
1099  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
1100  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
1101 
1102  // copy boundary adjusted values into the boundary expansion
1103  // field[1] and field[2]
1104  bcexp = m_fields[1]->GetBndCondExpansions()[bcRegion];
1105  Vmath::Vcopy(npts, &z[id2], 1, &(bcexp->UpdatePhys())[id1], 1);
1106  }
1107 }

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

Referenced by SetBoundaryConditionsContVariables().

◆ WallBoundaryForcing()

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

Definition at line 971 of file NonlinearPeregrine.cpp.

973 {
974 
975  // std::cout << " WallBoundaryForcing" << std::endl;
976 
977  int nTraceNumPoints = GetTraceTotPoints();
978  int nvariables = 2;
979 
980  // get physical values of f1 and f2 for the forward trace
981  Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
982  for (int i = 0; i < nvariables; ++i)
983  {
984  Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
985  m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
986  }
987 
988  // Adjust the physical values of the trace to take
989  // user defined boundaries into account
990  int e, id1, id2, npts;
992  m_fields[0]->GetBndCondExpansions()[bcRegion];
993  for (e = 0; e < bcexp->GetExpSize(); ++e)
994  {
995  npts = bcexp->GetExp(e)->GetTotPoints();
996  id1 = bcexp->GetPhys_Offset(e);
997  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
998  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
999 
1000  switch (m_expdim)
1001  {
1002  case 1:
1003  {
1004  ASSERTL0(false, "1D not yet implemented for Boussinesq");
1005  break;
1006  }
1007  case 2:
1008  {
1009  Array<OneD, NekDouble> tmp_n(npts);
1010  Array<OneD, NekDouble> tmp_t(npts);
1011 
1012  Vmath::Vmul(npts, &Fwd[0][id2], 1, &m_traceNormals[0][id2], 1,
1013  &tmp_n[0], 1);
1014  Vmath::Vvtvp(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
1015  &tmp_n[0], 1, &tmp_n[0], 1);
1016 
1017  Vmath::Vmul(npts, &Fwd[0][id2], 1, &m_traceNormals[1][id2], 1,
1018  &tmp_t[0], 1);
1019  Vmath::Vvtvm(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
1020  &tmp_t[0], 1, &tmp_t[0], 1);
1021 
1022  // negate the normal flux
1023  Vmath::Neg(npts, tmp_n, 1);
1024 
1025  // rotate back to Cartesian
1026  Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
1027  &Fwd[0][id2], 1);
1028  Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
1029  &Fwd[0][id2], 1, &Fwd[0][id2], 1);
1030 
1031  Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
1032  &Fwd[1][id2], 1);
1033  Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
1034  &Fwd[1][id2], 1, &Fwd[1][id2], 1);
1035  break;
1036  }
1037  case 3:
1038  ASSERTL0(false, "3D not implemented for Boussinesq equations");
1039  break;
1040  default:
1041  ASSERTL0(false, "Illegal expansion dimension");
1042  }
1043 
1044  // copy boundary adjusted values into the boundary expansion
1045  bcexp = m_fields[1]->GetBndCondExpansions()[bcRegion];
1046  Vmath::Vcopy(npts, &Fwd[0][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
1047 
1048  bcexp = m_fields[2]->GetBndCondExpansions()[bcRegion];
1049  Vmath::Vcopy(npts, &Fwd[1][id2], 1, &(bcexp->UpdatePhys())[id1], 1);
1050  }
1051 }

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 
)
private

Definition at line 883 of file NonlinearPeregrine.cpp.

884 {
885  int nq = GetTotPoints();
886 
888 
889  for (int j = 0; j < nq; j++)
890  {
891  (m_fields[3]->UpdatePhys())[j] = fce[j];
892  }
893 
894  m_fields[3]->SetPhysState(true);
895 
896  m_fields[3]->HelmSolve(m_fields[3]->GetPhys(), m_fields[3]->UpdateCoeffs(),
897  m_factors);
898 
899  m_fields[3]->BwdTrans(m_fields[3]->GetCoeffs(), m_fields[3]->UpdatePhys());
900 
901  m_fields[3]->SetPhysState(true);
902 
903  Vmath::Vcopy(nq, m_fields[3]->GetPhys(), 1, fce, 1);
904 }

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

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.
Definition: NekFactory.hpp:198
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 75 of file NonlinearPeregrine.h.

◆ m_const_depth

NekDouble Nektar::NonlinearPeregrine::m_const_depth
private

Definition at line 126 of file NonlinearPeregrine.h.

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

◆ m_dBwd

Array<OneD, NekDouble> Nektar::NonlinearPeregrine::m_dBwd
protected

Definition at line 92 of file NonlinearPeregrine.h.

Referenced by GetDepthBwd().

◆ m_dFwd

Array<OneD, NekDouble> Nektar::NonlinearPeregrine::m_dFwd
protected

Still water depth traces.

Definition at line 91 of file NonlinearPeregrine.h.

Referenced by GetDepthFwd().

◆ m_factors

StdRegions::ConstFactorMap Nektar::NonlinearPeregrine::m_factors
protected

Definition at line 83 of file NonlinearPeregrine.h.

Referenced by NonlinearPeregrine(), and WCESolve().

◆ m_problemType

ProblemType Nektar::NonlinearPeregrine::m_problemType

Definition at line 80 of file NonlinearPeregrine.h.

Referenced by v_InitObject(), and v_SetInitialConditions().