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

#include <VelocityCorrectionScheme.h>

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

Public Member Functions

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

Static Public Member Functions

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

Static Public Attributes

static std::string className
 Name of class.

Private Member Functions

virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
virtual void v_TransCoeffToPhys (void)
virtual void v_TransPhysToCoeff (void)
virtual void v_DoInitialise (void)
virtual Array< OneD, bool > v_GetSystemSingularChecks ()
virtual int v_GetForceDimension ()

Private Attributes

bool m_useHomo1DSpecVanVisc
 bool to identify if spectral vanishing viscosity is active.
bool m_useSpecVanVisc
 bool to identify if spectral vanishing viscosity is active.
NekDouble m_sVVCutoffRatio
 cutt off ratio from which to start decayhing modes
NekDouble m_sVVDiffCoeff
 Diffusion coefficient of SVV modes.

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1).
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous }
 Parameter for homogeneous expansions. More...
- Protected Member Functions inherited from Nektar::IncNavierStokes
 IncNavierStokes (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor.
EquationType GetEquationType (void)
void EvaluateAdvectionTerms (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
void WriteModalEnergy (void)
void SetBoundaryConditions (NekDouble time)
 time dependent boundary conditions updating
void SetRadiationBoundaryForcing (int fieldid)
 Set Radiation forcing term.
bool CalcSteadyState (void)
 evaluate steady state
virtual
MultiRegions::ExpListSharedPtr 
v_GetPressure ()
virtual void v_TransCoeffToPhys (void)
 Virtual function for transformation to physical space.
virtual void v_TransPhysToCoeff (void)
 Virtual function for transformation to coefficient space.
virtual int v_GetForceDimension ()=0
virtual bool v_PreIntegrate (int step)
virtual bool v_PostIntegrate (int step)
- Protected Attributes inherited from Nektar::IncNavierStokes
ExtrapolateSharedPtr m_extrapolation
std::ofstream m_mdlFile
 modal energy file
bool m_subSteppingScheme
 bool to identify if using a substepping scheme
bool m_SmoothAdvection
 bool to identify if advection term smoothing is requested
LibUtilities::TimeIntegrationWrapperSharedPtr m_subStepIntegrationScheme
AdvectionTermSharedPtr m_advObject
 Advection term.
std::vector
< SolverUtils::ForcingSharedPtr
m_forcing
 Forcing terms.
int m_nConvectiveFields
 Number of fields to be convected;.
Array< OneD, int > m_velocity
 int which identifies which components of m_fields contains the velocity (u,v,w);
MultiRegions::ExpListSharedPtr m_pressure
 Pointer to field holding pressure field.
NekDouble m_kinvis
 Kinematic viscosity.
int m_energysteps
 dump energy to file at steps time
int m_cflsteps
 dump cfl estimate
int m_steadyStateSteps
 Check for steady state at step interval.
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at.
EquationType m_equationType
 equation type;
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
 Mapping from BCs to Elmt IDs.
Array< OneD, Array< OneD, int > > m_fieldsBCToTraceID
 Mapping from BCs to Elmt Edge IDs.
Array< OneD, Array< OneD,
NekDouble > > 
m_fieldsRadiationFactor
 RHS Factor for Radiation Condition.
int m_intSteps
 Number of time integration steps AND Order of extrapolation for pressure boundary conditions.

Detailed Description

Definition at line 43 of file VelocityCorrectionScheme.h.

Constructor & Destructor Documentation

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

Constructor.

Constructor. Creates ...

Parameters
\param

Definition at line 55 of file VelocityCorrectionScheme.cpp.

:
IncNavierStokes(pSession)
{
}
Nektar::VelocityCorrectionScheme::~VelocityCorrectionScheme ( void  )
virtual

Destructor

Definition at line 158 of file VelocityCorrectionScheme.cpp.

{
}

Member Function Documentation

static SolverUtils::EquationSystemSharedPtr Nektar::VelocityCorrectionScheme::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 48 of file VelocityCorrectionScheme.h.

void Nektar::VelocityCorrectionScheme::EvaluateAdvection_SetPressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)

Explicit part of the method - Advection, Forcing + HOPBCs

Definition at line 274 of file VelocityCorrectionScheme.cpp.

References Nektar::IncNavierStokes::m_advObject, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_forcing, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::m_pressure, Nektar::IncNavierStokes::m_SmoothAdvection, Nektar::SolverUtils::EquationSystem::m_time, and Nektar::IncNavierStokes::m_velocity.

Referenced by v_InitObject().

{
// Evaluate convection terms
inarray, outarray, m_time);
// Smooth advection
{
for(int i = 0; i < m_nConvectiveFields; ++i)
{
m_pressure->SmoothField(outarray[i]);
}
}
// Add forcing terms
std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x;
for (x = m_forcing.begin(); x != m_forcing.end(); ++x)
{
(*x)->Apply(m_fields, inarray, outarray, time);
}
// Calculate High-Order pressure boundary conditions
m_extrapolation->EvaluatePressureBCs(inarray,outarray,m_kinvis);
}
void Nektar::VelocityCorrectionScheme::SetUpPressureForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
Array< OneD, Array< OneD, NekDouble > > &  Forcing,
const NekDouble  aii_Dt 
)

