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

#include <LinearSWE.h>

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

Public Member Functions

virtual ~LinearSWE ()
- Public Member Functions inherited from Nektar::ShallowWaterSystem
virtual ~ShallowWaterSystem ()
 Destructor.
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable.
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor.
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object.
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem.
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem.
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space.
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space.
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve.
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation.
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name.
template<class T >
boost::shared_ptr< T > as ()
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name.
SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr 
GetSession ()
 Get Session name.
SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
GetPressure ()
 Get pressure field if available.
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics.
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda.
SOLVER_UTILS_EXPORT void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name.
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow.
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields.
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution.
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution.
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields.
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf].
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $.
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection.
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection.
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form.
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields.
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields.
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields.
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename.
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file.
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file.
SOLVER_UTILS_EXPORT void ScanForHistoryPoints ()
 Builds map of which element holds each history point.
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to file.
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary.
SOLVER_UTILS_EXPORT Array
< OneD,
MultiRegions::ExpListSharedPtr > & 
UpdateFields ()
SOLVER_UTILS_EXPORT
LibUtilities::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated.
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time.
SOLVER_UTILS_EXPORT int GetNcoeffs ()
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
SOLVER_UTILS_EXPORT int GetNumExpModes ()
SOLVER_UTILS_EXPORT const
Array< OneD, int > 
GetNumExpModesPerExp ()
SOLVER_UTILS_EXPORT int GetNvariables ()
SOLVER_UTILS_EXPORT const
std::string 
GetVariable (unsigned int i)
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
SOLVER_UTILS_EXPORT int GetExpSize ()
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
SOLVER_UTILS_EXPORT int GetTotPoints ()
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
SOLVER_UTILS_EXPORT int GetNpoints ()
SOLVER_UTILS_EXPORT int GetNumElmVelocity ()
SOLVER_UTILS_EXPORT int GetSteps ()
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, Array< OneD, NekDouble > &output)
SOLVER_UTILS_EXPORT void SetStepsToOne ()
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
SOLVER_UTILS_EXPORT void FwdTransFields ()
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const string &s1, const string &s2)
 Perform a case-insensitive string comparison.

Static Public Member Functions

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

Static Public Attributes

static std::string className
 Name of class.
- Static Public Attributes inherited from Nektar::ShallowWaterSystem
static std::string className
 Name of class.

Protected Member Functions

 LinearSWE (const LibUtilities::SessionReaderSharedPtr &pSession)
virtual void v_InitObject ()
 Init object for UnsteadySystem class.
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)
 Print a summary of time stepping parameters.
virtual void v_PrimitiveToConservative ()
virtual void v_ConservativeToPrimitive ()
const Array< OneD, NekDouble > & GetDepthFwd ()
const Array< OneD, NekDouble > & GetDepthBwd ()
- Protected Member Functions inherited from Nektar::ShallowWaterSystem
 ShallowWaterSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members.
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)
 Initialises UnsteadySystem class members.
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control.
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file.
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Return the timestep to be used for the next step in the time-marching loop.
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members.
int nocase_cmp (const string &s1, const string &s2)
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time.
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution.
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution.
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space.
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space.
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
virtual SOLVER_UTILS_EXPORT void v_Output (void)
virtual SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
v_GetPressure (void)
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)

Protected Attributes

Array< OneD, NekDoublem_dFwd
 Still water depth traces.
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.
bool m_constantDepth
 Indicates if constant depth case.
NekDouble m_g
 Acceleration of gravity.
Array< OneD, NekDoublem_depth
 Still water depth.
Array< OneD, Array< OneD,
NekDouble > > 
m_bottomSlope
Array< OneD, NekDoublem_coriolis
 Coriolis force.
Array< OneD, Array< OneD,
NekDouble > > 
m_vecLocs
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status information.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme.
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use.
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
NekDouble m_epsilon
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used.
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used.
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used.
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state.
std::vector< int > m_intVariables
std::vector< FilterSharedPtrm_filters
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator.
LibUtilities::SessionReaderSharedPtr m_session
 The session reader.
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 Array holding all dependent variables.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object.
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh.
std::string m_sessionName
 Name of the session.
