Nektar++
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:
[legend]

Public Member Functions

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

Static Public Member Functions

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

Static Public Attributes

static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::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
void SetupFlowrate (NekDouble aii_dt)
 Set up the Stokes solution used to impose constant flowrate through a boundary. More...
 
NekDouble MeasureFlowrate (const Array< OneD, Array< OneD, NekDouble > > &inarray)
 Measure the volumetric flow rate through the volumetric flow rate reference surface. More...
 
virtual bool v_PostIntegrate (int step)
 
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)
 
void SetUpSVV (void)
 
void SetUpExtrapolation (void)
 
void SVVVarDiffCoeff (const NekDouble velmag, Array< OneD, NekDouble > &diffcoeff, const Array< OneD, Array< OneD, NekDouble > > &vel=NullNekDoubleArrayOfArray)
 
void AppendSVVFactors (StdRegions::ConstFactorMap &factors, MultiRegions::VarFactorsMap &varFactorsMap)
 
- Protected Member Functions inherited from Nektar::IncNavierStokes
 IncNavierStokes (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
EquationType GetEquationType (void)
 
void EvaluateAdvectionTerms (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
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 fldid, const int bndid, std::string womstr)
 Set Up Womersley details. More...
 
virtual MultiRegions::ExpListSharedPtr v_GetPressure ()
 
virtual Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor)
 
virtual bool v_PreIntegrate (int step)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 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 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 void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
virtual SOLVER_UTILS_EXPORT bool UpdateTimeStepCheck ()
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_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 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...
 
NekDouble m_sVVCutoffRatioHomo1D
 
NekDouble m_sVVDiffCoeffHomo1D
 Diffusion coefficient of SVV modes in homogeneous 1D Direction. More...
 
Array< OneD, NekDoublem_svvVarDiffCoeff
 Array of coefficient if power kernel is used in SVV. More...
 
bool m_IsSVVPowerKernel
 Identifier for Power Kernel otherwise DG kernel. 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...
 
NekDouble m_flowrate
 Desired volumetric flowrate. More...
 
NekDouble m_flowrateArea
 Area of the boundary through which we are measuring the flowrate. More...
 
bool m_homd1DFlowinPlane
 
NekDouble m_greenFlux
 Flux of the Stokes function solution. More...
 
NekDouble m_alpha
 Current flowrate correction. More...
 
int m_flowrateBndID
 Boundary ID of the flowrate reference surface. More...
 
int m_planeID
 Plane ID for cases with homogeneous expansion. More...
 
MultiRegions::ExpListSharedPtr m_flowrateBnd
 Flowrate reference surface. More...
 
Array< OneD, Array< OneD, NekDouble > > m_flowrateStokes
 Stokes solution used to impose flowrate. More...
 
std::ofstream m_flowrateStream
 Output stream to record flowrate. More...
 
int m_flowrateSteps
 Interval at which to record flowrate data. More...
 
NekDouble m_flowrateAiidt
 Value of aii_dt used to compute Stokes flowrate solution. 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::ForcingSharedPtrm_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...
 
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, 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_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
NekDouble m_TimeIntegLambda =0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm =-1.0
 
int m_TotLinItePerStep =0
 
int m_StagesPerStep =1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts =0
 
int m_TotLinIts =0
 
int m_TotImpStages =0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
bool m_root
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_timestepMax = -1.0
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_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...
 
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...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth More...
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

Detailed Description

Definition at line 43 of file VCSMapping.h.

Constructor & Destructor Documentation

◆ VCSMapping()

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

Constructor.

Constructor. Creates ...

Parameters

Definition at line 56 of file VCSMapping.cpp.

