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

#include <VCSMapping.h>

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

Public Member Functions

 VCSMapping (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor. More...
 
void ApplyIncNSMappingForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual ~VCSMapping ()
 
virtual void v_InitObject ()
 Init object for UnsteadySystem class. More...
 
- Public Member Functions inherited from Nektar::VelocityCorrectionScheme
 VelocityCorrectionScheme (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor. More...
 
virtual ~VelocityCorrectionScheme ()
 
void SetUpPressureForcing (const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
void SetUpViscousForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
void SolvePressure (const Array< OneD, NekDouble > &Forcing)
 
void SolveViscous (const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
 
void SolveUnsteadyStokesSystem (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
 
void EvaluateAdvection_SetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
- Public Member Functions inherited from Nektar::IncNavierStokes
virtual ~IncNavierStokes ()
 
virtual void v_GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
virtual void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
int GetNConvectiveFields (void)
 
Array< OneD, int > & GetVelocity (void)
 
Array< OneD, NekDoubleGetElmtCFLVals (void)
 
NekDouble GetCFLEstimate (int &elmtid)
 
void AddForcing (const SolverUtils::ForcingSharedPtr &pForce)
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual SOLVER_UTILS_EXPORT ~AdvectionSystem ()
 
AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
boost::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 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 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. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
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. More...
 
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $. More...
 
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) $. More...
 
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. More...
 
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. More...
 
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. 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 ScanForHistoryPoints ()
 Builds map of which element holds each history point. More...
 
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to 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::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const
Array< OneD, int > 
GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const
std::string 
GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int 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 SetSteps (const int steps)
 
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 std::string &s1, const std::string &s2)
 Perform a case-insensitive string comparison. More...
 
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 void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Static Public Member Functions

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

Static Public Attributes

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

Protected Member Functions

virtual void v_DoInitialise (void)
 Sets up initial conditions. More...
 
virtual void v_SetUpPressureForcing (const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
virtual void v_SetUpViscousForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
virtual void v_SolvePressure (const Array< OneD, NekDouble > &Forcing)
 
virtual void v_SolveViscous (const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
 
virtual void v_EvaluateAdvection_SetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
- Protected Member Functions inherited from Nektar::VelocityCorrectionScheme
virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
 Print a summary of time stepping parameters. More...
 
virtual void v_TransCoeffToPhys (void)
 Virtual function for transformation to physical space. More...
 
virtual void v_TransPhysToCoeff (void)
 Virtual function for transformation to coefficient space. More...
 
virtual Array< OneD, bool > v_GetSystemSingularChecks ()
 
virtual int v_GetForceDimension ()
 
virtual bool v_RequireFwdTrans ()
 
virtual std::string v_GetExtrapolateStr (void)
 
virtual std::string v_GetSubSteppingExtrapolateStr (const std::string &instr)
 
- Protected Member Functions inherited from Nektar::IncNavierStokes
 IncNavierStokes (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor. More...
 
EquationType GetEquationType (void)
 
void EvaluateAdvectionTerms (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void WriteModalEnergy (void)
 
void SetBoundaryConditions (NekDouble time)
 time dependent boundary conditions updating More...
 
void SetRadiationBoundaryForcing (int fieldid)
 Set Radiation forcing term. More...
 
void SetZeroNormalVelocity ()
 Set Normal Velocity Component to Zero. More...
 
void SetWomersleyBoundary (const int fldid, const int bndid)
 Set Womersley Profile if specified. More...
 
void SetUpWomersley (const int bndid, std::string womstr)
 Set Up Womersley details. More...
 
bool CalcSteadyState (void)
 evaluate steady state More...
 
virtual
MultiRegions::ExpListSharedPtr 
v_GetPressure ()
 
virtual bool v_PreIntegrate (int step)
 
virtual bool v_PostIntegrate (int step)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file. More...
 
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. More...
 
virtual SOLVER_UTILS_EXPORT bool v_SteadyStateCheck (int step)
 
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...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const std::string &s1, const std::string &s2)
 
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_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 EvaluateFunctionExp (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void EvaluateFunctionFld (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void EvaluateFunctionPts (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void LoadPts (std::string funcFilename, std::string filename, Nektar::LibUtilities::PtsFieldSharedPtr &outPts)
 
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 void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

GlobalMapping::MappingSharedPtr m_mapping
 
bool m_verbose
 
bool m_implicitPressure
 
bool m_implicitViscous
 
bool m_neglectViscous
 
NekDouble m_pressureTolerance
 
NekDouble m_viscousTolerance
 
NekDouble m_pressureRelaxation
 
NekDouble m_viscousRelaxation
 
Array< OneD, Array< OneD,
NekDouble > > 
m_gradP
 
- Protected Attributes inherited from Nektar::VelocityCorrectionScheme
bool m_useHomo1DSpecVanVisc
 bool to identify if spectral vanishing viscosity is active. More...
 
bool m_useSpecVanVisc
 bool to identify if spectral vanishing viscosity is active. More...
 
NekDouble m_sVVCutoffRatio
 cutt off ratio from which to start decayhing modes More...
 
NekDouble m_sVVDiffCoeff
 Diffusion coefficient of SVV modes. More...
 
Array< OneD, NekDoublem_diffCoeff
 Diffusion coefficients (will be kinvis for velocities) More...
 
StdRegions::VarCoeffMap m_varCoeffLap
 Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise. More...
 
Array< OneD, Array< OneD,
NekDouble > > 
m_F
 
- Protected Attributes inherited from Nektar::IncNavierStokes
ExtrapolateSharedPtr m_extrapolation
 
std::ofstream m_mdlFile
 modal energy file More...
 
bool m_SmoothAdvection
 bool to identify if advection term smoothing is requested More...
 
std::vector
< SolverUtils::ForcingSharedPtr
m_forcing
 Forcing terms. More...
 
int m_nConvectiveFields
 Number of fields to be convected;. More...
 
Array< OneD, int > m_velocity
 int which identifies which components of m_fields contains the velocity (u,v,w); More...
 
MultiRegions::ExpListSharedPtr m_pressure
 Pointer to field holding pressure field. More...
 
NekDouble m_kinvis
 Kinematic viscosity. More...
 
int m_energysteps
 dump energy to file at steps time More...
 
int m_cflsteps
 dump cfl estimate 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...
 
EquationType m_equationType
 equation type; More...
 
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
 Mapping from BCs to Elmt IDs. More...
 
Array< OneD, Array< OneD, int > > m_fieldsBCToTraceID
 Mapping from BCs to Elmt Edge IDs. More...
 
Array< OneD, Array< OneD,
NekDouble > > 
m_fieldsRadiationFactor
 RHS Factor for Radiation Condition. More...
 
int m_intSteps
 Number of time integration steps AND Order of extrapolation for pressure boundary conditions. More...
 
std::map< int,
WomersleyParamsSharedPtr
m_womersleyParams
 Womersley parameters if required. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem
SolverUtils::AdvectionSharedPtr m_advObject
 Advection term. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::vector< int > m_intVariables
 
std::vector< FilterSharedPtrm_filters
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
std::map< std::string,
FieldUtils::Interpolator
m_interpolators
 Map of interpolator objects. More...
 
std::map< std::string,
std::pair< std::string,
LibUtilities::PtsFieldSharedPtr > > 
m_loadedPtsFields
 pts fields we already read from disk: {funcFilename: (filename, ptsfield)} More...
 
std::map< std::string,
std::pair< std::string,
loadedFldField > > 
m_loadedFldFields
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 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...
 
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_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, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq 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...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 

Private Member Functions

void MappingAdvectionCorrection (const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void MappingAccelerationCorrection (const Array< OneD, const Array< OneD, NekDouble > > &vel, const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void MappingPressureCorrection (Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void MappingViscousCorrection (const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
 

Private Attributes

Array< OneD, Array< OneD,
NekDouble > > 
m_presForcingCorrection
 

Additional Inherited Members

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

Detailed Description

Definition at line 44 of file VCSMapping.h.

Constructor & Destructor Documentation

Nektar::VCSMapping::VCSMapping ( const LibUtilities::SessionReaderSharedPtr pSession)

Constructor.

Constructor. Creates ...

Parameters
param

Definition at line 57 of file VCSMapping.cpp.

59  : UnsteadySystem(pSession),
60  VelocityCorrectionScheme(pSession)
61  {
62 
63  }
VelocityCorrectionScheme(const LibUtilities::SessionReaderSharedPtr &pSession)
Constructor.
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises UnsteadySystem class members.
Nektar::VCSMapping::~VCSMapping ( void  )
virtual

Destructor

Definition at line 143 of file VCSMapping.cpp.

144  {
145  }

Member Function Documentation

void Nektar::VCSMapping::ApplyIncNSMappingForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)

Explicit terms of the mapping

Definition at line 750 of file VCSMapping.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields, m_implicitPressure, m_implicitViscous, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, m_neglectViscous, MappingAccelerationCorrection(), MappingAdvectionCorrection(), MappingPressureCorrection(), MappingViscousCorrection(), Vmath::Smul(), Vmath::Vadd(), and Vmath::Vcopy().

Referenced by v_EvaluateAdvection_SetPressureBCs().

753  {
754  int physTot = m_fields[0]->GetTotPoints();
755  Array<OneD, Array<OneD, NekDouble> > vel(m_nConvectiveFields);
756  Array<OneD, Array<OneD, NekDouble> > velPhys(m_nConvectiveFields);
757  Array<OneD, Array<OneD, NekDouble> > Forcing(m_nConvectiveFields);
758  Array<OneD, Array<OneD, NekDouble> > tmp(m_nConvectiveFields);
759  for (int i = 0; i < m_nConvectiveFields; ++i)
760  {
761  velPhys[i] = Array<OneD, NekDouble> (physTot, 0.0);
762  Forcing[i] = Array<OneD, NekDouble> (physTot, 0.0);
763  tmp[i] = Array<OneD, NekDouble> (physTot, 0.0);
764  }
765 
766  // Get fields and store velocity in wavespace and physical space
767  if(m_fields[0]->GetWaveSpace())
768  {
769  for (int i = 0; i < m_nConvectiveFields; ++i)
770  {
771  vel[i] = inarray[i];
772  m_fields[0]->HomogeneousBwdTrans(vel[i],velPhys[i]);
773  }
774  }
775  else
776  {
777  for (int i = 0; i < m_nConvectiveFields; ++i)
778  {
779  vel[i] = inarray[i];
780  Vmath::Vcopy(physTot, inarray[i], 1, velPhys[i], 1);
781  }
782  }
783 
784  //Advection contribution
785  MappingAdvectionCorrection(velPhys, Forcing);
786 
787  // Time-derivative contribution
788  if ( m_mapping->IsTimeDependent() )
789  {
790  MappingAccelerationCorrection(vel, velPhys, tmp);
791  for (int i = 0; i < m_nConvectiveFields; ++i)
792  {
793  Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
794  }
795  }
796 
797  // Pressure contribution
798  if (!m_implicitPressure)
799  {
801  for (int i = 0; i < m_nConvectiveFields; ++i)
802  {
803  Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
804  }
805  }
806  // Viscous contribution
807  if ( (!m_implicitViscous) && (!m_neglectViscous))
808  {
809  MappingViscousCorrection(velPhys, tmp);
810  for (int i = 0; i < m_nConvectiveFields; ++i)
811  {
812  Vmath::Smul(physTot, m_kinvis, tmp[i], 1, tmp[i], 1);
813  Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
814  }
815  }
816 
817  // If necessary, transform to wavespace
818  if(m_fields[0]->GetWaveSpace())
819  {
820  for (int i = 0; i < m_nConvectiveFields; ++i)
821  {
822  m_fields[0]->HomogeneousFwdTrans(Forcing[i],Forcing[i]);
823  }
824  }
825 
826  // Add to outarray
827  for (int i = 0; i < m_nConvectiveFields; ++i)
828  {
829  Vmath::Vadd(physTot, outarray[i], 1, Forcing[i], 1, outarray[i], 1);
830  }
831  }
void MappingAccelerationCorrection(const Array< OneD, const Array< OneD, NekDouble > > &vel, const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:858
NekDouble m_kinvis
Kinematic viscosity.
void MappingAdvectionCorrection(const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:833
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
int m_nConvectiveFields
Number of fields to be convected;.
void MappingViscousCorrection(const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:964
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:213
void MappingPressureCorrection(Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:937
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:299
static SolverUtils::EquationSystemSharedPtr Nektar::VCSMapping::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 49 of file VCSMapping.h.

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

50  {
53  AllocateSharedPtr(pSession);
54  p->InitObject();
55  return p;
56  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
void Nektar::VCSMapping::MappingAccelerationCorrection ( const Array< OneD, const Array< OneD, NekDouble > > &  vel,
const Array< OneD, const Array< OneD, NekDouble > > &  velPhys,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 858 of file VCSMapping.cpp.

References Nektar::MultiRegions::DirCartesianMap, Nektar::SolverUtils::EquationSystem::m_fields, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Vmath::Neg(), Vmath::Vadd(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by ApplyIncNSMappingForcing().

862  {
863  int physTot = m_fields[0]->GetTotPoints();
864  int nvel = m_nConvectiveFields;
865 
866  Array<OneD, Array<OneD, NekDouble> > wk(nvel*nvel);
867  Array<OneD, Array<OneD, NekDouble> > tmp(nvel);
868  Array<OneD, Array<OneD, NekDouble> > coordVel(nvel);
869  for (int i = 0; i< nvel; i++)
870  {
871  tmp[i] = Array<OneD, NekDouble> (physTot, 0.0);
872  coordVel[i] = Array<OneD, NekDouble> (physTot, 0.0);
873  }
874  // Get coordinates velocity in transformed system
875  m_mapping->GetCoordVelocity(tmp);
876  m_mapping->ContravarFromCartesian(tmp, coordVel);
877 
878  // Calculate first term: U^j u^i,j = U^j (du^i/dx^j + {i,kj}u^k)
879  m_mapping->ApplyChristoffelContravar(velPhys, wk);
880  for (int i=0; i< nvel; i++)
881  {
882  Vmath::Zero(physTot,outarray[i],1);
883 
884  m_fields[0]->PhysDeriv(velPhys[i], tmp[0], tmp[1]);
885  for (int j=0; j< nvel; j++)
886  {
887  if (j == 2)
888  {
889  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j],
890  vel[i], tmp[2]);
891  if (m_fields[0]->GetWaveSpace())
892  {
893  m_fields[0]->HomogeneousBwdTrans(tmp[2],tmp[2]);
894  }
895  }
896 
897  Vmath::Vadd(physTot,wk[i*nvel+j],1,tmp[j],1,
898  wk[i*nvel+j], 1);
899 
900  Vmath::Vvtvp(physTot, coordVel[j], 1, wk[i*nvel+j], 1,
901  outarray[i], 1, outarray[i], 1);
902  }
903  }
904 
905  // Set wavespace to false and store current value
906  bool wavespace = m_fields[0]->GetWaveSpace();
907  m_fields[0]->SetWaveSpace(false);
908 
909  // Add -u^j U^i,j
910  m_mapping->ApplyChristoffelContravar(coordVel, wk);
911  for (int i=0; i< nvel; i++)
912  {
913  if(nvel == 2)
914  {
915  m_fields[0]->PhysDeriv(coordVel[i], tmp[0], tmp[1]);
916  }
917  else
918  {
919  m_fields[0]->PhysDeriv(coordVel[i], tmp[0], tmp[1], tmp[2]);
920  }
921 
922  for (int j=0; j< nvel; j++)
923  {
924  Vmath::Vadd(physTot,wk[i*nvel+j],1,tmp[j],1,
925  wk[i*nvel+j], 1);
926  Vmath::Neg(physTot, wk[i*nvel+j], 1);
927 
928  Vmath::Vvtvp(physTot, velPhys[j], 1, wk[i*nvel+j], 1,
929  outarray[i], 1, outarray[i], 1);
930  }
931  }
932 
933  // Restore value of wavespace
934  m_fields[0]->SetWaveSpace(wavespace);
935  }
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:442
int m_nConvectiveFields
Number of fields to be convected;.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:396
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:299
void Nektar::VCSMapping::MappingAdvectionCorrection ( const Array< OneD, const Array< OneD, NekDouble > > &  velPhys,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 833 of file VCSMapping.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Vmath::Neg(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by ApplyIncNSMappingForcing().

836  {
837  int physTot = m_fields[0]->GetTotPoints();
838  int nvel = m_nConvectiveFields;
839 
840  Array<OneD, Array<OneD, NekDouble> > wk(nvel*nvel);
841 
842  // Apply Christoffel symbols to obtain {i,kj}vel(k)
843  m_mapping->ApplyChristoffelContravar(velPhys, wk);
844 
845  // Calculate correction -U^j*{i,kj}vel(k)
846  for (int i = 0; i< nvel; i++)
847  {
848  Vmath::Zero(physTot,outarray[i],1);
849  for (int j = 0; j< nvel; j++)
850  {
851  Vmath::Vvtvp(physTot,wk[i*nvel+j],1,velPhys[j],1,
852  outarray[i],1,outarray[i],1);
853  }
854  Vmath::Neg(physTot, outarray[i], 1);
855  }
856  }
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:442
int m_nConvectiveFields
Number of fields to be convected;.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:396
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
void Nektar::VCSMapping::MappingPressureCorrection ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
private

Definition at line 937 of file VCSMapping.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields, m_gradP, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Vmath::Vdiv(), and Vmath::Vsub().

Referenced by ApplyIncNSMappingForcing().

939  {
940  int physTot = m_fields[0]->GetTotPoints();
941  int nvel = m_nConvectiveFields;
942 
943  // Calculate g^(ij)p_(,j)
944  m_mapping->RaiseIndex(m_gradP, outarray);
945 
946  // Calculate correction = (nabla p)/J - g^(ij)p_,j
947  // (Jac is not required if it is constant)
948  if ( !m_mapping->HasConstantJacobian())
949  {
950  Array<OneD, NekDouble> Jac(physTot, 0.0);
951  m_mapping->GetJacobian(Jac);
952  for(int i = 0; i < nvel; ++i)
953  {
954  Vmath::Vdiv(physTot, m_gradP[i], 1, Jac, 1, m_gradP[i], 1);
955  }
956  }
957  for(int i = 0; i < nvel; ++i)
958  {
959  Vmath::Vsub(physTot, m_gradP[i], 1,outarray[i], 1,
960  outarray[i],1);
961  }
962  }
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.cpp:241
int m_nConvectiveFields
Number of fields to be convected;.
Array< OneD, Array< OneD, NekDouble > > m_gradP
Definition: VCSMapping.h:93
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:343
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Nektar::VCSMapping::MappingViscousCorrection ( const Array< OneD, const Array< OneD, NekDouble > > &  velPhys,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 964 of file VCSMapping.cpp.

References m_mapping.

Referenced by ApplyIncNSMappingForcing().

967  {
968  // L(U) - 1.0*d^2(u^i)/dx^jdx^j
969  m_mapping->VelocityLaplacian(velPhys, outarray, 1.0);
970  }
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
void Nektar::VCSMapping::v_DoInitialise ( void  )
protectedvirtual

Sets up initial conditions.

Sets the initial conditions.

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 147 of file VCSMapping.cpp.

References Nektar::MultiRegions::DirCartesianMap, Nektar::VelocityCorrectionScheme::m_F, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, m_gradP, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::m_pressure, Nektar::SolverUtils::EquationSystem::m_timestep, and Nektar::SolverUtils::UnsteadySystem::v_DoInitialise().

148  {
150 
151  // Set up Field Meta Data for output files
152  m_fieldMetaDataMap["Kinvis"] =
153  boost::lexical_cast<std::string>(m_kinvis);
154  m_fieldMetaDataMap["TimeStep"] =
155  boost::lexical_cast<std::string>(m_timestep);
156 
157  // Correct Dirichlet boundary conditions to account for mapping
158  m_mapping->UpdateBCs(0.0);
159  //
160  m_F = Array<OneD, Array< OneD, NekDouble> > (m_nConvectiveFields);
161  for(int i = 0; i < m_nConvectiveFields; ++i)
162  {
163  m_fields[i]->LocalToGlobal();
164  m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
165  m_fields[i]->GlobalToLocal();
166  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
167  m_fields[i]->UpdatePhys());
168  m_F[i] = Array< OneD, NekDouble> (m_fields[0]->GetTotPoints(), 0.0);
169  }
170 
171  // Initialise m_gradP
172  int physTot = m_fields[0]->GetTotPoints();
173  m_gradP = Array<OneD, Array<OneD, NekDouble> >(m_nConvectiveFields);
174  for(int i = 0; i < m_nConvectiveFields; ++i)
175  {
176  m_gradP[i] = Array<OneD, NekDouble>(physTot,0.0);
178  m_pressure->GetPhys(),
179  m_gradP[i]);
180  if(m_pressure->GetWaveSpace())
181  {
182  m_pressure->HomogeneousBwdTrans(m_gradP[i],m_gradP[i]);
183  }
184  }
185  }
NekDouble m_kinvis
Kinematic viscosity.
NekDouble m_timestep
Time step size.
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
int m_nConvectiveFields
Number of fields to be convected;.
Array< OneD, Array< OneD, NekDouble > > m_gradP
Definition: VCSMapping.h:93
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
Array< OneD, Array< OneD, NekDouble > > m_F
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise()
Sets up initial conditions.
void Nektar::VCSMapping::v_EvaluateAdvection_SetPressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protectedvirtual

Explicit part of the method - Advection, Forcing + HOPBCs

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 190 of file VCSMapping.cpp.

References ApplyIncNSMappingForcing(), Nektar::IncNavierStokes::EvaluateAdvectionTerms(), Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_forcing, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::m_pressure, Nektar::IncNavierStokes::m_SmoothAdvection, and Nektar::SolverUtils::EquationSystem::m_timestep.

194  {
195  EvaluateAdvectionTerms(inarray, outarray);
196 
197  // Smooth advection
199  {
200  for(int i = 0; i < m_nConvectiveFields; ++i)
201  {
202  m_pressure->SmoothField(outarray[i]);
203  }
204  }
205 
206  // Add forcing terms
207  std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x;
208  for (x = m_forcing.begin(); x != m_forcing.end(); ++x)
209  {
210  (*x)->Apply(m_fields, inarray, outarray, time);
211  }
212 
213  // Add mapping terms
214  ApplyIncNSMappingForcing( inarray, outarray);
215 
216  // Calculate High-Order pressure boundary conditions
217  m_extrapolation->EvaluatePressureBCs(inarray,outarray,m_kinvis);
218 
219  // Update mapping and deal with Dirichlet boundary conditions
220  if (m_mapping->IsTimeDependent())
221  {
222  if (m_mapping->IsFromFunction())
223  {
224  // If the transformation is explicitly defined, update it here
225  // Otherwise, it will be done somewhere else (ForcingMovingBody)
226  m_mapping->UpdateMapping(time+m_timestep);
227  }
228  m_mapping->UpdateBCs(time+m_timestep);
229  }
230  }
void ApplyIncNSMappingForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:750
NekDouble m_kinvis
Kinematic viscosity.
NekDouble m_timestep
Time step size.
ExtrapolateSharedPtr m_extrapolation
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
void EvaluateAdvectionTerms(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
int m_nConvectiveFields
Number of fields to be convected;.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
void Nektar::VCSMapping::v_InitObject ( )
virtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 65 of file VCSMapping.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::LibUtilities::eIMEXOrder1, Nektar::LibUtilities::eIMEXOrder2, Nektar::LibUtilities::eIMEXOrder3, Nektar::GetExtrapolateFactory(), Nektar::GlobalMapping::Mapping::Load(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, m_implicitPressure, m_implicitViscous, Nektar::SolverUtils::UnsteadySystem::m_intScheme, m_mapping, m_neglectViscous, m_presForcingCorrection, Nektar::IncNavierStokes::m_pressure, m_pressureRelaxation, m_pressureTolerance, Nektar::SolverUtils::EquationSystem::m_session, Nektar::IncNavierStokes::m_velocity, m_verbose, m_viscousRelaxation, m_viscousTolerance, and Nektar::VelocityCorrectionScheme::v_InitObject().

66  {
68 
71  "Could not create mapping in VCSMapping.");
72 
73  std::string vExtrapolation = "Mapping";
75  vExtrapolation,
76  m_session,
77  m_fields,
78  m_pressure,
79  m_velocity,
80  m_advObject);
81  m_extrapolation->SubSteppingTimeIntegration(
82  m_intScheme->GetIntegrationMethod(), m_intScheme);
83  m_extrapolation->GenerateHOPBCMap(m_session);
84 
85  // Storage to extrapolate pressure forcing
86  int physTot = m_fields[0]->GetTotPoints();
87  int intSteps = 1;
88  int intMethod = m_intScheme->GetIntegrationMethod();
89  switch(intMethod)
90  {
92  {
93  intSteps = 1;
94  }
95  break;
97  {
98  intSteps = 2;
99  }
100  break;
102  {
103  intSteps = 3;
104  }
105  break;
106  }
107  m_presForcingCorrection= Array<OneD, Array<OneD, NekDouble> >(intSteps);
108  for(int i = 0; i < m_presForcingCorrection.num_elements(); i++)
109  {
110  m_presForcingCorrection[i] = Array<OneD, NekDouble>(physTot,0.0);
111  }
112  m_verbose = (m_session->DefinesCmdLineArgument("verbose"))? true :false;
113 
114  // Load solve parameters related to the mapping
115  // Flags determining if pressure/viscous terms should be treated implicitly
116  m_session->MatchSolverInfo("MappingImplicitPressure","True",
117  m_implicitPressure,false);
118  m_session->MatchSolverInfo("MappingImplicitViscous","True",
119  m_implicitViscous,false);
120  m_session->MatchSolverInfo("MappingNeglectViscous","True",
121  m_neglectViscous,false);
122 
123  if (m_neglectViscous)
124  {
125  m_implicitViscous = false;
126  }
127 
128  // Tolerances and relaxation parameters for implicit terms
129  m_session->LoadParameter("MappingPressureTolerance",
130  m_pressureTolerance,1e-12);
131  m_session->LoadParameter("MappingViscousTolerance",
132  m_viscousTolerance,1e-12);
133  m_session->LoadParameter("MappingPressureRelaxation",
135  m_session->LoadParameter("MappingViscousRelaxation",
136  m_viscousRelaxation,1.0);
137 
138  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
NekDouble m_pressureTolerance
Definition: VCSMapping.h:87
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
ExtrapolateFactory & GetExtrapolateFactory()
Definition: Extrapolate.cpp:50
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
ExtrapolateSharedPtr m_extrapolation
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
NekDouble m_viscousRelaxation
Definition: VCSMapping.h:90
IMEX 2nd order scheme using Backward Different Formula & Extrapolation.
Array< OneD, Array< OneD, NekDouble > > m_presForcingCorrection
Definition: VCSMapping.h:122
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
IMEX 3rd order scheme using Backward Different Formula & Extrapolation.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
NekDouble m_viscousTolerance
Definition: VCSMapping.h:88
virtual void v_InitObject()
Init object for UnsteadySystem class.
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
IMEX 1st order scheme using Euler Backwards/Euler Forwards.
static GLOBAL_MAPPING_EXPORT MappingSharedPtr Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Return a pointer to the mapping, creating it on first call.
Definition: Mapping.cpp:265
NekDouble m_pressureRelaxation
Definition: VCSMapping.h:89
void Nektar::VCSMapping::v_SetUpPressureForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
Array< OneD, Array< OneD, NekDouble > > &  Forcing,
const NekDouble  aii_Dt 
)
protectedvirtual

Forcing term for Poisson solver solver

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 236 of file VCSMapping.cpp.

References Nektar::MultiRegions::DirCartesianMap, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, m_implicitViscous, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, m_presForcingCorrection, Vmath::Smul(), Nektar::VelocityCorrectionScheme::v_SetUpPressureForcing(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vdiv(), Vmath::Vmul(), Vmath::Vsub(), and Vmath::Zero().

240  {
241  if (m_mapping->HasConstantJacobian())
242  {
244  Forcing, aii_Dt);
245  }
246  else
247  {
248  int physTot = m_fields[0]->GetTotPoints();
249  int nvel = m_nConvectiveFields;
250  Array<OneD, NekDouble> wk(physTot, 0.0);
251 
252  Array<OneD, NekDouble> Jac(physTot,0.0);
253  m_mapping->GetJacobian(Jac);
254 
255  // Calculate div(J*u/Dt)
256  Vmath::Zero(physTot,Forcing[0],1);
257  for(int i = 0; i < nvel; ++i)
258  {
259  if (m_fields[i]->GetWaveSpace())
260  {
261  m_fields[i]->HomogeneousBwdTrans(fields[i],wk);
262  }
263  else
264  {
265  Vmath::Vcopy(physTot, fields[i], 1, wk, 1);
266  }
267  Vmath::Vmul(physTot,wk,1,Jac,1,wk,1);
268  if (m_fields[i]->GetWaveSpace())
269  {
270  m_fields[i]->HomogeneousFwdTrans(wk,wk);
271  }
272  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[i],wk, wk);
273  Vmath::Vadd(physTot,wk,1,Forcing[0],1,Forcing[0],1);
274  }
275  Vmath::Smul(physTot,1.0/aii_Dt,Forcing[0],1,Forcing[0],1);
276 
277  //
278  // If the mapping viscous terms are being treated explicitly
279  // we need to apply a correction to the forcing
280  if (!m_implicitViscous)
281  {
282  bool wavespace = m_fields[0]->GetWaveSpace();
283  m_fields[0]->SetWaveSpace(false);
284 
285  //
286  // Part 1: div(J*grad(U/J . grad(J)))
287  Array<OneD, Array<OneD, NekDouble> > tmp (nvel);
288  Array<OneD, Array<OneD, NekDouble> > velocity (nvel);
289  for(int i = 0; i < tmp.num_elements(); i++)
290  {
291  tmp[i] = Array<OneD, NekDouble>(physTot,0.0);
292  velocity[i] = Array<OneD, NekDouble>(physTot,0.0);
293  if (wavespace)
294  {
295  m_fields[0]->HomogeneousBwdTrans(m_fields[i]->GetPhys(),
296  velocity[i]);
297  }
298  else
299  {
300  Vmath::Vcopy(physTot, m_fields[i]->GetPhys(), 1,
301  velocity[i], 1);
302  }
303  }
304  // Calculate wk = U.grad(J)
305  m_mapping->DotGradJacobian(velocity, wk);
306  // Calculate wk = (U.grad(J))/J
307  Vmath::Vdiv(physTot, wk, 1, Jac, 1, wk, 1);
308  // J*grad[(U.grad(J))/J]
309  for(int i = 0; i < nvel; ++i)
310  {
311  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i],
312  wk, tmp[i]);
313  Vmath::Vmul(physTot, Jac, 1, tmp[i], 1, tmp[i], 1);
314  }
315  // div(J*grad[(U.grad(J))/J])
316  Vmath::Zero(physTot, wk, 1);
317  for(int i = 0; i < nvel; ++i)
318  {
319  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i],
320  tmp[i], tmp[i]);
321  Vmath::Vadd(physTot, wk, 1, tmp[i], 1, wk, 1);
322  }
323 
324  // Part 2: grad(J) . curl(curl(U))
325  m_mapping->CurlCurlField(velocity, tmp, m_implicitViscous);
326  // dont need velocity any more, so reuse it
327  m_mapping->DotGradJacobian(tmp, velocity[0]);
328 
329  // Add two parts
330  Vmath::Vadd(physTot, velocity[0], 1, wk, 1, wk, 1);
331 
332  // Multiply by kinvis and prepare to extrapolate
333  int nlevels = m_presForcingCorrection.num_elements();
334  Vmath::Smul(physTot, m_kinvis, wk, 1,
335  m_presForcingCorrection[nlevels-1], 1);
336 
337  // Extrapolate correction
338  m_extrapolation->ExtrapolateArray(m_presForcingCorrection);
339 
340  // Put in wavespace
341  if (wavespace)
342  {
343  m_fields[0]->HomogeneousFwdTrans(
344  m_presForcingCorrection[nlevels-1],wk);
345  }
346  else
347  {
348  Vmath::Vcopy(physTot, m_presForcingCorrection[nlevels-1], 1,
349  wk, 1);
350  }
351  // Apply correction: Forcing = Forcing - correction
352  Vmath::Vsub(physTot, Forcing[0], 1, wk, 1, Forcing[0], 1);
353 
354  m_fields[0]->SetWaveSpace(wavespace);
355  }
356  }
357  }
NekDouble m_kinvis
Kinematic viscosity.
ExtrapolateSharedPtr m_extrapolation
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.cpp:241
int m_nConvectiveFields
Number of fields to be convected;.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:213
Array< OneD, Array< OneD, NekDouble > > m_presForcingCorrection
Definition: VCSMapping.h:122
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:343
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
virtual void v_SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:299
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:183
void Nektar::VCSMapping::v_SetUpViscousForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  Forcing,
const NekDouble  aii_Dt 
)
protectedvirtual

Forcing term for Helmholtz solver

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 362 of file VCSMapping.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields, m_gradP, m_implicitPressure, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_pressure, Nektar::IncNavierStokes::m_velocity, Vmath::Vcopy(), and Vmath::Vdiv().

366  {
367  NekDouble aii_dtinv = 1.0/aii_Dt;
368  int physTot = m_fields[0]->GetTotPoints();
369 
370  // Grad p
371  m_pressure->BwdTrans(m_pressure->GetCoeffs(),m_pressure->UpdatePhys());
372 
373  int nvel = m_velocity.num_elements();
374  if(nvel == 2)
375  {
376  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1]);
377  }
378  else
379  {
380  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1],
381  Forcing[2]);
382  }
383 
384  // Copy grad p in physical space to m_gradP to reuse later
385  if (m_pressure->GetWaveSpace())
386  {
387  for (int i=0; i<nvel; i++)
388  {
389  m_pressure->HomogeneousBwdTrans(Forcing[i],m_gradP[i]);
390  }
391  }
392  else
393  {
394  for (int i=0; i<nvel; i++)
395  {
396  Vmath::Vcopy(physTot, Forcing[i], 1, m_gradP[i], 1);
397  }
398  }
399 
400  if ( (!m_mapping->HasConstantJacobian()) || m_implicitPressure)
401  {
402  // If pressure terms are treated explicitly, we need to divide by J
403  // if they are implicit, we need to calculate G(p)
404  if (m_implicitPressure)
405  {
406  m_mapping->RaiseIndex(m_gradP, Forcing);
407  }
408  else
409  {
410  Array<OneD, NekDouble> Jac(physTot,0.0);
411  m_mapping->GetJacobian(Jac);
412  for (int i=0; i<nvel; i++)
413  {
414  Vmath::Vdiv(physTot, m_gradP[i], 1, Jac, 1, Forcing[i], 1);
415  }
416  }
417  // Transform back to wavespace
418  if (m_pressure->GetWaveSpace())
419  {
420  for (int i=0; i<nvel; i++)
421  {
422  m_pressure->HomogeneousFwdTrans(Forcing[i],Forcing[i]);
423  }
424  }
425  }
426 
427  // Subtract inarray/(aii_dt) and divide by kinvis. Kinvis will
428  // need to be updated for the convected fields.
429  for(int i = 0; i < nvel; ++i)
430  {
431  Blas::Daxpy(physTot,-aii_dtinv,inarray[i],1,Forcing[i],1);
432  Blas::Dscal(physTot,1.0/m_kinvis,&(Forcing[i])[0],1);
433  }
434  }
NekDouble m_kinvis
Kinematic viscosity.
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.cpp:241
Array< OneD, Array< OneD, NekDouble > > m_gradP
Definition: VCSMapping.h:93
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Nektar::VCSMapping::v_SolvePressure ( const Array< OneD, NekDouble > &  Forcing)
protectedvirtual

Solve pressure system

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 439 of file VCSMapping.cpp.

References ASSERTL0, Nektar::MultiRegions::DirCartesianMap, Nektar::StdRegions::eFactorLambda, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, m_implicitPressure, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::m_pressure, m_pressureRelaxation, m_pressureTolerance, Nektar::SolverUtils::EquationSystem::m_session, m_verbose, Vmath::Neg(), Nektar::NullFlagList, Vmath::Smul(), Nektar::VelocityCorrectionScheme::v_SolvePressure(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vsub().

441  {
442  if (!m_implicitPressure)
443  {
445  }
446  else
447  {
448  int physTot = m_fields[0]->GetTotPoints();
449  int nvel = m_nConvectiveFields;
450  bool converged = false; // flag to mark if system converged
451  int s = 0; // iteration counter
452  NekDouble error; // L2 error at current iteration
453  NekDouble forcing_L2 = 0.0; // L2 norm of F
454 
455  int maxIter;
456  m_session->LoadParameter("MappingMaxIter",maxIter,5000);
457 
458  // rhs of the equation at current iteration
459  Array< OneD, NekDouble> F_corrected(physTot, 0.0);
460  // Pressure field at previous iteration
461  Array<OneD, NekDouble> previous_iter (physTot, 0.0);
462  // Temporary variables
463  Array<OneD, Array<OneD, NekDouble> > wk1(nvel);
464  Array<OneD, Array<OneD, NekDouble> > wk2(nvel);
465  Array<OneD, Array<OneD, NekDouble> > gradP(nvel);
466  for(int i = 0; i < nvel; ++i)
467  {
468  wk1[i] = Array<OneD, NekDouble> (physTot, 0.0);
469  wk2[i] = Array<OneD, NekDouble> (physTot, 0.0);
470  gradP[i] = Array<OneD, NekDouble> (physTot, 0.0);
471  }
472 
473  // Jacobian
474  Array<OneD, NekDouble> Jac(physTot, 0.0);
475  m_mapping->GetJacobian(Jac);
476 
477  // Factors for Laplacian system
479  factors[StdRegions::eFactorLambda] = 0.0;
480 
481  m_pressure->BwdTrans(m_pressure->GetCoeffs(),
482  m_pressure->UpdatePhys());
483  forcing_L2 = m_pressure->L2(Forcing, wk1[0]);
484  while (!converged)
485  {
486  // Update iteration counter and set previous iteration field
487  // (use previous timestep solution for first iteration)
488  s++;
489  ASSERTL0(s < maxIter,
490  "VCSMapping exceeded maximum number of iterations.");
491 
492  Vmath::Vcopy(physTot, m_pressure->GetPhys(), 1,
493  previous_iter, 1);
494 
495  // Correct pressure bc to account for iteration
496  m_extrapolation->CorrectPressureBCs(previous_iter);
497 
498  //
499  // Calculate forcing term for this iteration
500  //
501  for(int i = 0; i < nvel; ++i)
502  {
504  previous_iter, gradP[i]);
505  if(m_pressure->GetWaveSpace())
506  {
507  m_pressure->HomogeneousBwdTrans(gradP[i], wk1[i]);
508  }
509  else
510  {
511  Vmath::Vcopy(physTot, gradP[i], 1, wk1[i], 1);
512  }
513  }
514  m_mapping->RaiseIndex(wk1, wk2); // G(p)
515 
516  m_mapping->Divergence(wk2, F_corrected); // div(G(p))
517  if (!m_mapping->HasConstantJacobian())
518  {
519  Vmath::Vmul(physTot, F_corrected, 1,
520  Jac, 1,
521  F_corrected, 1);
522  }
523  // alpha*J*div(G(p))
524  Vmath::Smul(physTot, m_pressureRelaxation, F_corrected, 1,
525  F_corrected, 1);
526  if(m_pressure->GetWaveSpace())
527  {
528  m_pressure->HomogeneousFwdTrans(F_corrected, F_corrected);
529  }
530  // alpha*J*div(G(p)) - p_ii
531  for (int i = 0; i < m_nConvectiveFields; ++i)
532  {
534  gradP[i], wk1[0]);
535  Vmath::Vsub(physTot, F_corrected, 1, wk1[0], 1,
536  F_corrected, 1);
537  }
538  // p_i,i - J*div(G(p))
539  Vmath::Neg(physTot, F_corrected, 1);
540  // alpha*F - alpha*J*div(G(p)) + p_i,i
541  Vmath::Smul(physTot, m_pressureRelaxation, Forcing, 1,
542  wk1[0], 1);
543  Vmath::Vadd(physTot, wk1[0], 1, F_corrected, 1, F_corrected, 1);
544 
545  //
546  // Solve system
547  //
548  m_pressure->HelmSolve(F_corrected, m_pressure->UpdateCoeffs(),
549  NullFlagList, factors);
550  m_pressure->BwdTrans(m_pressure->GetCoeffs(),
551  m_pressure->UpdatePhys());
552 
553  //
554  // Test convergence
555  //
556  error = m_pressure->L2(m_pressure->GetPhys(), previous_iter);
557  if ( forcing_L2 != 0)
558  {
559  if ( (error/forcing_L2 < m_pressureTolerance))
560  {
561  converged = true;
562  }
563  }
564  else
565  {
566  if ( error < m_pressureTolerance)
567  {
568  converged = true;
569  }
570  }
571  }
572  if (m_verbose && m_session->GetComm()->GetRank()==0)
573  {
574  std::cout << " Pressure system (mapping) converged in " << s <<
575  " iterations with error = " << error << std::endl;
576  }
577  }
578  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
NekDouble m_pressureTolerance
Definition: VCSMapping.h:87
ExtrapolateSharedPtr m_extrapolation
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
int m_nConvectiveFields
Number of fields to be convected;.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:213
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:396
double NekDouble
virtual void v_SolvePressure(const Array< OneD, NekDouble > &Forcing)
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:343
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:299
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:183
NekDouble m_pressureRelaxation
Definition: VCSMapping.h:89
static FlagList NullFlagList
An empty flag list.
void Nektar::VCSMapping::v_SolveViscous ( const Array< OneD, const Array< OneD, NekDouble > > &  Forcing,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  aii_Dt 
)
protectedvirtual

Solve velocity system

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 583 of file VCSMapping.cpp.

References ASSERTL0, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::SolverUtils::EquationSystem::m_fields, m_implicitViscous, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::VelocityCorrectionScheme::m_sVVCutoffRatio, Nektar::VelocityCorrectionScheme::m_sVVDiffCoeff, Nektar::VelocityCorrectionScheme::m_useSpecVanVisc, m_verbose, m_viscousRelaxation, m_viscousTolerance, Nektar::NullFlagList, Vmath::Smul(), Nektar::VelocityCorrectionScheme::v_SolveViscous(), Vmath::Vadd(), and Vmath::Vcopy().

587  {
588  if(!m_implicitViscous)
589  {
590  VelocityCorrectionScheme::v_SolveViscous(Forcing, outarray, aii_Dt);
591  }
592  else
593  {
594  int physTot = m_fields[0]->GetTotPoints();
595  int nvel = m_nConvectiveFields;
596  bool converged = false; // flag to mark if system converged
597  int s = 0; // iteration counter
598  NekDouble error, max_error; // L2 error at current iteration
599 
600  int maxIter;
601  m_session->LoadParameter("MappingMaxIter",maxIter,5000);
602 
603  //L2 norm of F
604  Array<OneD, NekDouble> forcing_L2(m_nConvectiveFields,0.0);
605 
606  // rhs of the equation at current iteration
607  Array<OneD, Array<OneD, NekDouble> > F_corrected(nvel);
608  // Solution at previous iteration
609  Array<OneD, Array<OneD, NekDouble> > previous_iter(nvel);
610  // Working space
611  Array<OneD, Array<OneD, NekDouble> > wk(nvel);
612  for(int i = 0; i < nvel; ++i)
613  {
614  F_corrected[i] = Array<OneD, NekDouble> (physTot, 0.0);
615  previous_iter[i] = Array<OneD, NekDouble> (physTot, 0.0);
616  wk[i] = Array<OneD, NekDouble> (physTot, 0.0);
617  }
618 
619  // Factors for Helmholtz system
621  factors[StdRegions::eFactorLambda] =
622  1.0*m_viscousRelaxation/aii_Dt/m_kinvis;
623  if(m_useSpecVanVisc)
624  {
628  }
629 
630  // Calculate L2-norm of F and set initial solution for iteration
631  for(int i = 0; i < nvel; ++i)
632  {
633  forcing_L2[i] = m_fields[0]->L2(Forcing[i],wk[0]);
634  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
635  previous_iter[i]);
636  }
637 
638  while (!converged)
639  {
640  converged = true;
641  // Iteration counter
642  s++;
643  ASSERTL0(s < maxIter,
644  "VCSMapping exceeded maximum number of iterations.");
645 
646  max_error = 0.0;
647 
648  //
649  // Calculate forcing term for next iteration
650  //
651 
652  // Calculate L(U)- in this parts all components might be coupled
653  if(m_fields[0]->GetWaveSpace())
654  {
655  for (int i = 0; i < nvel; ++i)
656  {
657  m_fields[0]->HomogeneousBwdTrans(previous_iter[i],
658  wk[i]);
659  }
660  }
661  else
662  {
663  for (int i = 0; i < nvel; ++i)
664  {
665  Vmath::Vcopy(physTot, previous_iter[i], 1, wk[i], 1);
666  }
667  }
668 
669  // (L(U^i) - 1/alpha*U^i_jj)
670  m_mapping->VelocityLaplacian(wk, F_corrected,
671  1.0/m_viscousRelaxation);
672 
673  if(m_fields[0]->GetWaveSpace())
674  {
675  for (int i = 0; i < nvel; ++i)
676  {
677  m_fields[0]->HomogeneousFwdTrans(F_corrected[i],
678  F_corrected[i]);
679  }
680  }
681  else
682  {
683  for (int i = 0; i < nvel; ++i)
684  {
685  Vmath::Vcopy(physTot, F_corrected[i], 1,
686  F_corrected[i], 1);
687  }
688  }
689 
690  // Loop velocity components
691  for (int i = 0; i < nvel; ++i)
692  {
693  // (-alpha*L(U^i) + U^i_jj)
694  Vmath::Smul(physTot, -1.0*m_viscousRelaxation,
695  F_corrected[i], 1,
696  F_corrected[i], 1);
697  // F_corrected = alpha*F + (-alpha*L(U^i) + U^i_jj)
698  Vmath::Smul(physTot, m_viscousRelaxation, Forcing[i], 1,
699  wk[0], 1);
700  Vmath::Vadd(physTot, wk[0], 1, F_corrected[i], 1,
701  F_corrected[i], 1);
702 
703  //
704  // Solve System
705  //
706  m_fields[i]->HelmSolve(F_corrected[i],
707  m_fields[i]->UpdateCoeffs(),
708  NullFlagList, factors);
709  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),outarray[i]);
710 
711  //
712  // Test convergence
713  //
714  error = m_fields[i]->L2(outarray[i], previous_iter[i]);
715 
716  if ( forcing_L2[i] != 0)
717  {
718  if ( (error/forcing_L2[i] >= m_viscousTolerance))
719  {
720  converged = false;
721  }
722  }
723  else
724  {
725  if ( error >= m_viscousTolerance)
726  {
727  converged = false;
728  }
729  }
730  if (error > max_error)
731  {
732  max_error = error;
733  }
734 
735  // Copy field to previous_iter
736  Vmath::Vcopy(physTot, outarray[i], 1, previous_iter[i], 1);
737  }
738  }
739  if (m_verbose && m_session->GetComm()->GetRank()==0)
740  {
741  std::cout << " Velocity system (mapping) converged in " << s <<
742  " iterations with error = " << max_error << std::endl;
743  }
744  }
745  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
NekDouble m_kinvis
Kinematic viscosity.
NekDouble m_sVVDiffCoeff
Diffusion coefficient of SVV modes.
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:76
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
int m_nConvectiveFields
Number of fields to be convected;.
NekDouble m_viscousRelaxation
Definition: VCSMapping.h:90
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:213
double NekDouble
virtual void v_SolveViscous(const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
NekDouble m_sVVCutoffRatio
cutt off ratio from which to start decayhing modes
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
bool m_useSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
NekDouble m_viscousTolerance
Definition: VCSMapping.h:88
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:299
static FlagList NullFlagList
An empty flag list.

Member Data Documentation

string Nektar::VCSMapping::className
static
Initial value:

Name of class.

Definition at line 59 of file VCSMapping.h.

Array<OneD, Array<OneD, NekDouble> > Nektar::VCSMapping::m_gradP
protected

Definition at line 93 of file VCSMapping.h.

Referenced by MappingPressureCorrection(), v_DoInitialise(), and v_SetUpViscousForcing().

bool Nektar::VCSMapping::m_implicitPressure
protected
bool Nektar::VCSMapping::m_implicitViscous
protected
GlobalMapping::MappingSharedPtr Nektar::VCSMapping::m_mapping
protected
bool Nektar::VCSMapping::m_neglectViscous
protected

Definition at line 84 of file VCSMapping.h.

Referenced by ApplyIncNSMappingForcing(), and v_InitObject().

Array<OneD, Array<OneD, NekDouble> > Nektar::VCSMapping::m_presForcingCorrection
private

Definition at line 122 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SetUpPressureForcing().

NekDouble Nektar::VCSMapping::m_pressureRelaxation
protected

Definition at line 89 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

NekDouble Nektar::VCSMapping::m_pressureTolerance
protected

Definition at line 87 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

bool Nektar::VCSMapping::m_verbose
protected

Definition at line 78 of file VCSMapping.h.

Referenced by v_InitObject(), v_SolvePressure(), and v_SolveViscous().

NekDouble Nektar::VCSMapping::m_viscousRelaxation
protected

Definition at line 90 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().

NekDouble Nektar::VCSMapping::m_viscousTolerance
protected

Definition at line 88 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().