NekDouble m_time
 Current time of simulation.
NekDouble m_fintime
 Finish time of the simulation.
NekDouble m_timestep
 Time step size.
NekDouble m_lambda
 Lambda constant in real system if one required.
NekDouble m_checktime
 Time between checkpoints.
int m_steps
 Number of steps to take.
int m_checksteps
 Number of steps between checkpoints.
int m_spacedim
 Spatial dimension (>= expansion dim).
int m_expdim
 Expansion dimension.
bool m_SingleMode
 Flag to determine if single homogeneous mode is used.
bool m_HalfMode
 Flag to determine if half homogeneous mode is used.
bool m_MultipleModes
 Flag to determine if use multiple homogenenous modes are used.
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform.
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations.
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous.
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields.
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error.
enum HomogeneousType m_HomogeneousType
NekDouble m_LhomX
 physical length in X direction (if homogeneous)
NekDouble m_LhomY
 physical length in Y direction (if homogeneous)
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous)
int m_npointsX
 number of points in X direction (if homogeneous)
int m_npointsY
 number of points in Y direction (if homogeneous)
int m_npointsZ
 number of points in Z direction (if homogeneous)
int m_HomoDirec
 number of homogenous directions
int m_NumMode
 Mode to use in case of single mode analysis.