59  : UnsteadySystem(pSession, pGraph),
60  VelocityCorrectionScheme(pSession, pGraph)
61  {
62 
63  }
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
VelocityCorrectionScheme(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Constructor.

◆ ~VCSMapping()

Nektar::VCSMapping::~VCSMapping ( void  )
virtual

Destructor

Definition at line 135 of file VCSMapping.cpp.

136  {
137  }

Member Function Documentation

◆ ApplyIncNSMappingForcing()

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 738 of file VCSMapping.cpp.

741  {
742  int physTot = m_fields[0]->GetTotPoints();
743  Array<OneD, Array<OneD, NekDouble> > vel(m_nConvectiveFields);
744  Array<OneD, Array<OneD, NekDouble> > velPhys(m_nConvectiveFields);
745  Array<OneD, Array<OneD, NekDouble> > Forcing(m_nConvectiveFields);
746  Array<OneD, Array<OneD, NekDouble> > tmp(m_nConvectiveFields);
747  for (int i = 0; i < m_nConvectiveFields; ++i)
748  {
749  velPhys[i] = Array<OneD, NekDouble> (physTot, 0.0);
750  Forcing[i] = Array<OneD, NekDouble> (physTot, 0.0);
751  tmp[i] = Array<OneD, NekDouble> (physTot, 0.0);
752  }
753 
754  // Get fields and store velocity in wavespace and physical space
755  if(m_fields[0]->GetWaveSpace())
756  {
757  for (int i = 0; i < m_nConvectiveFields; ++i)
758  {
759  vel[i] = inarray[i];
760  m_fields[0]->HomogeneousBwdTrans(vel[i],velPhys[i]);
761  }
762  }
763  else
764  {
765  for (int i = 0; i < m_nConvectiveFields; ++i)
766  {
767  vel[i] = inarray[i];
768  Vmath::Vcopy(physTot, inarray[i], 1, velPhys[i], 1);
769  }
770  }
771 
772  //Advection contribution
773  MappingAdvectionCorrection(velPhys, Forcing);
774 
775  // Time-derivative contribution
776  if ( m_mapping->IsTimeDependent() )
777  {
778  MappingAccelerationCorrection(vel, velPhys, tmp);
779  for (int i = 0; i < m_nConvectiveFields; ++i)
780  {
781  Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
782  }
783  }
784 
785  // Pressure contribution
786  if (!m_implicitPressure)
787  {
789  for (int i = 0; i < m_nConvectiveFields; ++i)
790  {
791  Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
792  }
793  }
794  // Viscous contribution
795  if ( (!m_implicitViscous) && (!m_neglectViscous))
796  {
797  MappingViscousCorrection(velPhys, tmp);
798  for (int i = 0; i < m_nConvectiveFields; ++i)
799  {
800  Vmath::Smul(physTot, m_kinvis, tmp[i], 1, tmp[i], 1);
801  Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
802  }
803  }
804 
805  // If necessary, transform to wavespace
806  if(m_fields[0]->GetWaveSpace())
807  {
808  for (int i = 0; i < m_nConvectiveFields; ++i)
809  {
810  m_fields[0]->HomogeneousFwdTrans(Forcing[i],Forcing[i]);
811  }
812  }
813 
814  // Add to outarray
815  for (int i = 0; i < m_nConvectiveFields; ++i)
816  {
817  Vmath::Vadd(physTot, outarray[i], 1, Forcing[i], 1, outarray[i], 1);
818  }
819  }
NekDouble m_kinvis
Kinematic viscosity.
int m_nConvectiveFields
Number of fields to be convected;.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void MappingViscousCorrection(const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:952
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:77
void MappingPressureCorrection(Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:925
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:846
void MappingAdvectionCorrection(const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:821
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:322
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199

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

◆ create()

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

Creates an instance of this class.

Definition at line 48 of file VCSMapping.h.

51  {
54  p->InitObject();
55  return p;
56  }
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.

◆ MappingAccelerationCorrection()

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 846 of file VCSMapping.cpp.

850  {
851  int physTot = m_fields[0]->GetTotPoints();
852  int nvel = m_nConvectiveFields;
853 
854  Array<OneD, Array<OneD, NekDouble> > wk(nvel*nvel);
855  Array<OneD, Array<OneD, NekDouble> > tmp(nvel);
856  Array<OneD, Array<OneD, NekDouble> > coordVel(nvel);
857  for (int i = 0; i< nvel; i++)
858  {
859  tmp[i] = Array<OneD, NekDouble> (physTot, 0.0);
860  coordVel[i] = Array<OneD, NekDouble> (physTot, 0.0);
861  }
862  // Get coordinates velocity in transformed system
863  m_mapping->GetCoordVelocity(tmp);
864  m_mapping->ContravarFromCartesian(tmp, coordVel);
865 
866  // Calculate first term: U^j u^i,j = U^j (du^i/dx^j + {i,kj}u^k)
867  m_mapping->ApplyChristoffelContravar(velPhys, wk);
868  for (int i=0; i< nvel; i++)
869  {
870  Vmath::Zero(physTot,outarray[i],1);
871 
872  m_fields[0]->PhysDeriv(velPhys[i], tmp[0], tmp[1]);
873  for (int j=0; j< nvel; j++)
874  {
875  if (j == 2)
876  {
877  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j],
878  vel[i], tmp[2]);
879  if (m_fields[0]->GetWaveSpace())
880  {
881  m_fields[0]->HomogeneousBwdTrans(tmp[2],tmp[2]);
882  }
883  }
884 
885  Vmath::Vadd(physTot,wk[i*nvel+j],1,tmp[j],1,
886  wk[i*nvel+j], 1);
887 
888  Vmath::Vvtvp(physTot, coordVel[j], 1, wk[i*nvel+j], 1,
889  outarray[i], 1, outarray[i], 1);
890  }
891  }
892 
893  // Set wavespace to false and store current value
894  bool wavespace = m_fields[0]->GetWaveSpace();
895  m_fields[0]->SetWaveSpace(false);
896 
897  // Add -u^j U^i,j
898  m_mapping->ApplyChristoffelContravar(coordVel, wk);
899  for (int i=0; i< nvel; i++)
900  {
901  if(nvel == 2)
902  {
903  m_fields[0]->PhysDeriv(coordVel[i], tmp[0], tmp[1]);
904  }
905  else
906  {
907  m_fields[0]->PhysDeriv(coordVel[i], tmp[0], tmp[1], tmp[2]);
908  }
909 
910  for (int j=0; j< nvel; j++)
911  {
912  Vmath::Vadd(physTot,wk[i*nvel+j],1,tmp[j],1,
913  wk[i*nvel+j], 1);
914  Vmath::Neg(physTot, wk[i*nvel+j], 1);
915 
916  Vmath::Vvtvp(physTot, velPhys[j], 1, wk[i*nvel+j], 1,
917  outarray[i], 1, outarray[i], 1);
918  }
919  }
920 
921  // Restore value of wavespace
922  m_fields[0]->SetWaveSpace(wavespace);
923  }
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:90
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461
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:513
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436

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

◆ MappingAdvectionCorrection()

void Nektar::VCSMapping::MappingAdvectionCorrection ( const Array< OneD, const Array< OneD, NekDouble > > &  velPhys,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 821 of file VCSMapping.cpp.

824  {
825  int physTot = m_fields[0]->GetTotPoints();
826  int nvel = m_nConvectiveFields;
827 
828  Array<OneD, Array<OneD, NekDouble> > wk(nvel*nvel);
829 
830  // Apply Christoffel symbols to obtain {i,kj}vel(k)
831  m_mapping->ApplyChristoffelContravar(velPhys, wk);
832 
833  // Calculate correction -U^j*{i,kj}vel(k)
834  for (int i = 0; i< nvel; i++)
835  {
836  Vmath::Zero(physTot,outarray[i],1);
837  for (int j = 0; j< nvel; j++)
838  {
839  Vmath::Vvtvp(physTot,wk[i*nvel+j],1,velPhys[j],1,
840  outarray[i],1,outarray[i],1);
841  }
842  Vmath::Neg(physTot, outarray[i], 1);
843  }
844  }

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

Referenced by ApplyIncNSMappingForcing().

◆ MappingPressureCorrection()

void Nektar::VCSMapping::MappingPressureCorrection ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
private

Definition at line 925 of file VCSMapping.cpp.

927  {
928  int physTot = m_fields[0]->GetTotPoints();
929  int nvel = m_nConvectiveFields;
930 
931  // Calculate g^(ij)p_(,j)
932  m_mapping->RaiseIndex(m_gradP, outarray);
933 
934  // Calculate correction = (nabla p)/J - g^(ij)p_,j
935  // (Jac is not required if it is constant)
936  if ( !m_mapping->HasConstantJacobian())
937  {
938  Array<OneD, NekDouble> Jac(physTot, 0.0);
939  m_mapping->GetJacobian(Jac);
940  for(int i = 0; i < nvel; ++i)
941  {
942  Vmath::Vdiv(physTot, m_gradP[i], 1, Jac, 1, m_gradP[i], 1);
943  }
944  }
945  for(int i = 0; i < nvel; ++i)
946  {
947  Vmath::Vsub(physTot, m_gradP[i], 1,outarray[i], 1,
948  outarray[i],1);
949  }
950  }
Array< OneD, Array< OneD, NekDouble > > m_gradP
Definition: VCSMapping.h:94
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:257
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:372

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

Referenced by ApplyIncNSMappingForcing().

◆ MappingViscousCorrection()

void Nektar::VCSMapping::MappingViscousCorrection ( const Array< OneD, const Array< OneD, NekDouble > > &  velPhys,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 952 of file VCSMapping.cpp.

955  {
956  // L(U) - 1.0*d^2(u^i)/dx^jdx^j
957  m_mapping->VelocityLaplacian(velPhys, outarray, 1.0);
958  }

References m_mapping.

Referenced by ApplyIncNSMappingForcing().

◆ v_DoInitialise()

void Nektar::VCSMapping::v_DoInitialise ( void  )
protectedvirtual

Sets up initial conditions.

Sets the initial conditions.

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 139 of file VCSMapping.cpp.

140  {
142 
143  // Set up Field Meta Data for output files
144  m_fieldMetaDataMap["Kinvis"] =
145  boost::lexical_cast<std::string>(m_kinvis);
146  m_fieldMetaDataMap["TimeStep"] =
147  boost::lexical_cast<std::string>(m_timestep);
148 
149  // Correct Dirichlet boundary conditions to account for mapping
150  m_mapping->UpdateBCs(0.0);
151  //
152  m_F = Array<OneD, Array< OneD, NekDouble> > (m_nConvectiveFields);
153  for(int i = 0; i < m_nConvectiveFields; ++i)
154  {
155  m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
156  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
157  m_fields[i]->UpdatePhys());
158  m_F[i] = Array< OneD, NekDouble> (m_fields[0]->GetTotPoints(), 0.0);
159  }
160 
161  // Initialise m_gradP
162  int physTot = m_fields[0]->GetTotPoints();
163  m_gradP = Array<OneD, Array<OneD, NekDouble> >(m_nConvectiveFields);
164  for(int i = 0; i < m_nConvectiveFields; ++i)
165  {
166  m_gradP[i] = Array<OneD, NekDouble>(physTot,0.0);
168  m_pressure->GetPhys(),
169  m_gradP[i]);
170  if(m_pressure->GetWaveSpace())
171  {
172  m_pressure->HomogeneousBwdTrans(m_gradP[i],m_gradP[i]);
173  }
174  }
175  }
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
NekDouble m_timestep
Time step size.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
SOLVER_UTILS_EXPORT int GetTotPoints()
virtual SOLVER_UTILS_EXPORT void v_DoInitialise()
Sets up initial conditions.
Array< OneD, Array< OneD, NekDouble > > m_F

References Nektar::MultiRegions::DirCartesianMap, Nektar::SolverUtils::EquationSystem::GetTotPoints(), 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().

◆ v_EvaluateAdvection_SetPressureBCs()

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 180 of file VCSMapping.cpp.

184  {
185  EvaluateAdvectionTerms(inarray, outarray, time);
186 
187  // Smooth advection
189  {
190  for(int i = 0; i < m_nConvectiveFields; ++i)
191  {
192  m_pressure->SmoothField(outarray[i]);
193  }
194  }
195 
196  // Add forcing terms
197  for (auto &x : m_forcing)
198  {
199  x->Apply(m_fields, inarray, outarray, time);
200  }
201 
202  // Add mapping terms
203  ApplyIncNSMappingForcing( inarray, outarray);
204 
205  // Calculate High-Order pressure boundary conditions
206  m_extrapolation->EvaluatePressureBCs(inarray,outarray,m_kinvis);
207 
208  // Update mapping and deal with Dirichlet boundary conditions
209  if (m_mapping->IsTimeDependent())
210  {
211  if (m_mapping->IsFromFunction())
212  {
213  // If the transformation is explicitly defined, update it here
214  // Otherwise, it will be done somewhere else (ForcingMovingBody)
215  m_mapping->UpdateMapping(time+m_timestep);
216  }
217  m_mapping->UpdateBCs(time+m_timestep);
218  }
219  }
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested
ExtrapolateSharedPtr m_extrapolation
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
void EvaluateAdvectionTerms(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void ApplyIncNSMappingForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:738

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.

◆ v_InitObject()

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.

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(m_intScheme);
82  m_extrapolation->GenerateHOPBCMap(m_session);
83 
84  // Storage to extrapolate pressure forcing
85  int physTot = m_fields[0]->GetTotPoints();
86  int intSteps = 1;
87 
88  if ( m_intScheme->GetName() == "IMEX" ||
89  m_intScheme->GetName() == "IMEXGear" )
90  {
91  m_intSteps = m_intScheme->GetOrder();
92  }
93  else
94  {
95  NEKERROR(ErrorUtil::efatal, "Integration method not suitable: "
96  "Options include IMEXGear or IMEXOrder{1,2,3,4}");
97  }
98 
99  m_presForcingCorrection= Array<OneD, Array<OneD, NekDouble> >(intSteps);
100  for(int i = 0; i < m_presForcingCorrection.size(); i++)
101  {
102  m_presForcingCorrection[i] = Array<OneD, NekDouble>(physTot,0.0);
103  }
104  m_verbose = (m_session->DefinesCmdLineArgument("verbose"))? true :false;
105 
106  // Load solve parameters related to the mapping
107  // Flags determining if pressure/viscous terms should be treated implicitly
108  m_session->MatchSolverInfo("MappingImplicitPressure","True",
109  m_implicitPressure,false);
110  m_session->MatchSolverInfo("MappingImplicitViscous","True",
111  m_implicitViscous,false);
112  m_session->MatchSolverInfo("MappingNeglectViscous","True",
113  m_neglectViscous,false);
114 
115  if (m_neglectViscous)
116  {
117  m_implicitViscous = false;
118  }
119 
120  // Tolerances and relaxation parameters for implicit terms
121  m_session->LoadParameter("MappingPressureTolerance",
122  m_pressureTolerance,1e-12);
123  m_session->LoadParameter("MappingViscousTolerance",
124  m_viscousTolerance,1e-12);
125  m_session->LoadParameter("MappingPressureRelaxation",
127  m_session->LoadParameter("MappingViscousRelaxation",
128  m_viscousRelaxation,1.0);
129 
130  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
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:268
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
int m_intSteps
Number of time integration steps AND Order of extrapolation for pressure boundary conditions.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.
NekDouble m_pressureTolerance
Definition: VCSMapping.h:88
NekDouble m_viscousRelaxation
Definition: VCSMapping.h:91
NekDouble m_pressureRelaxation
Definition: VCSMapping.h:90
NekDouble m_viscousTolerance
Definition: VCSMapping.h:89
Array< OneD, Array< OneD, NekDouble > > m_presForcingCorrection
Definition: VCSMapping.h:123
virtual void v_InitObject()
Init object for UnsteadySystem class.
ExtrapolateFactory & GetExtrapolateFactory()
Definition: Extrapolate.cpp:49

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::ErrorUtil::efatal, 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, Nektar::IncNavierStokes::m_intSteps, 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, NEKERROR, and Nektar::VelocityCorrectionScheme::v_InitObject().

◆ v_SetUpPressureForcing()

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 225 of file VCSMapping.cpp.

229  {
230  if (m_mapping->HasConstantJacobian())
231  {
233  Forcing, aii_Dt);
234  }
235  else
236  {
237  int physTot = m_fields[0]->GetTotPoints();
238  int nvel = m_nConvectiveFields;
239  Array<OneD, NekDouble> wk(physTot, 0.0);
240 
241  Array<OneD, NekDouble> Jac(physTot,0.0);
242  m_mapping->GetJacobian(Jac);
243 
244  // Calculate div(J*u/Dt)
245  Vmath::Zero(physTot,Forcing[0],1);
246  for(int i = 0; i < nvel; ++i)
247  {
248  if (m_fields[i]->GetWaveSpace())
249  {
250  m_fields[i]->HomogeneousBwdTrans(fields[i],wk);
251  }
252  else
253  {
254  Vmath::Vcopy(physTot, fields[i], 1, wk, 1);
255  }
256  Vmath::Vmul(physTot,wk,1,Jac,1,wk,1);
257  if (m_fields[i]->GetWaveSpace())
258  {
259  m_fields[i]->HomogeneousFwdTrans(wk,wk);
260  }
261  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[i],wk, wk);
262  Vmath::Vadd(physTot,wk,1,Forcing[0],1,Forcing[0],1);
263  }
264  Vmath::Smul(physTot,1.0/aii_Dt,Forcing[0],1,Forcing[0],1);
265 
266  //
267  // If the mapping viscous terms are being treated explicitly
268  // we need to apply a correction to the forcing
269  if (!m_implicitViscous)
270  {
271  bool wavespace = m_fields[0]->GetWaveSpace();
272  m_fields[0]->SetWaveSpace(false);
273 
274  //
275  // Part 1: div(J*grad(U/J . grad(J)))
276  Array<OneD, Array<OneD, NekDouble> > tmp (nvel);
277  Array<OneD, Array<OneD, NekDouble> > velocity (nvel);
278  for(int i = 0; i < tmp.size(); i++)
279  {
280  tmp[i] = Array<OneD, NekDouble>(physTot,0.0);
281  velocity[i] = Array<OneD, NekDouble>(physTot,0.0);
282  if (wavespace)
283  {
284  m_fields[0]->HomogeneousBwdTrans(m_fields[i]->GetPhys(),
285  velocity[i]);
286  }
287  else
288  {
289  Vmath::Vcopy(physTot, m_fields[i]->GetPhys(), 1,
290  velocity[i], 1);
291  }
292  }
293  // Calculate wk = U.grad(J)
294  m_mapping->DotGradJacobian(velocity, wk);
295  // Calculate wk = (U.grad(J))/J
296  Vmath::Vdiv(physTot, wk, 1, Jac, 1, wk, 1);
297  // J*grad[(U.grad(J))/J]
298  for(int i = 0; i < nvel; ++i)
299  {
300  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i],
301  wk, tmp[i]);
302  Vmath::Vmul(physTot, Jac, 1, tmp[i], 1, tmp[i], 1);
303  }
304  // div(J*grad[(U.grad(J))/J])
305  Vmath::Zero(physTot, wk, 1);
306  for(int i = 0; i < nvel; ++i)
307  {
308  m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i],
309  tmp[i], tmp[i]);
310  Vmath::Vadd(physTot, wk, 1, tmp[i], 1, wk, 1);
311  }
312 
313  // Part 2: grad(J) . curl(curl(U))
314  m_mapping->CurlCurlField(velocity, tmp, m_implicitViscous);
315  // dont need velocity any more, so reuse it
316  m_mapping->DotGradJacobian(tmp, velocity[0]);
317 
318  // Add two parts
319  Vmath::Vadd(physTot, velocity[0], 1, wk, 1, wk, 1);
320 
321  // Multiply by kinvis and prepare to extrapolate
322  int nlevels = m_presForcingCorrection.size();
323  Vmath::Smul(physTot, m_kinvis, wk, 1,
324  m_presForcingCorrection[nlevels-1], 1);
325 
326  // Extrapolate correction
327  m_extrapolation->ExtrapolateArray(m_presForcingCorrection);
328 
329  // Put in wavespace
330  if (wavespace)
331  {
332  m_fields[0]->HomogeneousFwdTrans(
333  m_presForcingCorrection[nlevels-1],wk);
334  }
335  else
336  {
337  Vmath::Vcopy(physTot, m_presForcingCorrection[nlevels-1], 1,
338  wk, 1);
339  }
340  // Apply correction: Forcing = Forcing - correction
341  Vmath::Vsub(physTot, Forcing[0], 1, wk, 1, Forcing[0], 1);
342 
343  m_fields[0]->SetWaveSpace(wavespace);
344  }
345  }
346  }
virtual void v_SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
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:192

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

