40 #include <boost/algorithm/string.hpp> 
   46             "VelocityCorrectionScheme", 
 
   57         : UnsteadySystem(pSession),
 
   77             ASSERTL0(
false,
"Need to set up pressure field definition");
 
   83             std::string vExtrapolation = 
"Standard";
 
   85             if (
m_session->DefinesSolverInfo(
"Extrapolation"))
 
   87                 vExtrapolation = 
m_session->GetSolverInfo(
"Extrapolation");
 
  116             ASSERTL0(m_nConvectiveFields > 2,
"Expect to have three velcoity fields with homogenous expansion");
 
  128                 int num_planes = planes.num_elements();
 
  131                 int kmodes = 
m_fields[0]->GetHomogeneousBasis()->GetNumModes();
 
  136                 for(n = 0; n < num_planes; ++n)
 
  138                     if(planes[n] > pstart)
 
  140                         fac = (
NekDouble)((planes[n] - kmodes)*(planes[n] - kmodes))/
 
  141                             ((
NekDouble)((planes[n] - pstart)*(planes[n] - pstart)));
 
  146                 for(
int i = 0; i < 
m_velocity.num_elements(); ++i)
 
  159                      "Projection must be set to Mixed_CG_Discontinuous for " 
  186         UnsteadySystem::v_GenerateSummary(s);
 
  198             dealias += (dealias == 
"" ? 
"" : 
" + ") + 
string(
"spectral/hp");
 
  208             smoothing += (smoothing == 
"" ? 
"" : 
" + ") + 
string(
"Homogeneous1D");
 
  213                 s, 
"Smoothing", 
"SVV (" + smoothing + 
" SVV (cut-off = " 
  226         UnsteadySystem::v_DoInitialise();
 
  248         int nfields = 
m_fields.num_elements() - 1;
 
  249         for (
int k=0 ; k < nfields; ++k)
 
  263         int nfields = 
m_fields.num_elements() - 1;
 
  264         for (
int k=0 ; k < nfields; ++k)
 
  276         int vVar = 
m_session->GetVariables().size();
 
  278         vChecks[vVar-1] = 
true;
 
  287         return m_session->GetVariables().size() - 1;
 
  310         std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x;
 
  313             (*x)->Apply(
m_fields, inarray, outarray, time);
 
  330         int phystot = 
m_fields[0]->GetTotPoints();
 
  381         int physTot = 
m_fields[0]->GetTotPoints();
 
  387         for(i = 0; i < nvel; ++i)
 
  390             Vmath::Vadd(physTot,wk,1,Forcing[0],1,Forcing[0],1);
 
  393         Vmath::Smul(physTot,1.0/aii_Dt,Forcing[0],1,Forcing[0],1);        
 
  405         int phystot = 
m_fields[0]->GetTotPoints();
 
  423         for(
int i = 0; i < nvel; ++i)
 
  425             Blas::Daxpy(phystot,-aii_dtinv,inarray[i],1,Forcing[i],1);
 
  426             Blas::Dscal(phystot,1.0/
m_kinvis,&(Forcing[i])[0],1);
 
EquationType m_equationType
equation type; 
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class. 
virtual Array< OneD, bool > v_GetSystemSingularChecks()
virtual void v_DoInitialise(void)
Sets up initial conditions. 
#define ASSERTL0(condition, msg)
virtual ~VelocityCorrectionScheme()
void SetBoundaryConditions(NekDouble time)
time dependent boundary conditions updating 
void EvaluateAdvection_SetPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
SolverUtils::AdvectionSharedPtr m_advObject
Advection term. 
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey. 
void SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
NekDouble m_kinvis
Kinematic viscosity. 
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
ExtrapolateFactory & GetExtrapolateFactory()
VelocityCorrectionScheme(const LibUtilities::SessionReaderSharedPtr &pSession)
Constructor. 
bool m_subSteppingScheme
bool to identify if using a substepping scheme 
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use. 
NekDouble m_timestep
Time step size. 
std::vector< std::pair< std::string, std::string > > SummaryList
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print a summary of time stepping parameters. 
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w); 
bool m_useHomo1DSpecVanVisc
bool to identify if spectral vanishing viscosity is active. 
NekDouble m_sVVDiffCoeff
Diffusion coefficient of SVV modes. 
ExtrapolateSharedPtr m_extrapolation
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous. 
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation. 
std::map< ConstFactorType, NekDouble > ConstFactorMap
static std::string className
Name of class. 
const char *const TimeIntegrationMethodMap[]
int m_nConvectiveFields
Number of fields to be convected;. 
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
virtual int v_GetForceDimension()
void SetUpViscousForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y. 
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms. 
void SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested 
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list. 
virtual void v_TransCoeffToPhys(void)
Virtual function for transformation to physical space. 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields. 
bool m_homogen_dealiasing
Flag to determine if dealiasing is used for homogeneous simulations. 
EquationSystemFactory & GetEquationSystemFactory()
MultiRegions::Direction const DirCartesianMap[]
void EvaluateAdvectionTerms(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
NekDouble m_sVVCutoffRatio
cutt off ratio from which to start decayhing modes 
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables. 
bool m_useSpecVanVisc
bool to identify if spectral vanishing viscosity is active. 
LibUtilities::SessionReaderSharedPtr m_session
The session reader. 
virtual void v_TransPhysToCoeff(void)
Virtual function for transformation to coefficient space. 
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme. 
This class is the base class for Navier Stokes problems. 
LibUtilities::TimeIntegrationWrapperSharedPtr m_subStepIntegrationScheme
virtual void v_InitObject()
Init object for UnsteadySystem class. 
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field. 
void Zero(int n, T *x, const int incx)
Zero vector. 
virtual void v_InitObject()
Init object for UnsteadySystem class. 
enum HomogeneousType m_HomogeneousType
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. 
static FlagList NullFlagList
An empty flag list. 
std::vector< int > m_intVariables
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.