Private Member Functions

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 SetBoundaryConditions (Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
void WallBoundary2D (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &physarray)
void WallBoundary (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &physarray)
 Wall boundary condition.
void AddCoriolis (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} $.

Friends

class MemoryManager< LinearSWE >

Additional Inherited Members

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

Detailed Description

Definition at line 50 of file LinearSWE.h.

Constructor & Destructor Documentation

Nektar::LinearSWE::~LinearSWE ( )
virtual

Definition at line 166 of file LinearSWE.cpp.

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

Definition at line 50 of file LinearSWE.cpp.

: ShallowWaterSystem(pSession)
{
}

Member Function Documentation

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

Definition at line 172 of file LinearSWE.cpp.

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

{
int ncoeffs = GetNcoeffs();
int nq = GetTotPoints();
Array<OneD, NekDouble> tmp(nq);
Array<OneD, NekDouble> mod(ncoeffs);
{
{
// add to u equation
Vmath::Vmul(nq,m_coriolis,1,physarray[2],1,tmp,1);
m_fields[0]->IProductWRTBase(tmp,mod);
m_fields[0]->MultiplyByElmtInvMass(mod,mod);
m_fields[0]->BwdTrans(mod,tmp);
Vmath::Vadd(nq,tmp,1,outarray[1],1,outarray[1],1);
// add to v equation
Vmath::Vmul(nq,m_coriolis,1,physarray[1],1,tmp,1);
Vmath::Neg(nq,tmp,1);
m_fields[0]->IProductWRTBase(tmp,mod);
m_fields[0]->MultiplyByElmtInvMass(mod,mod);
m_fields[0]->BwdTrans(mod,tmp);
Vmath::Vadd(nq,tmp,1,outarray[2],1,outarray[2],1);
}
break;
{
// add to u equation
Vmath::Vmul(nq,m_coriolis,1,physarray[2],1,tmp,1);
Vmath::Vadd(nq,tmp,1,outarray[1],1,outarray[1],1);
// add to v equation
Vmath::Vmul(nq,m_coriolis,1,physarray[1],1,tmp,1);
Vmath::Neg(nq,tmp,1);
Vmath::Vadd(nq,tmp,1,outarray[2],1,outarray[2],1);
}
break;
default:
ASSERTL0(false,"Unknown projection scheme for the NonlinearSWE");
break;
}
}
void Nektar::LinearSWE::ConservativeToPrimitive ( const Array< OneD, const Array< OneD, NekDouble > > &  physin,
Array< OneD, Array< OneD, NekDouble > > &  physout 
)
private

Definition at line 584 of file LinearSWE.cpp.

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

{
int nq = GetTotPoints();
if(physin.get() == physout.get())
{
// copy indata and work with tmp array
Array<OneD, Array<OneD, NekDouble> >tmp(3);
for (int i = 0; i < 3; ++i)
{
// deep copy
tmp[i] = Array<OneD, NekDouble>(nq);
Vmath::Vcopy(nq,physin[i],1,tmp[i],1);
}
// \eta = h - d
Vmath::Vsub(nq,tmp[0],1,m_depth,1,physout[0],1);
// u = hu/h
Vmath::Vdiv(nq,tmp[1],1,tmp[0],1,physout[1],1);
// v = hv/ v
Vmath::Vdiv(nq,tmp[2],1,tmp[0],1,physout[2],1);
}
else
{
// \eta = h - d
Vmath::Vsub(nq,physin[0],1,m_depth,1,physout[0],1);
// u = hu/h
Vmath::Vdiv(nq,physin[1],1,physin[0],1,physout[1],1);
// v = hv/ v
Vmath::Vdiv(nq,physin[2],1,physin[0],1,physout[2],1);
}
}
static SolverUtils::EquationSystemSharedPtr Nektar::LinearSWE::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 56 of file LinearSWE.h.

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

Definition at line 326 of file LinearSWE.cpp.

References ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, SetBoundaryConditions(), and Vmath::Vcopy().

Referenced by v_InitObject().

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

Definition at line 223 of file LinearSWE.cpp.

References AddCoriolis(), ASSERTL0, Nektar::MultiRegions::DirCartesianMap, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, GetFluxVector(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::ShallowWaterSystem::m_advection, Nektar::ShallowWaterSystem::m_coriolis, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Neg(), and Vmath::Vadd().

Referenced by v_InitObject().

{
int i, j;
int ndim = m_spacedim;
int nvariables = inarray.num_elements();
int nq = GetTotPoints();
{
{
//-------------------------------------------------------
// Compute the DG advection including the numerical flux
// by using SolverUtils/Advection
// Input and output in physical space
Array<OneD, Array<OneD, NekDouble> > advVel;
m_advection->Advect(nvariables, m_fields, advVel, inarray,
outarray, time);
//-------------------------------------------------------
//-------------------------------------------------------
// negate the outarray since moving terms to the rhs
for(i = 0; i < nvariables; ++i)
{
Vmath::Neg(nq,outarray[i],1);
}
//-------------------------------------------------------
//-------------------------------------------------
// Add "source terms"
// Input and output in physical space
// Coriolis forcing
if (m_coriolis.num_elements() != 0)
{
AddCoriolis(inarray,outarray);
}
//-------------------------------------------------
}
break;
{
//-------------------------------------------------------
// Compute the fluxvector in physical space
Array<OneD, Array<OneD, Array<OneD, NekDouble> > >
fluxvector(nvariables);
for (i = 0; i < nvariables; ++i)
{
fluxvector[i] = Array<OneD, Array<OneD, NekDouble> >(ndim);
for(j = 0; j < ndim; ++j)
{
fluxvector[i][j] = Array<OneD, NekDouble>(nq);
}
}
LinearSWE::GetFluxVector(inarray, fluxvector);
//-------------------------------------------------------
//-------------------------------------------------------
// Take the derivative of the flux terms
// and negate the outarray since moving terms to the rhs
Array<OneD,NekDouble> tmp(nq);
Array<OneD, NekDouble>tmp1(nq);
for(i = 0; i < nvariables; ++i)
{
m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[0],fluxvector[i][0],tmp);
m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[1],fluxvector[i][1],tmp1);
Vmath::Vadd(nq,tmp,1,tmp1,1,outarray[i],1);
Vmath::Neg(nq,outarray[i],1);
}
//-------------------------------------------------
// Add "source terms"
// Input and output in physical space
// Coriolis forcing
if (m_coriolis.num_elements() != 0)
{
AddCoriolis(inarray,outarray);
}
//-------------------------------------------------
}
break;
default:
ASSERTL0(false,"Unknown projection scheme for the NonlinearSWE");
break;
}
}
const Array<OneD, NekDouble>& Nektar::LinearSWE::GetDepthBwd ( )
inlineprotected