◆ v_SetUpViscousForcing()

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 351 of file VCSMapping.cpp.

355  {
356  NekDouble aii_dtinv = 1.0/aii_Dt;
357  int physTot = m_fields[0]->GetTotPoints();
358 
359  // Grad p
360  m_pressure->BwdTrans(m_pressure->GetCoeffs(),m_pressure->UpdatePhys());
361 
362  int nvel = m_velocity.size();
363  if(nvel == 2)
364  {
365  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1]);
366  }
367  else
368  {
369  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1],
370  Forcing[2]);
371  }
372 
373  // Copy grad p in physical space to m_gradP to reuse later
374  if (m_pressure->GetWaveSpace())
375  {
376  for (int i=0; i<nvel; i++)
377  {
378  m_pressure->HomogeneousBwdTrans(Forcing[i],m_gradP[i]);
379  }
380  }
381  else
382  {
383  for (int i=0; i<nvel; i++)
384  {
385  Vmath::Vcopy(physTot, Forcing[i], 1, m_gradP[i], 1);
386  }
387  }
388 
389  if ( (!m_mapping->HasConstantJacobian()) || m_implicitPressure)
390  {
391  // If pressure terms are treated explicitly, we need to divide by J
392  // if they are implicit, we need to calculate G(p)
393  if (m_implicitPressure)
394  {
395  m_mapping->RaiseIndex(m_gradP, Forcing);
396  }
397  else
398  {
399  Array<OneD, NekDouble> Jac(physTot,0.0);
400  m_mapping->GetJacobian(Jac);
401  for (int i=0; i<nvel; i++)
402  {
403  Vmath::Vdiv(physTot, m_gradP[i], 1, Jac, 1, Forcing[i], 1);
404  }
405  }
406  // Transform back to wavespace
407  if (m_pressure->GetWaveSpace())
408  {
409  for (int i=0; i<nvel; i++)
410  {
411  m_pressure->HomogeneousFwdTrans(Forcing[i],Forcing[i]);
412  }
413  }
414  }
415 
416  // Subtract inarray/(aii_dt) and divide by kinvis. Kinvis will
417  // need to be updated for the convected fields.
418  for(int i = 0; i < nvel; ++i)
419  {
420  Blas::Daxpy(physTot,-aii_dtinv,inarray[i],1,Forcing[i],1);
421  Blas::Dscal(physTot,1.0/m_kinvis,&(Forcing[i])[0],1);
422  }
423  }
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:182
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:167
double NekDouble

