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

#include <MMFAdvection.h>

Inheritance diagram for Nektar::MMFAdvection:
[legend]

Public Member Functions

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

Static Public Member Functions

static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Public Attributes

TestType m_TestType
 
- Public Attributes inherited from Nektar::SolverUtils::MMFSystem
NekDouble m_pi
 
int m_shapedim
 
SurfaceType m_surfaceType
 
UpwindType m_upwindType
 
TestMaxwellType m_TestMaxwellType
 
PolType m_PolType
 
IncType m_IncType
 
Array< OneD, NekDoublem_MMFfactors
 

Static Public Attributes

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

 MMFAdvection (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Session reader. More...
 
void WeakDGDirectionalAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField)
 
void GetFluxVector (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 Evaluate the flux at each solution point. More...
 
void GetFluxVectorDeAlias (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 Evaluate the flux at each solution point using dealiasing. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the RHS. More...
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the projection. More...
 
void EvaluateAdvectionVelocity (Array< OneD, Array< OneD, NekDouble > > &velocity)
 
NekDouble ComputeCirculatingArclength (const NekDouble zlevel, const NekDouble Rhs)
 
Array< OneD, NekDouble > & GetNormalVelocity ()
 Get the normal velocity. More...
 
void ComputeNablaCdotVelocity (Array< OneD, NekDouble > &vellc)
 
void ComputeveldotMF (Array< OneD, Array< OneD, NekDouble > > &veldotMF)
 
void AdvectionBellPlane (Array< OneD, NekDouble > &outfield)
 
void AdvectionBellSphere (Array< OneD, NekDouble > &outfield)
 
void Test2Dproblem (const NekDouble time, Array< OneD, NekDouble > &outfield)
 
void Test3Dproblem (const NekDouble time, Array< OneD, NekDouble > &outfield)
 
void v_InitObject (bool DeclareFields=true) override
 Initialise the object. More...
 
void v_DoSolve () override
 Solves an unsteady problem. More...
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print Summary. More...
 
void v_SetInitialConditions (const NekDouble initialtime, bool dumpInitialConditions, const int domain) override
 
void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time) override
 
- Protected Member Functions inherited from Nektar::SolverUtils::MMFSystem
void SetUpMovingFrames (const Array< OneD, const Array< OneD, NekDouble > > &Anisotropy, const int TangentXelem)
 
void CheckMovingFrames (const Array< OneD, const Array< OneD, NekDouble > > &movingframes)
 
SOLVER_UTILS_EXPORT void ComputencdotMF ()
 
SOLVER_UTILS_EXPORT void ComputeDivCurlMF ()
 
SOLVER_UTILS_EXPORT void ComputeMFtrace ()
 
SOLVER_UTILS_EXPORT void VectorDotProd (const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, NekDouble > &v3)
 
SOLVER_UTILS_EXPORT void VectorCrossProd (const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &v3)
 
SOLVER_UTILS_EXPORT void VectorCrossProd (const Array< OneD, NekDouble > &v1, const Array< OneD, NekDouble > &v2, Array< OneD, NekDouble > &v3)
 
SOLVER_UTILS_EXPORT void ComputeCurl (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleCartesianToMovingframes (const Array< OneD, const Array< OneD, NekDouble > > &uvec, unsigned int field)
 
SOLVER_UTILS_EXPORT void DeriveCrossProductMF (Array< OneD, Array< OneD, NekDouble > > &CrossProductMF)
 
SOLVER_UTILS_EXPORT void ComputeNtimesMF ()
 
SOLVER_UTILS_EXPORT void ComputeNtimesFz (const int dir, const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &imFwd, const Array< OneD, const NekDouble > &imBwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
 
SOLVER_UTILS_EXPORT void ComputeNtimesF12 (const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &im1Fwd, const Array< OneD, const NekDouble > &im1Bwd, const Array< OneD, const NekDouble > &im2Fwd, const Array< OneD, const NekDouble > &im2Bwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
 
SOLVER_UTILS_EXPORT void ComputeNtimestimesdFz (const int dir, const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &imFwd, const Array< OneD, const NekDouble > &imBwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
 
SOLVER_UTILS_EXPORT void ComputeNtimestimesdF12 (const Array< OneD, Array< OneD, NekDouble > > &Fwd, const Array< OneD, Array< OneD, NekDouble > > &Bwd, const Array< OneD, const NekDouble > &im1Fwd, const Array< OneD, const NekDouble > &im1Bwd, const Array< OneD, const NekDouble > &im2Fwd, const Array< OneD, const NekDouble > &im2Bwd, Array< OneD, NekDouble > &outarrayFwd, Array< OneD, NekDouble > &outarrayBwd)
 
SOLVER_UTILS_EXPORT void CartesianToSpherical (const NekDouble x0j, const NekDouble x1j, const NekDouble x2j, NekDouble &sin_varphi, NekDouble &cos_varphi, NekDouble &sin_theta, NekDouble &cos_theta)
 
SOLVER_UTILS_EXPORT void ComputeZimYim (Array< OneD, Array< OneD, NekDouble > > &epsvec, Array< OneD, Array< OneD, NekDouble > > &muvec)
 
SOLVER_UTILS_EXPORT void AdddedtMaxwell (const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void GetMaxwellFluxVector (const int var, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetMaxwellFlux1D (const int var, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetMaxwellFlux2D (const int var, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void LaxFriedrichMaxwellFlux1D (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd)
 
SOLVER_UTILS_EXPORT void UpwindMaxwellFlux1D (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd)
 
SOLVER_UTILS_EXPORT void AverageMaxwellFlux1D (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd)
 
SOLVER_UTILS_EXPORT void NumericalMaxwellFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd, const NekDouble time=0.0)
 
SOLVER_UTILS_EXPORT void NumericalMaxwellFluxTM (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd, const NekDouble time)
 
SOLVER_UTILS_EXPORT void NumericalMaxwellFluxTE (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxFwd, Array< OneD, Array< OneD, NekDouble > > &numfluxBwd, const NekDouble time)
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetIncidentField (const int var, const NekDouble time)
 
SOLVER_UTILS_EXPORT void Computedemdxicdote ()
 
SOLVER_UTILS_EXPORT NekDouble AvgInt (const Array< OneD, const NekDouble > &inarray)
 
SOLVER_UTILS_EXPORT NekDouble AvgAbsInt (const Array< OneD, const NekDouble > &inarray)
 
SOLVER_UTILS_EXPORT NekDouble AbsIntegral (const Array< OneD, const NekDouble > &inarray)
 
SOLVER_UTILS_EXPORT NekDouble RootMeanSquare (const Array< OneD, const NekDouble > &inarray)
 
SOLVER_UTILS_EXPORT NekDouble VectorAvgMagnitude (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 
SOLVER_UTILS_EXPORT void GramSchumitz (const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &outarray, bool KeepTheMagnitude=true)
 
SOLVER_UTILS_EXPORT void BubbleSort (Array< OneD, NekDouble > &refarray, Array< OneD, NekDouble > &sortarray)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics. More...
 
SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble > > vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform dummy projection. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
- Protected Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step) override
 
virtual SOLVER_UTILS_EXPORT Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor=1.0)
 

Protected Attributes

SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
NekDouble m_advx
 
NekDouble m_advy
 
NekDouble m_advz
 
NekDouble m_waveFreq
 
NekDouble m_RotAngle
 
NekDouble m_Mass0
 
int m_VelProjection
 
Array< OneD, Array< OneD, NekDouble > > m_velocity
 Advection velocity. More...
 
Array< OneD, NekDoublem_traceVn
 
Array< OneD, Array< OneD, NekDouble > > m_veldotMF
 
Array< OneD, NekDoublem_vellc
 
- Protected Attributes inherited from Nektar::SolverUtils::MMFSystem
NekDouble m_alpha
 
NekDouble m_Incfreq
 
int m_SmoothFactor
 
NekDouble m_SFinit
 
Array< OneD, Array< OneD, NekDouble > > m_movingframes
 
Array< OneD, Array< OneD, NekDouble > > m_surfaceNormal
 
Array< OneD, Array< OneD, NekDouble > > m_ncdotMFFwd
 
Array< OneD, Array< OneD, NekDouble > > m_ncdotMFBwd
 
Array< OneD, Array< OneD, NekDouble > > m_nperpcdotMFFwd
 
Array< OneD, Array< OneD, NekDouble > > m_nperpcdotMFBwd
 
Array< OneD, Array< OneD, NekDouble > > m_DivMF
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_CurlMF
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_MFtraceFwd
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_MFtraceBwd
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimesMFFwd
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimesMFBwd
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimes_ntimesMFFwd
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_ntimes_ntimesMFBwd
 
Array< OneD, Array< OneD, NekDouble > > m_ZimFwd
 
Array< OneD, Array< OneD, NekDouble > > m_ZimBwd
 
Array< OneD, Array< OneD, NekDouble > > m_YimFwd
 
Array< OneD, Array< OneD, NekDouble > > m_YimBwd
 
Array< OneD, Array< OneD, NekDouble > > m_epsvec
 
Array< OneD, Array< OneD, NekDouble > > m_muvec
 
Array< OneD, Array< OneD, NekDouble > > m_negepsvecminus1
 
Array< OneD, Array< OneD, NekDouble > > m_negmuvecminus1
 
Array< OneD, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > > m_dedxi_cdot_e
 
SpatialDomains::GeomMMF m_MMFdir
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
Array< OneD, NekDoublem_movingFrameData
 Moving reference frame status in the inertial frame X, Y, Z, Theta_x, Theta_y, Theta_z, U, V, W, Omega_x, Omega_y, Omega_z, A_x, A_y, A_z, DOmega_x, DOmega_y, DOmega_z, pivot_x, pivot_y, pivot_z. More...
 
std::vector< std::string > m_strFrameData
 variable name in m_movingFrameData More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 
- Protected Attributes inherited from Nektar::SolverUtils::ALEHelper
Array< OneD, MultiRegions::ExpListSharedPtrm_fieldsALE
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 
std::vector< ALEBaseShPtrm_ALEs
 
bool m_ALESolver = false
 
bool m_ImplicitALESolver = false
 
NekDouble m_prevStageTime = 0.0
 
int m_spaceDim
 
- Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem
SolverUtils::AdvectionSharedPtr m_advObject
 Advection term. More...
 

Friends

class MemoryManager< MMFAdvection >
 

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 []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

Definition at line 61 of file MMFAdvection.h.

Constructor & Destructor Documentation

◆ ~MMFAdvection()

Nektar::MMFAdvection::~MMFAdvection ( )
overridedefault

Destructor.

◆ MMFAdvection()

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

Session reader.

Definition at line 50 of file MMFAdvection.cpp.

52 : UnsteadySystem(pSession, pGraph), MMFSystem(pSession, pGraph),
53 AdvectionSystem(pSession, pGraph)
54{
55}
SOLVER_UTILS_EXPORT AdvectionSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
SOLVER_UTILS_EXPORT MMFSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: MMFSystem.cpp:41
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

Member Function Documentation

◆ AdvectionBellPlane()

void Nektar::MMFAdvection::AdvectionBellPlane ( Array< OneD, NekDouble > &  outfield)
protected

Definition at line 901 of file MMFAdvection.cpp.

902{
903 int nq = m_fields[0]->GetNpoints();
904
905 NekDouble dist, m_radius_limit;
906
907 Array<OneD, NekDouble> x(nq);
908 Array<OneD, NekDouble> y(nq);
909 Array<OneD, NekDouble> z(nq);
910
911 m_fields[0]->GetCoords(x, y, z);
912
913 // Sets of parameters
914 m_radius_limit = 0.5;
915
916 NekDouble x0j, x1j;
917 outfield = Array<OneD, NekDouble>(nq, 0.0);
918 for (int j = 0; j < nq; ++j)
919 {
920 x0j = x[j];
921 x1j = y[j];
922
923 dist = sqrt(x0j * x0j + x1j * x1j);
924
925 if (dist < m_radius_limit)
926 {
927 outfield[j] = 0.5 * (1.0 + cos(m_pi * dist / m_radius_limit));
928 }
929 else
930 {
931 outfield[j] = 0.0;
932 }
933 }
934}
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
std::vector< double > z(NPUPPER)
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::MMFSystem::m_pi, tinysimd::sqrt(), and Nektar::UnitTests::z().

Referenced by v_EvaluateExactSolution(), and v_SetInitialConditions().

◆ AdvectionBellSphere()

void Nektar::MMFAdvection::AdvectionBellSphere ( Array< OneD, NekDouble > &  outfield)
protected

Definition at line 936 of file MMFAdvection.cpp.

937{
938 int nq = m_fields[0]->GetNpoints();
939
940 NekDouble dist, radius, cosdiff, sin_theta, cos_theta, sin_varphi,
941 cos_varphi;
942 NekDouble m_theta_c, m_varphi_c, m_radius_limit, m_c0;
943
944 Array<OneD, NekDouble> x(nq);
945 Array<OneD, NekDouble> y(nq);
946 Array<OneD, NekDouble> z(nq);
947
948 m_fields[0]->GetCoords(x, y, z);
949
950 // Sets of parameters
951 m_theta_c = 0.0;
952 m_varphi_c = 3.0 * m_pi / 2.0;
953 m_radius_limit = 7.0 * m_pi / 64.0;
954 m_c0 = 0.0;
955
956 NekDouble x0j, x1j, x2j;
957 outfield = Array<OneD, NekDouble>(nq, 0.0);
958 for (int j = 0; j < nq; ++j)
959 {
960 x0j = x[j];
961 x1j = y[j];
962 x2j = z[j];
963
964 radius = sqrt(x0j * x0j + x1j * x1j + x2j * x2j);
965
966 sin_varphi = x1j / sqrt(x0j * x0j + x1j * x1j);
967 cos_varphi = x0j / sqrt(x0j * x0j + x1j * x1j);
968
969 sin_theta = x2j / radius;
970 cos_theta = sqrt(x0j * x0j + x1j * x1j) / radius;
971
972 cosdiff = cos_varphi * cos(m_varphi_c) + sin_varphi * sin(m_varphi_c);
973 dist = radius * acos(sin(m_theta_c) * sin_theta +
974 cos(m_theta_c) * cos_theta * cosdiff);
975
976 if (dist < m_radius_limit)
977 {
978 outfield[j] =
979 0.5 * (1.0 + cos(m_pi * dist / m_radius_limit)) + m_c0;
980 }
981 else
982 {
983 outfield[j] = m_c0;
984 }
985 }
986}

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::MMFSystem::m_pi, tinysimd::sqrt(), and Nektar::UnitTests::z().

Referenced by v_EvaluateExactSolution(), and v_SetInitialConditions().

◆ ComputeCirculatingArclength()

NekDouble Nektar::MMFAdvection::ComputeCirculatingArclength ( const NekDouble  zlevel,
const NekDouble  Rhs 
)
protected

Definition at line 678 of file MMFAdvection.cpp.

680{
681
682 NekDouble Tol = 0.0001, Maxiter = 1000, N = 100;
683 NekDouble newy, F = 0.0, dF = 1.0, y0, tmp;
684
685 Array<OneD, NekDouble> xp(N + 1);
686 Array<OneD, NekDouble> yp(N + 1);
687
688 NekDouble intval = 0.0;
689 switch (m_surfaceType)
690 {
693 {
694 intval = sqrt(Rhs - zlevel * zlevel);
695 }
696 break;
697
699 {
700 intval = sqrt(0.5 * (Rhs - zlevel * zlevel * zlevel * zlevel -
701 zlevel * zlevel));
702 }
703 break;
704
706 {
707 tmp = 0.5 *
708 (Rhs - zlevel * zlevel * zlevel * zlevel - zlevel * zlevel);
709 intval = sqrt(0.5 * (1.0 + sqrt(1.0 + 4.0 * tmp)));
710 }
711 break;
712
713 default:
714 break;
715 }
716
717 switch (m_surfaceType)
718 {
719 // Find the half of all the xp and yp on zlevel ....
723 {
724 for (int j = 0; j < N + 1; ++j)
725 {
726 xp[j] = j * 2.0 * intval / N - intval;
727
728 y0 = 1.0;
729 for (int i = 0; i < Maxiter; ++i)
730 {
731 switch (m_surfaceType)
732 {
733 // Find the half of all the xp and yp on zlevel ....
736 {
737 F = xp[j] * xp[j] + y0 * y0 + zlevel * zlevel - Rhs;
738 dF = 2.0 * y0;
739 }
740 break;
741
743 {
744 F = 2.0 * xp[j] * xp[j] + y0 * y0 * y0 * y0 +
745 y0 * y0 + zlevel * zlevel * zlevel * zlevel +
746 zlevel * zlevel - Rhs;
747 dF = 4.0 * y0 * y0 * y0 + 2.0 * y0;
748 }
749 break;
750
751 default:
752 break;
753 }
754
755 newy = y0 - F / dF;
756
757 if (fabs(F / dF) < Tol)
758 {
759 yp[j] = newy;
760 break;
761 }
762
763 else
764 {
765 y0 = newy;
766 }
767
768 ASSERTL0(i < Maxiter,
769 "Advection Velocity convergence fails");
770
771 } // i-loop
772 }
773 }
774 break;
775
777 {
778 for (int j = 0; j < N + 1; ++j)
779 {
780 xp[j] = j * 2.0 * intval / N - intval;
781 tmp = 0.5 * Rhs -
782 0.5 * (zlevel * zlevel * zlevel * zlevel +
783 zlevel * zlevel) -
784 (xp[j] * xp[j] * xp[j] * xp[j] - xp[j] * xp[j]);
785 if (tmp < 0)
786 {
787 tmp = -1.0 * tmp;
788 }
789 yp[j] = sqrt(tmp);
790 } // j-loop
791 }
792 break;
793
794 default:
795 break;
796
797 } // switch-loop
798
799 NekDouble pi = 3.14159265358979323846;
800 NekDouble arclength = 0.0;
801 for (int j = 0; j < N; ++j)
802 {
803 arclength =
804 arclength + sqrt((yp[j + 1] - yp[j]) * (yp[j + 1] - yp[j]) +
805 (xp[j + 1] - xp[j]) * (xp[j + 1] - xp[j])) /
806 pi;
807 }
808
809 return arclength;
810}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208

References ASSERTL0, Nektar::SolverUtils::eIrregular, Nektar::SolverUtils::eNonconvex, Nektar::SolverUtils::eSphere, Nektar::SolverUtils::eTRSphere, Nektar::SolverUtils::MMFSystem::m_surfaceType, and tinysimd::sqrt().

◆ ComputeNablaCdotVelocity()

void Nektar::MMFAdvection::ComputeNablaCdotVelocity ( Array< OneD, NekDouble > &  vellc)
protected

Definition at line 1033 of file MMFAdvection.cpp.

1034{
1035 int nq = m_fields[0]->GetNpoints();
1036
1037 Array<OneD, NekDouble> velcoeff(nq, 0.0);
1038
1039 Array<OneD, NekDouble> Dtmp0(nq);
1040 Array<OneD, NekDouble> Dtmp1(nq);
1041 Array<OneD, NekDouble> Dtmp2(nq);
1042 Array<OneD, NekDouble> Drv(nq);
1043
1044 vellc = Array<OneD, NekDouble>(nq, 0.0);
1045
1046 // m_vellc = \nabla m_vel \cdot tan_i
1047 Array<OneD, NekDouble> tmp(nq);
1048 Array<OneD, NekDouble> vessel(nq);
1049
1050 for (int j = 0; j < m_shapedim; ++j)
1051 {
1052 Vmath::Zero(nq, velcoeff, 1);
1053 for (int k = 0; k < m_spacedim; ++k)
1054 {
1055 // a_j = tan_j cdot m_vel
1056 Vmath::Vvtvp(nq, &m_movingframes[j][k * nq], 1, &m_velocity[k][0],
1057 1, &velcoeff[0], 1, &velcoeff[0], 1);
1058 }
1059
1060 // d a_j / d x^k
1061 m_fields[0]->PhysDeriv(velcoeff, Dtmp0, Dtmp1, Dtmp2);
1062
1063 for (int k = 0; k < m_spacedim; ++k)
1064 {
1065 // tan_j^k ( d a_j / d x^k )
1066 switch (k)
1067 {
1068 case (0):
1069 {
1070 Vmath::Vvtvp(nq, &Dtmp0[0], 1, &m_movingframes[j][k * nq],
1071 1, &vellc[0], 1, &vellc[0], 1);
1072 }
1073 break;
1074
1075 case (1):
1076 {
1077 Vmath::Vvtvp(nq, &Dtmp1[0], 1, &m_movingframes[j][k * nq],
1078 1, &vellc[0], 1, &vellc[0], 1);
1079 }
1080 break;
1081
1082 case (2):
1083 {
1084 Vmath::Vvtvp(nq, &Dtmp2[0], 1, &m_movingframes[j][k * nq],
1085 1, &vellc[0], 1, &vellc[0], 1);
1086 }
1087 break;
1088 }
1089 }
1090 }
1091}
Array< OneD, Array< OneD, NekDouble > > m_velocity
Advection velocity.
Definition: MMFAdvection.h:95
int m_spacedim
Spatial dimension (>= expansion dim).
Array< OneD, Array< OneD, NekDouble > > m_movingframes
Definition: MMFSystem.h:183
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::MMFSystem::m_movingframes, Nektar::SolverUtils::MMFSystem::m_shapedim, Nektar::SolverUtils::EquationSystem::m_spacedim, m_velocity, Vmath::Vvtvp(), and Vmath::Zero().

Referenced by v_InitObject().

◆ ComputeveldotMF()

void Nektar::MMFAdvection::ComputeveldotMF ( Array< OneD, Array< OneD, NekDouble > > &  veldotMF)
protected

Definition at line 1093 of file MMFAdvection.cpp.

1095{
1096 int nq = m_fields[0]->GetNpoints();
1097
1098 veldotMF = Array<OneD, Array<OneD, NekDouble>>(m_shapedim);
1099
1100 Array<OneD, NekDouble> magMF(nq);
1101 for (int j = 0; j < m_shapedim; ++j)
1102 {
1103 veldotMF[j] = Array<OneD, NekDouble>(nq, 0.0);
1104
1105 for (int k = 0; k < m_spacedim; ++k)
1106 {
1107 Vmath::Vvtvp(nq, &m_movingframes[j][k * nq], 1, &m_velocity[k][0],
1108 1, &veldotMF[j][0], 1, &veldotMF[j][0], 1);
1109 }
1110 }
1111}

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::MMFSystem::m_movingframes, Nektar::SolverUtils::MMFSystem::m_shapedim, Nektar::SolverUtils::EquationSystem::m_spacedim, m_velocity, and Vmath::Vvtvp().

Referenced by v_InitObject().

◆ create()

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

Creates an instance of this class.

Definition at line 68 of file MMFAdvection.h.

71 {
74 p->InitObject();
75 return p;
76 }
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::MMFAdvection::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Compute the projection.

Compute the projection for the linear advection equation.

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 447 of file MMFAdvection.cpp.

450{
451 // Counter variable
452 int i;
453
454 // Number of fields (variables of the problem)
455 int nVariables = inarray.size();
456
457 // Set the boundary conditions
459
460 // Switch on the projection type (Discontinuous or Continuous)
461 switch (m_projectionType)
462 {
463 // Discontinuous projection
465 {
466 // Number of quadrature points
467 int nQuadraturePts = GetNpoints();
468
469 // Just copy over array
470 if (inarray != outarray)
471 {
472 for (i = 0; i < nVariables; ++i)
473 {
474 Vmath::Vcopy(nQuadraturePts, inarray[i], 1, outarray[i], 1);
475 }
476 }
477 break;
478 }
479
480 // Continuous projection
483 {
484 Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
485 for (i = 0; i < nVariables; ++i)
486 {
487 m_fields[i]->FwdTrans(inarray[i], coeffs);
488 m_fields[i]->BwdTrans(coeffs, outarray[i]);
489 }
490 break;
491 }
492
493 default:
494 ASSERTL0(false, "Unknown projection scheme");
495 break;
496 }
497}
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

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

Referenced by v_InitObject().

◆ DoOdeRhs()

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

Compute the RHS.

Compute the right-hand side for the linear advection equation.

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 384 of file MMFAdvection.cpp.

388{
389 int i;
390 int nvariables = inarray.size();
391 int npoints = GetNpoints();
392
393 switch (m_projectionType)
394 {
396 {
397 int ncoeffs = inarray[0].size();
398
399 if (m_spacedim == 3)
400 {
401 Array<OneD, Array<OneD, NekDouble>> WeakAdv(nvariables);
402
403 WeakAdv[0] = Array<OneD, NekDouble>(ncoeffs * nvariables);
404 for (i = 1; i < nvariables; ++i)
405 {
406 WeakAdv[i] = WeakAdv[i - 1] + ncoeffs;
407 }
408
409 // Compute \nabla \cdot \vel u according to MMF scheme
410 WeakDGDirectionalAdvection(inarray, WeakAdv);
411
412 for (i = 0; i < nvariables; ++i)
413 {
414 m_fields[i]->MultiplyByElmtInvMass(WeakAdv[i], WeakAdv[i]);
415 m_fields[i]->BwdTrans(WeakAdv[i], outarray[i]);
416 Vmath::Neg(npoints, outarray[i], 1);
417 }
418 }
419 else
420 {
421 m_advObject->Advect(2, m_fields, m_velocity, inarray, outarray,
422 0.0);
423
424 for (i = 0; i < nvariables; ++i)
425 {
426 Vmath::Neg(npoints, outarray[i], 1);
427 }
428 }
429 }
430 break;
431
432 default:
433 {
434 ASSERTL0(false, "Unknown projection scheme");
435 }
436 break;
437 }
438}
void WeakDGDirectionalAdvection(const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField)
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292

References ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_spacedim, m_velocity, Vmath::Neg(), and WeakDGDirectionalAdvection().

Referenced by v_InitObject().

◆ EvaluateAdvectionVelocity()

void Nektar::MMFAdvection::EvaluateAdvectionVelocity ( Array< OneD, Array< OneD, NekDouble > > &  velocity)
protected

Definition at line 596 of file MMFAdvection.cpp.

598{
599 int nq = m_fields[0]->GetNpoints();
600
601 NekDouble vel_phi, vel_theta, sin_varphi, cos_varphi, sin_theta, cos_theta;
602 NekDouble x0j, x1j, x2j;
603
604 Array<OneD, NekDouble> x0(nq);
605 Array<OneD, NekDouble> x1(nq);
606 Array<OneD, NekDouble> x2(nq);
607
608 m_fields[0]->GetCoords(x0, x1, x2);
609
610 // theta = a*sin(z/r), phi = a*tan(y/x);
611 for (int j = 0; j < nq; j++)
612 {
613 x0j = x0[j];
614 x1j = x1[j];
615 x2j = x2[j];
616
617 CartesianToSpherical(x0j, x1j, x2j, sin_varphi, cos_varphi, sin_theta,
618 cos_theta);
619
620 vel_phi = m_waveFreq * (cos_theta * cos(m_RotAngle) +
621 sin_theta * cos_varphi * sin(m_RotAngle));
622 vel_theta = -1.0 * m_waveFreq * sin_theta * sin(m_RotAngle);
623
624 velocity[0][j] =
625 -1.0 * vel_phi * sin_varphi - vel_theta * sin_theta * cos_varphi;
626 velocity[1][j] =
627 vel_phi * cos_varphi - vel_theta * sin_theta * sin_varphi;
628 velocity[2][j] = vel_theta * cos_theta;
629 }
630
631 // Project the veloicty on the tangent plane
632
633 if (m_VelProjection)
634 {
635 // Check MovingFrames \cdot SurfaceNormals = 0
636 Array<OneD, Array<OneD, NekDouble>> newvelocity(m_spacedim);
637
638 Array<OneD, Array<OneD, NekDouble>> MF1(m_spacedim);
639 Array<OneD, Array<OneD, NekDouble>> MF2(m_spacedim);
640 Array<OneD, Array<OneD, NekDouble>> SN(m_spacedim);
641
642 for (int k = 0; k < m_spacedim; ++k)
643 {
644 newvelocity[k] = Array<OneD, NekDouble>(nq);
645 MF1[k] = Array<OneD, NekDouble>(nq);
646 MF2[k] = Array<OneD, NekDouble>(nq);
647 SN[k] = Array<OneD, NekDouble>(nq);
648
649 Vmath::Vcopy(nq, &m_movingframes[0][k * nq], 1, &MF1[k][0], 1);
650 Vmath::Vcopy(nq, &m_movingframes[1][k * nq], 1, &MF2[k][0], 1);
651 }
652
653 VectorCrossProd(MF1, MF2, SN);
654 GramSchumitz(SN, m_velocity, newvelocity, true);
655
656 Array<OneD, NekDouble> tmp(nq, 0.0);
657 Array<OneD, NekDouble> Projection(nq, 0.0);
658
659 for (int k = 0; k < m_spacedim; ++k)
660 {
661 Vmath::Vsub(nq, &m_velocity[0][0], 1, &newvelocity[0][0], 1,
662 &tmp[0], 1);
663 Vmath::Vmul(nq, &tmp[0], 1, &tmp[0], 1, &tmp[0], 1);
664 Vmath::Vadd(nq, tmp, 1, Projection, 1, Projection, 1);
665 }
666
667 std::cout
668 << "Velocity vector is projected onto the tangent plane: Diff = "
669 << RootMeanSquare(Projection) << std::endl;
670
671 for (int k = 0; k < m_spacedim; ++k)
672 {
673 Vmath::Vcopy(nq, &newvelocity[k][0], 1, &m_velocity[k][0], 1);
674 }
675 }
676}
SOLVER_UTILS_EXPORT void VectorCrossProd(const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &v3)
Definition: MMFSystem.cpp:696
SOLVER_UTILS_EXPORT void GramSchumitz(const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &outarray, bool KeepTheMagnitude=true)
Definition: MMFSystem.cpp:2399
SOLVER_UTILS_EXPORT void CartesianToSpherical(const NekDouble x0j, const NekDouble x1j, const NekDouble x2j, NekDouble &sin_varphi, NekDouble &cos_varphi, NekDouble &sin_theta, NekDouble &cos_theta)
Definition: MMFSystem.cpp:792
SOLVER_UTILS_EXPORT NekDouble RootMeanSquare(const Array< OneD, const NekDouble > &inarray)
Definition: MMFSystem.cpp:2338
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void 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.hpp:220

References Nektar::SolverUtils::MMFSystem::CartesianToSpherical(), Nektar::SolverUtils::MMFSystem::GramSchumitz(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::MMFSystem::m_movingframes, m_RotAngle, Nektar::SolverUtils::EquationSystem::m_spacedim, m_velocity, m_VelProjection, m_waveFreq, Nektar::SolverUtils::MMFSystem::RootMeanSquare(), Vmath::Vadd(), Vmath::Vcopy(), Nektar::SolverUtils::MMFSystem::VectorCrossProd(), Vmath::Vmul(), and Vmath::Vsub().

Referenced by v_InitObject().

◆ GetFluxVector()

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

Evaluate the flux at each solution point.

Return the flux vector for the linear advection equation.

Parameters
iComponent of the flux vector to calculate.
physfieldFields.
fluxResulting flux.

Definition at line 506 of file MMFAdvection.cpp.

509{
510 ASSERTL1(flux[0].size() == m_velocity.size(),
511 "Dimension of flux array and velocity array do not match");
512
513 int i, j;
514 int nq = physfield[0].size();
515
516 for (i = 0; i < flux.size(); ++i)
517 {
518 for (j = 0; j < flux[0].size(); ++j)
519 {
520 Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
521 }
522 }
523}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242

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

Referenced by v_InitObject().

◆ GetFluxVectorDeAlias()

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

Evaluate the flux at each solution point using dealiasing.

◆ GetNormalVelocity()

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

Get the normal velocity.

Get the normal velocity for the linear advection equation.

Definition at line 354 of file MMFAdvection.cpp.

355{
356 // Number of trace (interface) points
357 int i;
358 int nTracePts = GetTraceNpoints();
359
360 // Auxiliary variable to compute the normal velocity
361 Array<OneD, NekDouble> tmp(nTracePts);
362
363 // Reset the normal velocity
364 Vmath::Zero(nTracePts, m_traceVn, 1);
365
366 for (i = 0; i < m_velocity.size(); ++i)
367 {
368 m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
369
370 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
371 m_traceVn, 1);
372 }
373
374 return m_traceVn;
375}
Array< OneD, NekDouble > m_traceVn
Definition: MMFAdvection.h:96
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.

References Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, m_traceVn, m_velocity, Vmath::Vvtvp(), and Vmath::Zero().

Referenced by v_InitObject().

◆ Test2Dproblem()

void Nektar::MMFAdvection::Test2Dproblem ( const NekDouble  time,
Array< OneD, NekDouble > &  outfield 
)
protected

Definition at line 988 of file MMFAdvection.cpp.

990{
991 int nq = m_fields[0]->GetNpoints();
992
993 Array<OneD, NekDouble> x0(nq);
994 Array<OneD, NekDouble> x1(nq);
995 Array<OneD, NekDouble> x2(nq);
996
997 m_fields[0]->GetCoords(x0, x1, x2);
998
999 Array<OneD, NekDouble> u(nq);
1000
1001 for (int i = 0; i < nq; ++i)
1002 {
1003 u[i] = cos(m_pi * (x0[i] - m_advx * time)) *
1004 cos(m_pi * (x1[i] - m_advy * time));
1005 }
1006
1007 outfield = u;
1008}

References m_advx, m_advy, Nektar::SolverUtils::EquationSystem::m_fields, and Nektar::SolverUtils::MMFSystem::m_pi.

Referenced by v_EvaluateExactSolution(), and v_SetInitialConditions().

◆ Test3Dproblem()

void Nektar::MMFAdvection::Test3Dproblem ( const NekDouble  time,
Array< OneD, NekDouble > &  outfield 
)
protected

Definition at line 1010 of file MMFAdvection.cpp.

1012{
1013 int nq = m_fields[0]->GetNpoints();
1014
1015 Array<OneD, NekDouble> x0(nq);
1016 Array<OneD, NekDouble> x1(nq);
1017 Array<OneD, NekDouble> x2(nq);
1018
1019 m_fields[0]->GetCoords(x0, x1, x2);
1020
1021 Array<OneD, NekDouble> u(nq);
1022
1023 for (int i = 0; i < nq; ++i)
1024 {
1025 u[i] = cos(m_pi * (x0[i] - m_advx * time)) *
1026 cos(m_pi * (x1[i] - m_advy * time)) *
1027 cos(m_pi * (x2[i] - m_advz * time));
1028 }
1029
1030 outfield = u;
1031}

References m_advx, m_advy, m_advz, Nektar::SolverUtils::EquationSystem::m_fields, and Nektar::SolverUtils::MMFSystem::m_pi.

Referenced by v_EvaluateExactSolution(), and v_SetInitialConditions().

◆ v_DoSolve()

void Nektar::MMFAdvection::v_DoSolve ( void  )
overrideprotectedvirtual

Solves an unsteady problem.

Initialises the time integration scheme (as specified in the session file), and perform the time integration.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 196 of file MMFAdvection.cpp.

197{
198 ASSERTL0(m_intScheme != nullptr, "No time integration scheme.");
199
200 int i, nchk = 1;
201 int nvariables = 0;
202 int nfields = m_fields.size();
203 int nq = m_fields[0]->GetNpoints();
204
205 if (m_intVariables.empty())
206 {
207 for (i = 0; i < nfields; ++i)
208 {
209 m_intVariables.push_back(i);
210 }
211 nvariables = nfields;
212 }
213 else
214 {
215 nvariables = m_intVariables.size();
216 }
217
218 // Set up wrapper to fields data storage.
219 Array<OneD, Array<OneD, NekDouble>> fields(nvariables);
220 Array<OneD, Array<OneD, NekDouble>> tmp(nvariables);
221
222 // Order storage to list time-integrated fields first.
223 for (i = 0; i < nvariables; ++i)
224 {
225 fields[i] = m_fields[m_intVariables[i]]->GetPhys();
226 m_fields[m_intVariables[i]]->SetPhysState(false);
227 }
228
229 // Initialise time integration scheme
230 m_intScheme->InitializeScheme(m_timestep, fields, m_time, m_ode);
231
232 // Check uniqueness of checkpoint output
233 ASSERTL0((m_checktime == 0.0 && m_checksteps == 0) ||
234 (m_checktime > 0.0 && m_checksteps == 0) ||
235 (m_checktime == 0.0 && m_checksteps > 0),
236 "Only one of IO_CheckTime and IO_CheckSteps "
237 "should be set!");
238
239 LibUtilities::Timer timer;
240 bool doCheckTime = false;
241 int step = 0;
242 NekDouble intTime = 0.0;
243 NekDouble cpuTime = 0.0;
244 NekDouble elapsed = 0.0;
245
246 int Ntot, indx;
247 // Perform integration in time.
248 Ntot = m_checksteps ? m_steps / m_checksteps + 1 : 0;
249
250 Array<OneD, NekDouble> dMass(Ntot ? Ntot : 1);
251
252 Array<OneD, NekDouble> zeta(nq);
253 Array<OneD, Array<OneD, NekDouble>> fieldsprimitive(nvariables);
254 for (int i = 0; i < nvariables; ++i)
255 {
256 fieldsprimitive[i] = Array<OneD, NekDouble>(nq);
257 }
258
260 {
261 timer.Start();
262 fields = m_intScheme->TimeIntegrate(step, m_timestep);
263 timer.Stop();
264
266 elapsed = timer.TimePerTest(1);
267 intTime += elapsed;
268 cpuTime += elapsed;
269
270 // Write out status information
271 if (m_infosteps && !((step + 1) % m_infosteps) &&
272 m_session->GetComm()->GetRank() == 0)
273 {
274 std::cout << "Steps: " << std::setw(8) << std::left << step + 1
275 << " "
276 << "Time: " << std::setw(12) << std::left << m_time;
277
278 std::stringstream ss;
279 ss << cpuTime << "s";
280 std::cout << " CPU Time: " << std::setw(8) << std::left << ss.str()
281 << std::endl;
282
283 // Masss = h^*
284 indx = m_checksteps ? (step + 1) / m_checksteps : 0;
285 dMass[indx] =
286 (m_fields[0]->Integral(fields[0]) - m_Mass0) / m_Mass0;
287
288 std::cout << "dMass = " << std::setw(8) << std::left << dMass[indx]
289 << std::endl;
290
291 cpuTime = 0.0;
292 }
293
294 // Transform data into coefficient space
295 for (i = 0; i < nvariables; ++i)
296 {
297 m_fields[m_intVariables[i]]->SetPhys(fields[i]);
298 m_fields[m_intVariables[i]]->FwdTransLocalElmt(
299 fields[i], m_fields[m_intVariables[i]]->UpdateCoeffs());
300 m_fields[m_intVariables[i]]->SetPhysState(false);
301 }
302
303 // Write out checkpoint files
304 if ((m_checksteps && step && !((step + 1) % m_checksteps)) ||
305 doCheckTime)
306 {
307 Checkpoint_Output(nchk++);
308 doCheckTime = false;
309 }
310
311 // Step advance
312 ++step;
313 }
314
315 std::cout << "dMass = ";
316 for (i = 0; i < Ntot; ++i)
317 {
318 std::cout << dMass[i] << " , ";
319 }
320 std::cout << std::endl << std::endl;
321
322 // Print out summary statistics
323 if (m_session->GetComm()->GetRank() == 0)
324 {
325 if (m_cflSafetyFactor > 0.0)
326 {
327 std::cout << "CFL safety factor : " << m_cflSafetyFactor
328 << std::endl
329 << "CFL time-step : " << m_timestep << std::endl;
330 }
331
332 if (m_session->GetSolverInfo("Driver") != "SteadyState")
333 {
334 std::cout << "Time-integration : " << intTime << "s" << std::endl;
335 }
336 }
337
338 for (i = 0; i < nvariables; ++i)
339 {
340 m_fields[m_intVariables[i]]->SetPhys(fields[i]);
341 m_fields[m_intVariables[i]]->SetPhysState(true);
342 }
343
344 for (i = 0; i < nvariables; ++i)
345 {
346 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
347 m_fields[i]->UpdateCoeffs());
348 }
349}
NekDouble m_timestep
Time step size.
int m_infosteps
Number of time steps between outputting status information.
NekDouble m_time
Current time of simulation.
NekDouble m_fintime
Finish time of the simulation.
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
NekDouble m_checktime
Time between checkpoints.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
int m_steps
Number of steps to take.
int m_checksteps
Number of steps between checkpoints.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.
static const NekDouble kNekZeroTol

References ASSERTL0, Nektar::SolverUtils::EquationSystem::Checkpoint_Output(), Nektar::NekConstants::kNekZeroTol, Nektar::SolverUtils::UnsteadySystem::m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_checksteps, Nektar::SolverUtils::EquationSystem::m_checktime, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_fintime, Nektar::SolverUtils::EquationSystem::m_infosteps, Nektar::SolverUtils::UnsteadySystem::m_intScheme, Nektar::SolverUtils::UnsteadySystem::m_intVariables, m_Mass0, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_steps, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_timestep, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Nektar::LibUtilities::Timer::TimePerTest().

◆ v_EvaluateExactSolution()

void Nektar::MMFAdvection::v_EvaluateExactSolution ( unsigned int  field,
Array< OneD, NekDouble > &  outfield,
const NekDouble  time 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 1113 of file MMFAdvection.cpp.

1116{
1117 switch (m_TestType)
1118 {
1119 case eAdvectionBell:
1120 {
1121 AdvectionBellSphere(outfield);
1122 }
1123 break;
1124
1125 case eTestPlane:
1126 {
1127 Test2Dproblem(time, outfield);
1128 }
1129 break;
1130
1132 {
1133 AdvectionBellPlane(outfield);
1134 }
1135 break;
1136
1137 case eTestCube:
1138 {
1139 Test3Dproblem(time, outfield);
1140 }
1141 break;
1142
1143 default:
1144 break;
1145 }
1146}
@ eTestPlaneMassConsv
Definition: MMFAdvection.h:46
@ eAdvectionBell
Definition: MMFAdvection.h:48
void Test2Dproblem(const NekDouble time, Array< OneD, NekDouble > &outfield)
void Test3Dproblem(const NekDouble time, Array< OneD, NekDouble > &outfield)
void AdvectionBellSphere(Array< OneD, NekDouble > &outfield)
void AdvectionBellPlane(Array< OneD, NekDouble > &outfield)
@ eTestPlane
Definition: MMFDiffusion.h:48
@ eTestCube
Definition: MMFDiffusion.h:49

References AdvectionBellPlane(), AdvectionBellSphere(), eAdvectionBell, Nektar::eTestCube, Nektar::eTestPlane, eTestPlaneMassConsv, m_TestType, Test2Dproblem(), and Test3Dproblem().

◆ v_GenerateSummary()

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

Print Summary.

Reimplemented from Nektar::SolverUtils::MMFSystem.

Definition at line 1148 of file MMFAdvection.cpp.

1149{
1152 SolverUtils::AddSummaryItem(s, "Rotation Angle", m_RotAngle);
1153}
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Virtual function for generating summary information.
Definition: MMFSystem.cpp:2463
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47
const char *const TestTypeMap[]
Definition: MMFDiffusion.h:58

References Nektar::SolverUtils::AddSummaryItem(), m_RotAngle, m_TestType, Nektar::TestTypeMap, and Nektar::SolverUtils::MMFSystem::v_GenerateSummary().

◆ v_InitObject()

void Nektar::MMFAdvection::v_InitObject ( bool  DeclareFields = true)
overrideprotectedvirtual

Initialise the object.

Initialisation object for the unsteady linear advection equation.

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Definition at line 60 of file MMFAdvection.cpp.

61{
62 // Call to the initialisation object
63 UnsteadySystem::v_InitObject(DeclareFields);
64
65 int nq = m_fields[0]->GetNpoints();
66 int shapedim = m_fields[0]->GetShapeDimension();
67 Array<OneD, Array<OneD, NekDouble>> Anisotropy(shapedim);
68 for (int j = 0; j < shapedim; ++j)
69 {
70 Anisotropy[j] = Array<OneD, NekDouble>(nq, 1.0);
71 }
72
73 MMFSystem::MMFInitObject(Anisotropy);
74
75 // Define TestType
76 ASSERTL0(m_session->DefinesSolverInfo("TESTTYPE"),
77 "No TESTTYPE defined in session.");
78 std::string TestTypeStr = m_session->GetSolverInfo("TESTTYPE");
79 for (int i = 0; i < (int)SIZE_TestType; ++i)
80 {
81 if (boost::iequals(TestTypeMap[i], TestTypeStr))
82 {
84 break;
85 }
86 }
87
88 m_session->LoadParameter("Angular Frequency", m_waveFreq, m_pi);
89 m_session->LoadParameter("Rotational Angle", m_RotAngle, 0.0);
90 m_session->LoadParameter("Velocity Projection", m_VelProjection, 0);
91
92 // Read the advection velocities from session file
93 m_session->LoadParameter("advx", m_advx, 1.0);
94 m_session->LoadParameter("advy", m_advy, 1.0);
95 m_session->LoadParameter("advz", m_advz, 1.0);
96
97 std::vector<std::string> vel;
98 vel.push_back("Vx");
99 vel.push_back("Vy");
100 vel.push_back("Vz");
101
102 // Resize the advection velocities vector to dimension of the problem
103 vel.resize(m_spacedim);
104
105 // Store in the global variable m_velocity the advection velocities
106 m_velocity = Array<OneD, Array<OneD, NekDouble>>(m_spacedim);
107 for (int k = 0; k < m_spacedim; ++k)
108 {
109 m_velocity[k] = Array<OneD, NekDouble>(nq);
110 }
111
112 switch (m_surfaceType)
113 {
118 {
119 // true = project velocity onto the tangent plane
121 }
122 break;
123
126 {
127 GetFunction("AdvectionVelocity")->Evaluate(vel, m_velocity);
128 }
129 break;
130
131 default:
132 break;
133 }
134
135 std::cout << "|Velocity vector| = ( " << RootMeanSquare(m_velocity[0])
136 << " , " << RootMeanSquare(m_velocity[1]) << " , "
137 << RootMeanSquare(m_velocity[2]) << " ) " << std::endl;
138
139 // Define the normal velocity fields
140 if (m_fields[0]->GetTrace())
141 {
142 m_traceVn = Array<OneD, NekDouble>(GetTraceNpoints());
143 }
144
145 std::string advName;
146 std::string riemName;
147 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
150 m_advObject->SetFluxVector(&MMFAdvection::GetFluxVector, this);
151 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
153 riemName, m_session);
154 m_riemannSolver->SetScalar("Vn", &MMFAdvection::GetNormalVelocity, this);
155
156 m_advObject->SetRiemannSolver(m_riemannSolver);
157 m_advObject->InitObject(m_session, m_fields);
158
159 // Compute m_traceVn = n \cdot v
161
162 // Compute m_vellc = nabal a^j \cdot m_vel
164 std::cout << "m_vellc is generated with mag = " << AvgInt(m_vellc)
165 << std::endl;
166
167 // Compute vel \cdot MF
169
170 // Modify e^i as v^i e^i
171 for (int j = 0; j < m_shapedim; j++)
172 {
173 for (int k = 0; k < m_spacedim; k++)
174 {
175 Vmath::Vmul(nq, &m_veldotMF[j][0], 1, &m_movingframes[j][k * nq], 1,
176 &m_movingframes[j][k * nq], 1);
177 }
178 }
179
180 // Reflect it into m_ncdotMFFwd and Bwd
182
183 // If explicit it computes RHS and PROJECTION for the time integration
185 {
188 }
189 // Otherwise it gives an error (no implicit integration)
190 else
191 {
192 ASSERTL0(false, "Implicit unsteady Advection not set up.");
193 }
194}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point.
Array< OneD, Array< OneD, NekDouble > > m_veldotMF
Definition: MMFAdvection.h:98
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
void ComputeNablaCdotVelocity(Array< OneD, NekDouble > &vellc)
void ComputeveldotMF(Array< OneD, Array< OneD, NekDouble > > &veldotMF)
void EvaluateAdvectionVelocity(Array< OneD, Array< OneD, NekDouble > > &velocity)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Definition: MMFAdvection.h:86
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
Array< OneD, NekDouble > m_vellc
Definition: MMFAdvection.h:99
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
SOLVER_UTILS_EXPORT NekDouble AvgInt(const Array< OneD, const NekDouble > &inarray)
Definition: MMFSystem.cpp:2292
SOLVER_UTILS_EXPORT void ComputencdotMF()
Definition: MMFSystem.cpp:317
SOLVER_UTILS_EXPORT void MMFInitObject(const Array< OneD, const Array< OneD, NekDouble > > &Anisotropy, const int TangentXelem=-1)
Definition: MMFSystem.cpp:51
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
RiemannSolverFactory & GetRiemannSolverFactory()
@ SIZE_TestType
Length of enum list.
Definition: MMFDiffusion.h:55

References ASSERTL0, Nektar::SolverUtils::MMFSystem::AvgInt(), ComputeNablaCdotVelocity(), Nektar::SolverUtils::MMFSystem::ComputencdotMF(), ComputeveldotMF(), Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoOdeProjection(), DoOdeRhs(), Nektar::SolverUtils::eCube, Nektar::SolverUtils::eIrregular, Nektar::SolverUtils::eNonconvex, Nektar::SolverUtils::ePlane, Nektar::SolverUtils::eSphere, Nektar::SolverUtils::eTRSphere, EvaluateAdvectionVelocity(), Nektar::SolverUtils::GetAdvectionFactory(), GetFluxVector(), Nektar::SolverUtils::EquationSystem::GetFunction(), GetNormalVelocity(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::AdvectionSystem::m_advObject, m_advx, m_advy, m_advz, Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::MMFSystem::m_movingframes, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::MMFSystem::m_pi, m_riemannSolver, m_RotAngle, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::MMFSystem::m_shapedim, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::MMFSystem::m_surfaceType, m_TestType, m_traceVn, m_veldotMF, m_vellc, m_velocity, m_VelProjection, m_waveFreq, Nektar::SolverUtils::MMFSystem::MMFInitObject(), Nektar::SolverUtils::MMFSystem::RootMeanSquare(), Nektar::SIZE_TestType, Nektar::TestTypeMap, Nektar::SolverUtils::UnsteadySystem::v_InitObject(), and Vmath::Vmul().

◆ v_SetInitialConditions()

void Nektar::MMFAdvection::v_SetInitialConditions ( const NekDouble  initialtime,
bool  dumpInitialConditions,
const int  domain 
)
overrideprotectedvirtual

Set the physical fields based on a restart file, or a function describing the initial condition given in the session.

Parameters
initialtimeTime at which to evaluate the function.
dumpInitialConditionsWrite the initial condition to file?

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 812 of file MMFAdvection.cpp.

815{
816 int nq = m_fields[0]->GetNpoints();
817
818 Array<OneD, NekDouble> u(nq);
819
820 switch (m_TestType)
821 {
822 case eAdvectionBell:
823 {
825 m_fields[0]->SetPhys(u);
826
827 m_Mass0 = m_fields[0]->Integral(u);
828
829 // forward transform to fill the modal coeffs
830 for (int i = 0; i < m_fields.size(); ++i)
831 {
832 m_fields[i]->SetPhysState(true);
833 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
834 m_fields[i]->UpdateCoeffs());
835 }
836 }
837 break;
838
839 case eTestPlane:
840 {
841 Test2Dproblem(initialtime, u);
842 m_fields[0]->SetPhys(u);
843
844 m_Mass0 = m_fields[0]->Integral(u);
845
846 // forward transform to fill the modal coeffs
847 for (int i = 0; i < m_fields.size(); ++i)
848 {
849 m_fields[i]->SetPhysState(true);
850 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
851 m_fields[i]->UpdateCoeffs());
852 }
853 }
854 break;
855
857 {
859 m_fields[0]->SetPhys(u);
860
861 m_Mass0 = m_fields[0]->Integral(u);
862 std::cout << "m_Mass0 = " << m_Mass0 << std::endl;
863
864 // forward transform to fill the modal coeffs
865 for (int i = 0; i < m_fields.size(); ++i)
866 {
867 m_fields[i]->SetPhysState(true);
868 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
869 m_fields[i]->UpdateCoeffs());
870 }
871 }
872 break;
873
874 case eTestCube:
875 {
876 Test3Dproblem(initialtime, u);
877 m_fields[0]->SetPhys(u);
878
879 // forward transform to fill the modal coeffs
880 for (int i = 0; i < m_fields.size(); ++i)
881 {
882 m_fields[i]->SetPhysState(true);
883 m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
884 m_fields[i]->UpdateCoeffs());
885 }
886 }
887 break;
888
889 default:
890 break;
891 }
892
893 if (dumpInitialConditions)
894 {
895 // dump initial conditions to file
896 std::string outname = m_sessionName + "_initial.chk";
897 WriteFld(outname);
898 }
899}
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
std::string m_sessionName
Name of the session.

References AdvectionBellPlane(), AdvectionBellSphere(), eAdvectionBell, Nektar::eTestCube, Nektar::eTestPlane, eTestPlaneMassConsv, Nektar::SolverUtils::EquationSystem::m_fields, m_Mass0, Nektar::SolverUtils::EquationSystem::m_sessionName, m_TestType, Test2Dproblem(), Test3Dproblem(), and Nektar::SolverUtils::EquationSystem::WriteFld().

◆ WeakDGDirectionalAdvection()

void Nektar::MMFAdvection::WeakDGDirectionalAdvection ( const Array< OneD, Array< OneD, NekDouble > > &  InField,
Array< OneD, Array< OneD, NekDouble > > &  OutField 
)
protected

Definition at line 525 of file MMFAdvection.cpp.

528{
529 int i, j;
530 int ncoeffs = GetNcoeffs();
531 int nTracePointsTot = GetTraceNpoints();
532 int nvariables = m_fields.size();
533
534 Array<OneD, Array<OneD, NekDouble>> physfield(nvariables);
535
536 // Get the variables in physical space
537 // already in physical space
538 for (i = 0; i < nvariables; ++i)
539 {
540 physfield[i] = InField[i];
541 }
542
543 Array<OneD, Array<OneD, NekDouble>> WeakDeriv(m_shapedim);
544 for (i = 0; i < nvariables; ++i)
545 {
546 for (j = 0; j < m_shapedim; ++j)
547 {
548 WeakDeriv[j] = Array<OneD, NekDouble>(ncoeffs, 0.0);
549
550 // Directional derivation with respect to the j'th moving frame
551 // tmp[j] = \nabla \physfield[i] \cdot \mathbf{e}^j
552 // Implemented at TriExp::v_IProductWRTDirectionalDerivBase_SumFa
553 m_fields[0]->IProductWRTDirectionalDerivBase(
554 m_movingframes[j], physfield[i], WeakDeriv[j]);
555 }
556
557 // Get the numerical flux and add to the modal coeffs
558 // if the NumericalFluxs function already includes the
559 // normal in the output
560
561 Array<OneD, NekDouble> Fwd(nTracePointsTot);
562 Array<OneD, NekDouble> Bwd(nTracePointsTot);
563
564 Array<OneD, NekDouble> flux(nTracePointsTot, 0.0);
565 Array<OneD, NekDouble> fluxFwd(nTracePointsTot);
566 Array<OneD, NekDouble> fluxBwd(nTracePointsTot);
567
568 // Evaluate numerical flux in physical space which may in
569 // general couple all component of vectors
570 m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
571
572 // evaulate upwinded m_fields[i]
573 m_fields[i]->GetTrace()->Upwind(m_traceVn, Fwd, Bwd, flux);
574
575 OutField[i] = Array<OneD, NekDouble>(ncoeffs, 0.0);
576 for (j = 0; j < m_shapedim; ++j)
577 {
578 // calculate numflux = (n \cdot MF)*flux
579 Vmath::Vmul(nTracePointsTot, &flux[0], 1, &m_ncdotMFFwd[j][0], 1,
580 &fluxFwd[0], 1);
581 Vmath::Vmul(nTracePointsTot, &flux[0], 1, &m_ncdotMFBwd[j][0], 1,
582 &fluxBwd[0], 1);
583 Vmath::Neg(ncoeffs, WeakDeriv[j], 1);
584
585 // FwdBwdtegral because generallize (N \cdot MF)_{FWD} \neq -(N
586 // \cdot MF)_{BWD}
587 m_fields[i]->AddFwdBwdTraceIntegral(fluxFwd, fluxBwd, WeakDeriv[j]);
588 m_fields[i]->SetPhysState(false);
589
590 Vmath::Vadd(ncoeffs, &WeakDeriv[j][0], 1, &OutField[i][0], 1,
591 &OutField[i][0], 1);
592 }
593 }
594}
Array< OneD, Array< OneD, NekDouble > > m_ncdotMFBwd
Definition: MMFSystem.h:187
Array< OneD, Array< OneD, NekDouble > > m_ncdotMFFwd
Definition: MMFSystem.h:186

References Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::MMFSystem::m_movingframes, Nektar::SolverUtils::MMFSystem::m_ncdotMFBwd, Nektar::SolverUtils::MMFSystem::m_ncdotMFFwd, Nektar::SolverUtils::MMFSystem::m_shapedim, m_traceVn, Vmath::Neg(), Vmath::Vadd(), and Vmath::Vmul().

Referenced by DoOdeRhs().

Friends And Related Function Documentation

◆ MemoryManager< MMFAdvection >

friend class MemoryManager< MMFAdvection >
friend

Definition at line 52 of file MMFAdvection.h.

Member Data Documentation

◆ className

std::string Nektar::MMFAdvection::className
static
Initial value:
=
"MMFAdvection", MMFAdvection::create, "MMFAdvection equation.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: MMFAdvection.h:68
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 78 of file MMFAdvection.h.

◆ m_advx

NekDouble Nektar::MMFAdvection::m_advx
protected

Definition at line 88 of file MMFAdvection.h.

Referenced by Test2Dproblem(), Test3Dproblem(), and v_InitObject().

◆ m_advy

NekDouble Nektar::MMFAdvection::m_advy
protected

Definition at line 88 of file MMFAdvection.h.

Referenced by Test2Dproblem(), Test3Dproblem(), and v_InitObject().

◆ m_advz

NekDouble Nektar::MMFAdvection::m_advz
protected

Definition at line 88 of file MMFAdvection.h.

Referenced by Test3Dproblem(), and v_InitObject().

◆ m_Mass0

NekDouble Nektar::MMFAdvection::m_Mass0
protected

Definition at line 91 of file MMFAdvection.h.

Referenced by v_DoSolve(), and v_SetInitialConditions().

◆ m_riemannSolver

SolverUtils::RiemannSolverSharedPtr Nektar::MMFAdvection::m_riemannSolver
protected

Definition at line 86 of file MMFAdvection.h.

Referenced by v_InitObject().

◆ m_RotAngle

NekDouble Nektar::MMFAdvection::m_RotAngle
protected

Definition at line 89 of file MMFAdvection.h.

Referenced by EvaluateAdvectionVelocity(), v_GenerateSummary(), and v_InitObject().

◆ m_TestType

TestType Nektar::MMFAdvection::m_TestType

◆ m_traceVn

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

Definition at line 96 of file MMFAdvection.h.

Referenced by GetNormalVelocity(), v_InitObject(), and WeakDGDirectionalAdvection().

◆ m_veldotMF

Array<OneD, Array<OneD, NekDouble> > Nektar::MMFAdvection::m_veldotMF
protected

Definition at line 98 of file MMFAdvection.h.

Referenced by v_InitObject().

◆ m_vellc

Array<OneD, NekDouble> Nektar::MMFAdvection::m_vellc
protected

Definition at line 99 of file MMFAdvection.h.

Referenced by v_InitObject().

◆ m_velocity

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

◆ m_VelProjection

int Nektar::MMFAdvection::m_VelProjection
protected

Definition at line 92 of file MMFAdvection.h.

Referenced by EvaluateAdvectionVelocity(), and v_InitObject().

◆ m_waveFreq

NekDouble Nektar::MMFAdvection::m_waveFreq
protected

Definition at line 89 of file MMFAdvection.h.

Referenced by EvaluateAdvectionVelocity(), and v_InitObject().