Definition at line 100 of file LinearSWE.h.

References m_dBwd.

Referenced by v_InitObject().

{
return m_dBwd;
}
const Array<OneD, NekDouble>& Nektar::LinearSWE::GetDepthFwd ( )
inlineprotected

Definition at line 96 of file LinearSWE.h.

References m_dFwd.

Referenced by v_InitObject().

{
return m_dFwd;
}
void Nektar::LinearSWE::GetFluxVector ( const Array< OneD, const Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Definition at line 550 of file LinearSWE.cpp.

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

Referenced by DoOdeRhs(), and v_InitObject().

{
int i, j;
int nq = m_fields[0]->GetTotPoints();
// Flux vector for the mass equation
for (i = 0; i < m_spacedim; ++i)
{
Vmath::Vmul(nq, m_depth, 1, physfield[i+1], 1, flux[0][i], 1);
}
// Put (g eta) in tmp
Array<OneD, NekDouble> tmp(nq);
Vmath::Smul(nq, g, physfield[0], 1, tmp, 1);
// Flux vector for the momentum equations
for (i = 0; i < m_spacedim; ++i)
{
for (j = 0; j < m_spacedim; ++j)
{
// must zero fluxes as not initialised to zero in AdvectionWeakDG ...
Vmath::Zero(nq, flux[i+1][j], 1);
}
// Add (g eta) to appropriate field
Vmath::Vadd(nq, flux[i+1][i], 1, tmp, 1, flux[i+1][i], 1);
}
}
void Nektar::LinearSWE::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
physfieldVelocity field.
velocityVelocity field.

Definition at line 701 of file LinearSWE.cpp.

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

{
const int npts = physfield[0].num_elements();
for (int i = 0; i < m_spacedim; ++i)
{
Vmath::Vcopy(npts, physfield[1+i], 1, velocity[i], 1);
}
}
void Nektar::LinearSWE::NumericalFlux1D ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numfluxX 
)
private
void Nektar::LinearSWE::NumericalFlux2D ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numfluxX,
Array< OneD, Array< OneD, NekDouble > > &  numfluxY 
)
private
void Nektar::LinearSWE::PrimitiveToConservative ( const Array< OneD, const Array< OneD, NekDouble > > &  physin,
Array< OneD, Array< OneD, NekDouble > > &  physout 
)
private

Definition at line 637 of file LinearSWE.cpp.

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

{
int nq = GetTotPoints();
if(physin.get() == physout.get())
{
// copy indata and work with tmp array
Array<OneD, Array<OneD, NekDouble> >tmp(3);
for (int i = 0; i < 3; ++i)
{
// deep copy
tmp[i] = Array<OneD, NekDouble>(nq);
Vmath::Vcopy(nq,physin[i],1,tmp[i],1);
}
// h = \eta + d
Vmath::Vadd(nq,tmp[0],1,m_depth,1,physout[0],1);
// hu = h * u
Vmath::Vmul(nq,physout[0],1,tmp[1],1,physout[1],1);
// hv = h * v
Vmath::Vmul(nq,physout[0],1,tmp[2],1,physout[2],1);
}
else
{
// h = \eta + d
Vmath::Vadd(nq,physin[0],1,m_depth,1,physout[0],1);
// hu = h * u
Vmath::Vmul(nq,physout[0],1,physin[1],1,physout[1],1);
// hv = h * v
Vmath::Vmul(nq,physout[0],1,physin[2],1,physout[2],1);
}
}
void Nektar::LinearSWE::SetBoundaryConditions ( Array< OneD, Array< OneD, NekDouble > > &  physarray,
NekDouble  time 
)
private

Definition at line 371 of file LinearSWE.cpp.

References Nektar::SpatialDomains::eTimeDependent, Nektar::SpatialDomains::eWall, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, and WallBoundary2D().

Referenced by DoOdeProjection().