References Blas::Daxpy(), Blas::Dscal(), 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().

◆ v_SolvePressure()

void Nektar::VCSMapping::v_SolvePressure ( const Array< OneD, NekDouble > &  Forcing)
protectedvirtual

Solve pressure system

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 428 of file VCSMapping.cpp.

430  {
431  if (!m_implicitPressure)
432  {
434  }
435  else
436  {
437  int physTot = m_fields[0]->GetTotPoints();
438  int nvel = m_nConvectiveFields;
439  bool converged = false; // flag to mark if system converged
440  int s = 0; // iteration counter
441  NekDouble error; // L2 error at current iteration
442  NekDouble forcing_L2 = 0.0; // L2 norm of F
443 
444  int maxIter;
445  m_session->LoadParameter("MappingMaxIter",maxIter,5000);
446 
447  // rhs of the equation at current iteration
448  Array< OneD, NekDouble> F_corrected(physTot, 0.0);
449  // Pressure field at previous iteration
450  Array<OneD, NekDouble> previous_iter (physTot, 0.0);
451  // Temporary variables
452  Array<OneD, Array<OneD, NekDouble> > wk1(nvel);
453  Array<OneD, Array<OneD, NekDouble> > wk2(nvel);
454  Array<OneD, Array<OneD, NekDouble> > gradP(nvel);
455  for(int i = 0; i < nvel; ++i)
456  {
457  wk1[i] = Array<OneD, NekDouble> (physTot, 0.0);
458  wk2[i] = Array<OneD, NekDouble> (physTot, 0.0);
459  gradP[i] = Array<OneD, NekDouble> (physTot, 0.0);
460  }
461 
462  // Jacobian
463  Array<OneD, NekDouble> Jac(physTot, 0.0);
464  m_mapping->GetJacobian(Jac);
465 
466  // Factors for Laplacian system
468  factors[StdRegions::eFactorLambda] = 0.0;
469 
470  m_pressure->BwdTrans(m_pressure->GetCoeffs(),
471  m_pressure->UpdatePhys());
472  forcing_L2 = m_pressure->L2(Forcing, wk1[0]);
473  while (!converged)
474  {
475  // Update iteration counter and set previous iteration field
476  // (use previous timestep solution for first iteration)
477  s++;
478  ASSERTL0(s < maxIter,
479  "VCSMapping exceeded maximum number of iterations.");
480 
481  Vmath::Vcopy(physTot, m_pressure->GetPhys(), 1,
482  previous_iter, 1);
483 
484  // Correct pressure bc to account for iteration
485  m_extrapolation->CorrectPressureBCs(previous_iter);
486 
487  //
488  // Calculate forcing term for this iteration
489  //
490  for(int i = 0; i < nvel; ++i)
491  {
493  previous_iter, gradP[i]);
494  if(m_pressure->GetWaveSpace())
495  {
496  m_pressure->HomogeneousBwdTrans(gradP[i], wk1[i]);
497  }
498  else
499  {
500  Vmath::Vcopy(physTot, gradP[i], 1, wk1[i], 1);
501  }
502  }
503  m_mapping->RaiseIndex(wk1, wk2); // G(p)
504 
505  m_mapping->Divergence(wk2, F_corrected); // div(G(p))
506  if (!m_mapping->HasConstantJacobian())
507  {
508  Vmath::Vmul(physTot, F_corrected, 1,
509  Jac, 1,
510  F_corrected, 1);
511  }
512  // alpha*J*div(G(p))
513  Vmath::Smul(physTot, m_pressureRelaxation, F_corrected, 1,
514  F_corrected, 1);
515  if(m_pressure->GetWaveSpace())
516  {
517  m_pressure->HomogeneousFwdTrans(F_corrected, F_corrected);
518  }
519  // alpha*J*div(G(p)) - p_ii
520  for (int i = 0; i < m_nConvectiveFields; ++i)
521  {
523  gradP[i], wk1[0]);
524  Vmath::Vsub(physTot, F_corrected, 1, wk1[0], 1,
525  F_corrected, 1);
526  }
527  // p_i,i - J*div(G(p))
528  Vmath::Neg(physTot, F_corrected, 1);
529  // alpha*F - alpha*J*div(G(p)) + p_i,i
530  Vmath::Smul(physTot, m_pressureRelaxation, Forcing, 1,
531  wk1[0], 1);
532  Vmath::Vadd(physTot, wk1[0], 1, F_corrected, 1, F_corrected, 1);
533 
534  //
535  // Solve system
536  //
537  m_pressure->HelmSolve(F_corrected, m_pressure->UpdateCoeffs(),
538  factors);
539  m_pressure->BwdTrans(m_pressure->GetCoeffs(),
540  m_pressure->UpdatePhys());
541 
542  //
543  // Test convergence
544  //
545  error = m_pressure->L2(m_pressure->GetPhys(), previous_iter);
546  if ( forcing_L2 != 0)
547  {
548  if ( (error/forcing_L2 < m_pressureTolerance))
549  {
550  converged = true;
551  }
552  }
553  else
554  {
555  if ( error < m_pressureTolerance)
556  {
557  converged = true;
558  }
559  }
560  }
561  if (m_verbose && m_session->GetComm()->GetRank()==0)
562  {
563  std::cout << " Pressure system (mapping) converged in " << s <<
564  " iterations with error = " << error << std::endl;
565  }
566  }
567  }
virtual void v_SolvePressure(const Array< OneD, NekDouble > &Forcing)
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314

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(), Vmath::Smul(), Nektar::VelocityCorrectionScheme::v_SolvePressure(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vsub().

◆ v_SolveViscous()

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 572 of file VCSMapping.cpp.

576  {
577  if(!m_implicitViscous)
578  {
579  VelocityCorrectionScheme::v_SolveViscous(Forcing, outarray, aii_Dt);
580  }
581  else
582  {
583  int physTot = m_fields[0]->GetTotPoints();
584  int nvel = m_nConvectiveFields;
585  bool converged = false; // flag to mark if system converged
586  int s = 0; // iteration counter
587  NekDouble error, max_error; // L2 error at current iteration
588 
589  int maxIter;
590  m_session->LoadParameter("MappingMaxIter",maxIter,5000);
591 
592  //L2 norm of F
593  Array<OneD, NekDouble> forcing_L2(m_nConvectiveFields,0.0);
594 
595  // rhs of the equation at current iteration
596  Array<OneD, Array<OneD, NekDouble> > F_corrected(nvel);
597  // Solution at previous iteration
598  Array<OneD, Array<OneD, NekDouble> > previous_iter(nvel);
599  // Working space
600  Array<OneD, Array<OneD, NekDouble> > wk(nvel);
601  for(int i = 0; i < nvel; ++i)
602  {
603  F_corrected[i] = Array<OneD, NekDouble> (physTot, 0.0);
604  previous_iter[i] = Array<OneD, NekDouble> (physTot, 0.0);
605  wk[i] = Array<OneD, NekDouble> (physTot, 0.0);
606  }
607 
608  // Factors for Helmholtz system
610  factors[StdRegions::eFactorLambda] =
611  1.0*m_viscousRelaxation/aii_Dt/m_kinvis;
612  if(m_useSpecVanVisc)
613  {
617  }
618 
619  // Calculate L2-norm of F and set initial solution for iteration
620  for(int i = 0; i < nvel; ++i)
621  {
622  forcing_L2[i] = m_fields[0]->L2(Forcing[i],wk[0]);
623  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
624  previous_iter[i]);
625  }
626 
627  while (!converged)
628  {
629  converged = true;
630  // Iteration counter
631  s++;
632  ASSERTL0(s < maxIter,
633  "VCSMapping exceeded maximum number of iterations.");
634 
635  max_error = 0.0;
636 
637  //
638  // Calculate forcing term for next iteration
639  //
640 
641  // Calculate L(U)- in this parts all components might be coupled
642  if(m_fields[0]->GetWaveSpace())
643  {
644  for (int i = 0; i < nvel; ++i)
645  {
646  m_fields[0]->HomogeneousBwdTrans(previous_iter[i],
647  wk[i]);
648  }
649  }
650  else
651  {
652  for (int i = 0; i < nvel; ++i)
653  {
654  Vmath::Vcopy(physTot, previous_iter[i], 1, wk[i], 1);
655  }
656  }
657 
658  // (L(U^i) - 1/alpha*U^i_jj)
659  m_mapping->VelocityLaplacian(wk, F_corrected,
660  1.0/m_viscousRelaxation);
661 
662  if(m_fields[0]->GetWaveSpace())
663  {
664  for (int i = 0; i < nvel; ++i)
665  {
666  m_fields[0]->HomogeneousFwdTrans(F_corrected[i],
667  F_corrected[i]);
668  }
669  }
670  else
671  {
672  for (int i = 0; i < nvel; ++i)
673  {
674  Vmath::Vcopy(physTot, F_corrected[i], 1,
675  F_corrected[i], 1);
676  }
677  }
678 
679  // Loop velocity components
680  for (int i = 0; i < nvel; ++i)
681  {
682  // (-alpha*L(U^i) + U^i_jj)
683  Vmath::Smul(physTot, -1.0*m_viscousRelaxation,
684  F_corrected[i], 1,
685  F_corrected[i], 1);
686  // F_corrected = alpha*F + (-alpha*L(U^i) + U^i_jj)
687  Vmath::Smul(physTot, m_viscousRelaxation, Forcing[i], 1,
688  wk[0], 1);
689  Vmath::Vadd(physTot, wk[0], 1, F_corrected[i], 1,
690  F_corrected[i], 1);
691 
692  //
693  // Solve System
694  //
695  m_fields[i]->HelmSolve(F_corrected[i],
696  m_fields[i]->UpdateCoeffs(),factors);
697  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),outarray[i]);
698 
699  //
700  // Test convergence
701  //
702  error = m_fields[i]->L2(outarray[i], previous_iter[i]);
703 
704  if ( forcing_L2[i] != 0)
705  {
706  if ( (error/forcing_L2[i] >= m_viscousTolerance))
707  {
708  converged = false;
709  }
710  }
711  else
712  {
713  if ( error >= m_viscousTolerance)
714  {
715  converged = false;
716  }
717  }
718  if (error > max_error)
719  {
720  max_error = error;
721  }
722 
723  // Copy field to previous_iter
724  Vmath::Vcopy(physTot, outarray[i], 1, previous_iter[i], 1);
725  }
726  }
727  if (m_verbose && m_session->GetComm()->GetRank()==0)
728  {
729  std::cout << " Velocity system (mapping) converged in " << s <<
730  " iterations with error = " << max_error << std::endl;
731  }
732  }
733  }
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
NekDouble m_sVVDiffCoeff
Diffusion coefficient of SVV modes.
bool m_useSpecVanVisc
bool to identify if spectral vanishing viscosity is active.

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, Vmath::Smul(), Nektar::VelocityCorrectionScheme::v_SolveViscous(), Vmath::Vadd(), and Vmath::Vcopy().