Forcing term for Poisson solver solver

Definition at line 358 of file VelocityCorrectionScheme.cpp.

References Nektar::MultiRegions::DirCartesianMap, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_velocity, Vmath::Smul(), Vmath::Vadd(), and Vmath::Zero().

Referenced by SolveUnsteadyStokesSystem().

{
int i;
int physTot = m_fields[0]->GetTotPoints();
int nvel = m_velocity.num_elements();
Array<OneD, NekDouble> wk(physTot, 0.0);
Vmath::Zero(physTot,Forcing[0],1);
for(i = 0; i < nvel; ++i)
{
m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[i],fields[i], wk);
Vmath::Vadd(physTot,wk,1,Forcing[0],1,Forcing[0],1);
}
Vmath::Smul(physTot,1.0/aii_Dt,Forcing[0],1,Forcing[0],1);
}
void Nektar::VelocityCorrectionScheme::SetUpViscousForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  Forcing,
const NekDouble  aii_Dt 
)

Forcing term for Helmholtz solver

Definition at line 381 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_pressure, and Nektar::IncNavierStokes::m_velocity.

Referenced by SolveUnsteadyStokesSystem().

{
NekDouble aii_dtinv = 1.0/aii_Dt;
int phystot = m_fields[0]->GetTotPoints();
// Grad p
m_pressure->BwdTrans(m_pressure->GetCoeffs(),m_pressure->UpdatePhys());
int nvel = m_velocity.num_elements();
if(nvel == 2)
{
m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1]);
}
else
{
m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1],
Forcing[2]);
}
// Subtract inarray/(aii_dt) and divide by kinvis. Kinvis will
// need to be updated for the convected fields.
for(int i = 0; i < nvel; ++i)
{
Blas::Daxpy(phystot,-aii_dtinv,inarray[i],1,Forcing[i],1);
Blas::Dscal(phystot,1.0/m_kinvis,&(Forcing[i])[0],1);
}
}
void Nektar::VelocityCorrectionScheme::SolveUnsteadyStokesSystem ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time,
const NekDouble  aii_Dt 
)

Implicit part of the method - Poisson + nConv*Helmholtz

Definition at line 306 of file VelocityCorrectionScheme.cpp.

References Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::m_pressure, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useSpecVanVisc, Nektar::NullFlagList, Nektar::IncNavierStokes::SetBoundaryConditions(), SetUpPressureForcing(), and SetUpViscousForcing().

Referenced by v_InitObject().

{
int i;
int phystot = m_fields[0]->GetTotPoints();
Array<OneD, Array< OneD, NekDouble> > F(m_nConvectiveFields);
for(i = 0; i < m_nConvectiveFields; ++i)
{
F[i] = Array<OneD, NekDouble> (phystot);
}
// Enforcing boundary conditions on all fields
// Substep the pressure boundary condition
m_extrapolation->SubStepSetPressureBCs(inarray,aii_Dt,m_kinvis);
// Set up forcing term and coefficients for pressure Poisson equation
SetUpPressureForcing(inarray, F, aii_Dt);
factors[StdRegions::eFactorLambda] = 0.0;
// Solver Pressure Poisson Equation
m_pressure->HelmSolve(F[0], m_pressure->UpdateCoeffs(), NullFlagList,
factors);
// Set up forcing term and coefficients for Helmholtz problems
SetUpViscousForcing(inarray, F, aii_Dt);
factors[StdRegions::eFactorLambda] = 1.0/aii_Dt/m_kinvis;
{
}
// Solve Helmholtz system and put in Physical space
for(i = 0; i < m_nConvectiveFields; ++i)
{
m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(),
NullFlagList, factors);
m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),outarray[i]);
}
}
void Nektar::VelocityCorrectionScheme::v_DoInitialise ( void  )
privatevirtual