{
std::string varName;
int nvariables = m_fields.num_elements();
int cnt = 0;
// loop over Boundary Regions
for(int n = 0; n < m_fields[0]->GetBndConditions().num_elements(); ++n)
{
// Wall Boundary Condition
if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
{
WallBoundary2D(n, cnt, inarray);
}
// Time Dependent Boundary Condition (specified in meshfile)
if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
{
for (int i = 0; i < nvariables; ++i)
{
varName = m_session->GetVariable(i);
m_fields[i]->EvaluateBoundaryConditions(time, varName);
}
}
cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
}
}
void Nektar::LinearSWE::v_ConservativeToPrimitive ( )
protectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 623 of file LinearSWE.cpp.

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

{
int nq = GetTotPoints();
// u = hu/h
Vmath::Vdiv(nq,m_fields[1]->GetPhys(),1,m_fields[0]->GetPhys(),1,m_fields[1]->UpdatePhys(),1);
// v = hv/ v
Vmath::Vdiv(nq,m_fields[2]->GetPhys(),1,m_fields[0]->GetPhys(),1,m_fields[2]->UpdatePhys(),1);
// \eta = h - d
Vmath::Vsub(nq,m_fields[0]->GetPhys(),1,m_depth,1,m_fields[0]->UpdatePhys(),1);
}
void Nektar::LinearSWE::v_GenerateSummary ( SolverUtils::SummaryList s)
protectedvirtual

Print a summary of time stepping parameters.

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 714 of file LinearSWE.cpp.

References Nektar::SolverUtils::AddSummaryItem(), and Nektar::SolverUtils::EquationSystem::m_session.

{
if (m_session->DefinesSolverInfo("UpwindType"))
{
std::string UpwindType;
UpwindType = m_session->GetSolverInfo("UpwindType");
if (UpwindType == "LinearAverage")
{
SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear Average");
}
if (UpwindType == "LinearHLL")
{
SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear HLL");
}
}
SolverUtils::AddSummaryItem(s, "Variables", "eta should be in field[0]");
SolverUtils::AddSummaryItem(s, "", "u should be in field[1]");
SolverUtils::AddSummaryItem(s, "", "v should be in field[2]");
}
void Nektar::LinearSWE::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 56 of file LinearSWE.cpp.