Member Data Documentation

◆ className

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

Name of class.

Definition at line 59 of file VCSMapping.h.

◆ m_gradP

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

Definition at line 94 of file VCSMapping.h.

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

◆ m_implicitPressure

bool Nektar::VCSMapping::m_implicitPressure
protected

◆ m_implicitViscous

bool Nektar::VCSMapping::m_implicitViscous
protected

◆ m_mapping

GlobalMapping::MappingSharedPtr Nektar::VCSMapping::m_mapping
protected

◆ m_neglectViscous

bool Nektar::VCSMapping::m_neglectViscous
protected

Definition at line 85 of file VCSMapping.h.

Referenced by ApplyIncNSMappingForcing(), and v_InitObject().

◆ m_presForcingCorrection

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

Definition at line 123 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SetUpPressureForcing().

◆ m_pressureRelaxation

NekDouble Nektar::VCSMapping::m_pressureRelaxation
protected

Definition at line 90 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

◆ m_pressureTolerance

NekDouble Nektar::VCSMapping::m_pressureTolerance
protected

Definition at line 88 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

◆ m_verbose

bool Nektar::VCSMapping::m_verbose
protected

Definition at line 79 of file VCSMapping.h.

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

◆ m_viscousRelaxation

NekDouble Nektar::VCSMapping::m_viscousRelaxation
protected

Definition at line 91 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().

◆ m_viscousTolerance

NekDouble Nektar::VCSMapping::m_viscousTolerance
protected

Definition at line 89 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().