Definition at line 204 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_nConvectiveFields, and Nektar::SolverUtils::EquationSystem::m_timestep.

{
// Set up Field Meta Data for output files
m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<std::string>(m_timestep);
for(int i = 0; i < m_nConvectiveFields; ++i)
{
m_fields[i]->LocalToGlobal();
m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
m_fields[i]->GlobalToLocal();
m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
m_fields[i]->UpdatePhys());
}
}
void Nektar::VelocityCorrectionScheme::v_GenerateSummary ( SolverUtils::SummaryList s)
privatevirtual

Definition at line 165 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::AddSummaryItem(), Nektar::IncNavierStokes::m_advObject, Nektar::SolverUtils::EquationSystem::m_homogen_dealiasing, Nektar::IncNavierStokes::m_subStepIntegrationScheme, Nektar::IncNavierStokes::m_subSteppingScheme, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useHomo1DSpecVanVisc, m_useSpecVanVisc, and Nektar::LibUtilities::TimeIntegrationMethodMap.

{
{
m_subStepIntegrationScheme->GetIntegrationMethod()]);
}
string dealias = m_homogen_dealiasing ? "Homogeneous1D" : "";
if (m_advObject->GetSpecHPDealiasing())
{
dealias += (dealias == "" ? "" : " + ") + string("spectral/hp");
}
if (dealias != "")
{
SolverUtils::AddSummaryItem(s, "Dealiasing", dealias);
}
string smoothing = m_useSpecVanVisc ? "spectral/hp" : "";
{
smoothing += (smoothing == "" ? "" : " + ") + string("Homogeneous1D");
}
if (smoothing != "")
{
s, "Smoothing", "SVV (" + smoothing + " SVV (cut-off = "
+ boost::lexical_cast<string>(m_sVVCutoffRatio)
+ ", diff coeff = "
+ boost::lexical_cast<string>(m_sVVDiffCoeff)+")");
}
}
int Nektar::VelocityCorrectionScheme::v_GetForceDimension ( void  )
privatevirtual

Definition at line 266 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::EquationSystem::m_session.

{
return m_session->GetVariables().size() - 1;
}
Array< OneD, bool > Nektar::VelocityCorrectionScheme::v_GetSystemSingularChecks ( )
privatevirtual

Definition at line 255 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::EquationSystem::m_session.

{
int vVar = m_session->GetVariables().size();
Array<OneD, bool> vChecks(vVar, false);
vChecks[vVar-1] = true;
return vChecks;
}
void Nektar::VelocityCorrectionScheme::v_InitObject ( )
virtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::IncNavierStokes.

Definition at line 62 of file VelocityCorrectionScheme.cpp.

References ASSERTL0, Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::MultiRegions::eMixed_CG_Discontinuous, EvaluateAdvection_SetPressureBCs(), Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, Nektar::SolverUtils::UnsteadySystem::m_intScheme, Nektar::SolverUtils::UnsteadySystem::m_intVariables, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::IncNavierStokes::m_pressure, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_session, Nektar::IncNavierStokes::m_SmoothAdvection, Nektar::IncNavierStokes::m_subSteppingScheme, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useHomo1DSpecVanVisc, m_useSpecVanVisc, Nektar::IncNavierStokes::m_velocity, and SolveUnsteadyStokesSystem().