References ASSERTL0, Nektar::ShallowWaterSystem::CopyBoundaryTrace(), Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoOdeProjection(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::SolverUtils::GetAdvectionFactory(), GetDepthBwd(), GetDepthFwd(), GetFluxVector(), Nektar::ShallowWaterSystem::GetGravity(), Nektar::ShallowWaterSystem::GetNormals(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::ShallowWaterSystem::GetVecLocs(), Nektar::ShallowWaterSystem::m_advection, Nektar::ShallowWaterSystem::m_constantDepth, m_dBwd, Nektar::ShallowWaterSystem::m_depth, m_dFwd, Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::ShallowWaterSystem::m_riemannSolver, and Nektar::SolverUtils::EquationSystem::m_session.

{
{
}
else
{
ASSERTL0(false, "Implicit SWE not set up.");
}
// Type of advection class to be used
{
// Continuous field
{
// Do nothing
break;
}
// Discontinuous field
{
string advName;
string diffName;
string riemName;
//---------------------------------------------------------------
// Setting up advection and diffusion operators
// NB: diffusion not set up for SWE at the moment
// but kept here for future use ...
m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
// m_session->LoadSolverInfo("DiffusionType", diffName, "LDGEddy");
.CreateInstance(advName, advName);
// m_diffusion = SolverUtils::GetDiffusionFactory()
// .CreateInstance(diffName, diffName);
m_advection->SetFluxVector(&LinearSWE::GetFluxVector, this);
// m_diffusion->SetFluxVectorNS(&ShallowWaterSystem::
// GetEddyViscosityFluxVector, this);
// Setting up Riemann solver for advection operator
m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
if ((riemName == "LinearHLL") && (m_constantDepth != true))
{
ASSERTL0(false,"LinearHLL only valid for constant depth");
}
.CreateInstance(riemName);
// Setting up upwind solver for diffusion operator
// m_riemannSolverLDG = SolverUtils::GetRiemannSolverFactory()
// .CreateInstance("UpwindLDG");
// Setting up parameters for advection operator Riemann solver
m_riemannSolver->SetParam (
"gravity",
m_riemannSolver->SetAuxVec(
"vecLocs",
m_riemannSolver->SetVector(
"N",
// The numerical flux for linear SWE requires depth information
int nTracePointsTot = m_fields[0]->GetTrace()->GetTotPoints();
m_dFwd = Array<OneD, NekDouble>(nTracePointsTot);
m_dBwd = Array<OneD, NekDouble>(nTracePointsTot);
m_fields[0]->GetFwdBwdTracePhys(m_depth, m_dFwd, m_dBwd);
m_riemannSolver->SetScalar(
"depthFwd",
m_riemannSolver->SetScalar(
"depthBwd",
// Setting up parameters for diffusion operator Riemann solver
// m_riemannSolverLDG->AddParam (
// "gravity",
// &NonlinearSWE::GetGravity, this);
// m_riemannSolverLDG->SetAuxVec(
// "vecLocs",
// &NonlinearSWE::GetVecLocs, this);
// m_riemannSolverLDG->AddVector(
// "N",
// &NonlinearSWE::GetNormals, this);
// Concluding initialisation of advection / diffusion operators
m_advection->SetRiemannSolver (m_riemannSolver);
//m_diffusion->SetRiemannSolver (m_riemannSolverLDG);
m_advection->InitObject (m_session, m_fields);
//m_diffusion->InitObject (m_session, m_fields);
break;
}
default:
{
ASSERTL0(false, "Unsupported projection type.");
break;
}
}
}
void Nektar::LinearSWE::v_PrimitiveToConservative ( )
protectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 679 of file LinearSWE.cpp.

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

{
int nq = GetTotPoints();
// h = \eta + d
Vmath::Vadd(nq,m_fields[0]->GetPhys(),1,m_depth,1,m_fields[0]->UpdatePhys(),1);
// hu = h * u
Vmath::Vmul(nq,m_fields[0]->GetPhys(),1,m_fields[1]->GetPhys(),1,m_fields[1]->UpdatePhys(),1);
// hv = h * v
Vmath::Vmul(nq,m_fields[0]->GetPhys(),1,m_fields[2]->GetPhys(),1,m_fields[2]->UpdatePhys(),1);
}
void Nektar::LinearSWE::WallBoundary ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  physarray 
)
private

Wall boundary condition.

Definition at line 407 of file LinearSWE.cpp.

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

{
int i;
int nTracePts = GetTraceTotPoints();
int nvariables = physarray.num_elements();
// get physical values of the forward trace
Array<OneD, Array<OneD, NekDouble> > Fwd(nvariables);
for (i = 0; i < nvariables; ++i)
{
Fwd[i] = Array<OneD, NekDouble>(nTracePts);
m_fields[i]->ExtractTracePhys(physarray[i], Fwd[i]);
}
// Adjust the physical values of the trace to take
// user defined boundaries into account
int e, id1, id2, npts;
for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]
->GetExpSize(); ++e)
{
npts = m_fields[0]->GetBndCondExpansions()[bcRegion]->
GetExp(e)->GetTotPoints();
id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->
id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
m_fields[0]->GetTraceMap()->
GetBndCondCoeffsToGlobalCoeffsMap(cnt+e));
// For 2D/3D, define: v* = v - 2(v.n)n
Array<OneD, NekDouble> tmp(npts, 0.0);
// Calculate (v.n)
for (i = 0; i < m_spacedim; ++i)
{
&Fwd[1+i][id2], 1,
&m_traceNormals[i][id2], 1,
&tmp[0], 1,
&tmp[0], 1);
}
// Calculate 2.0(v.n)
Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
// Calculate v* = v - 2.0(v.n)n
for (i = 0; i < m_spacedim; ++i)
{
&tmp[0], 1,
&m_traceNormals[i][id2], 1,
&Fwd[1+i][id2], 1,
&Fwd[1+i][id2], 1);
}
// copy boundary adjusted values into the boundary expansion
for (i = 0; i < nvariables; ++i)
{
Vmath::Vcopy(npts, &Fwd[i][id2], 1,
&(m_fields[i]->GetBndCondExpansions()[bcRegion]->
UpdatePhys())[id1], 1);
}
}
}
void Nektar::LinearSWE::WallBoundary2D ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  physarray 
)
private

Definition at line 476 of file LinearSWE.cpp.

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

Referenced by SetBoundaryConditions().

{
int i;
int nTraceNumPoints = GetTraceTotPoints();
int nvariables = physarray.num_elements();
// get physical values of the forward trace
Array<OneD, Array<OneD, NekDouble> > Fwd(nvariables);
for (i = 0; i < nvariables; ++i)
{
Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
m_fields[i]->ExtractTracePhys(physarray[i],Fwd[i]);
}
// Adjust the physical values of the trace to take
// user defined boundaries into account
int e, id1, id2, npts;
for(e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize(); ++e)
{
npts = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExp(e)->GetNumPoints(0);
id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e) ;
id2 = m_fields[0]->GetTrace()->GetPhys_Offset(m_fields[0]->GetTraceMap()->GetBndCondCoeffsToGlobalCoeffsMap(cnt+e));
switch(m_expdim)
{
case 1:
{
// negate the forward flux
Vmath::Neg(npts,&Fwd[1][id2],1);
}
break;
case 2:
{
Array<OneD, NekDouble> tmp_n(npts);
Array<OneD, NekDouble> tmp_t(npts);
Vmath::Vmul(npts,&Fwd[1][id2],1,&m_traceNormals[0][id2],1,&tmp_n[0],1);
Vmath::Vvtvp(npts,&Fwd[2][id2],1,&m_traceNormals[1][id2],1,&tmp_n[0],1,&tmp_n[0],1);
Vmath::Vmul(npts,&Fwd[1][id2],1,&m_traceNormals[1][id2],1,&tmp_t[0],1);
Vmath::Vvtvm(npts,&Fwd[2][id2],1,&m_traceNormals[0][id2],1,&tmp_t[0],1,&tmp_t[0],1);
// negate the normal flux
Vmath::Neg(npts,tmp_n,1);
// rotate back to Cartesian
Vmath::Vmul(npts,&tmp_t[0],1,&m_traceNormals[1][id2],1,&Fwd[1][id2],1);
Vmath::Vvtvm(npts,&tmp_n[0],1,&m_traceNormals[0][id2],1,&Fwd[1][id2],1,&Fwd[1][id2],1);
Vmath::Vmul(npts,&tmp_t[0],1,&m_traceNormals[0][id2],1,&Fwd[2][id2],1);
Vmath::Vvtvp(npts,&tmp_n[0],1,&m_traceNormals[1][id2],1,&Fwd[2][id2],1,&Fwd[2][id2],1);
}
break;
case 3:
ASSERTL0(false,"3D not implemented for Shallow Water Equations");
break;
default:
ASSERTL0(false,"Illegal expansion dimension");
}
// copy boundary adjusted values into the boundary expansion
for (i = 0; i < nvariables; ++i)
{
Vmath::Vcopy(npts,&Fwd[i][id2], 1,&(m_fields[i]->GetBndCondExpansions()[bcRegion]->UpdatePhys())[id1],1);
}
}
}

Friends And Related Function Documentation

friend class MemoryManager< LinearSWE >
friend

Definition at line 53 of file LinearSWE.h.

Member Data Documentation

string Nektar::LinearSWE::className
static
Initial value:
"LinearSWE", LinearSWE::create,
"Linear shallow water equation in primitive variables.")

Name of class.

Definition at line 64 of file LinearSWE.h.

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

Definition at line 76 of file LinearSWE.h.

Referenced by GetDepthBwd(), and v_InitObject().

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

Still water depth traces.

Definition at line 75 of file LinearSWE.h.

Referenced by GetDepthFwd(), and v_InitObject().