{
int n;
// Set m_pressure to point to last field of m_fields;
if (boost::iequals(m_session->GetVariable(m_fields.num_elements()-1), "p"))
{
m_nConvectiveFields = m_fields.num_elements()-1;
}
else
{
ASSERTL0(false,"Need to set up pressure field definition");
}
// Integrate only the convective fields
for (n = 0; n < m_nConvectiveFields; ++n)
{
m_intVariables.push_back(n);
}
// Load parameters for Spectral Vanishing Viscosity
m_session->MatchSolverInfo("SpectralVanishingViscosity","True",m_useSpecVanVisc,false);
m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
// Needs to be set outside of next if so that it is turned off by default
m_session->MatchSolverInfo("SpectralVanishingViscosityHomo1D","True",m_useHomo1DSpecVanVisc,false);
{
ASSERTL0(m_nConvectiveFields > 2,"Expect to have three velcoity fields with homogenous expansion");
{
m_session->MatchSolverInfo("SpectralVanishingViscosity","True",m_useHomo1DSpecVanVisc,false);
}
{
Array<OneD, unsigned int> planes;
planes = m_fields[0]->GetZIDs();
int num_planes = planes.num_elements();
Array<OneD, NekDouble> SVV(num_planes,0.0);
NekDouble fac;
int kmodes = m_fields[0]->GetHomogeneousBasis()->GetNumModes();
int pstart;
pstart = m_sVVCutoffRatio*kmodes;
for(n = 0; n < num_planes; ++n)
{
if(planes[n] > pstart)
{
fac = (NekDouble)((planes[n] - kmodes)*(planes[n] - kmodes))/
((NekDouble)((planes[n] - pstart)*(planes[n] - pstart)));
SVV[n] = m_sVVDiffCoeff*exp(-fac)/m_kinvis;
}
}
for(int i = 0; i < m_velocity.num_elements(); ++i)
{
m_fields[m_velocity[i]]->SetHomo1DSpecVanVisc(SVV);
}
}
}
m_session->MatchSolverInfo("SmoothAdvection", "True", m_SmoothAdvection, false);
if(m_subSteppingScheme) // Substepping
{
"Projection must be set to Mixed_CG_Discontinuous for "
"substepping");
}
else // Standard velocity correction scheme
{
// set explicit time-intregration class operators
}
m_extrapolation->SubSteppingTimeIntegration(m_intScheme->GetIntegrationMethod(), m_intScheme);
m_extrapolation->GenerateHOPBCMap();
// set implicit time-intregration class operators
}
void Nektar::VelocityCorrectionScheme::v_TransCoeffToPhys ( void  )
privatevirtual

Definition at line 227 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields.

{
int nfields = m_fields.num_elements() - 1;
for (int k=0 ; k < nfields; ++k)
{
//Backward Transformation in physical space for time evolution
m_fields[k]->BwdTrans_IterPerExp(m_fields[k]->GetCoeffs(),
m_fields[k]->UpdatePhys());
}
}
void Nektar::VelocityCorrectionScheme::v_TransPhysToCoeff ( void  )
privatevirtual

Definition at line 241 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields.

{
int nfields = m_fields.num_elements() - 1;
for (int k=0 ; k < nfields; ++k)
{
//Forward Transformation in physical space for time evolution
m_fields[k]->FwdTrans_IterPerExp(m_fields[k]->GetPhys(),m_fields[k]->UpdateCoeffs());
}
}

Member Data Documentation

string Nektar::VelocityCorrectionScheme::className
static
Initial value:

Name of class.

Definition at line 58 of file VelocityCorrectionScheme.h.

NekDouble Nektar::VelocityCorrectionScheme::m_sVVCutoffRatio
private

cutt off ratio from which to start decayhing modes

Definition at line 93 of file VelocityCorrectionScheme.h.

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

NekDouble Nektar::VelocityCorrectionScheme::m_sVVDiffCoeff
private

Diffusion coefficient of SVV modes.

Definition at line 95 of file VelocityCorrectionScheme.h.

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

bool Nektar::VelocityCorrectionScheme::m_useHomo1DSpecVanVisc
private

bool to identify if spectral vanishing viscosity is active.

Definition at line 89 of file VelocityCorrectionScheme.h.

Referenced by v_GenerateSummary(), and v_InitObject().

bool Nektar::VelocityCorrectionScheme::m_useSpecVanVisc
private

bool to identify if spectral vanishing viscosity is active.

Definition at line 91 of file VelocityCorrectionScheme.h.

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