45 "UnsteadyViscousBurgers",
63 AdvectionSystem::v_InitObject();
89 ASSERTL0(
false,
"Need to implement for DG");
96 m_session->LoadSolverInfo(
"AdvectionType", advName,
"WeakDG");
98 CreateInstance(advName, advName);
101 m_session->LoadSolverInfo(
"UpwindType", riemName,
"Upwind");
103 CreateInstance(riemName);
108 std::string diffName;
109 m_session->LoadSolverInfo(
"DiffusionType", diffName,
"LDG");
111 CreateInstance(diffName, diffName);
123 m_session->LoadSolverInfo(
"AdvectionType", advName,
126 CreateInstance(advName, advName);
133 for(
int i = 0; i <
m_fields.num_elements(); ++i)
143 ASSERTL0(
false,
"Explicit Galerkin diffusion not set up.");
150 ASSERTL0(
false,
"Unsupported projection type.");
194 for (i = 0; i < inarray.num_elements(); ++i)
196 m_fields[0]->ExtractTracePhys(inarray[i], tmp);
222 int nVariables = inarray.num_elements();
229 for (
int i = 0; i < nVariables; ++i)
236 inarray, outarray, time);
239 for (
int i = 0; i < nVariables; ++i)
249 for (
int i = 0; i < nVariables; ++i)
252 &outarrayDiff[i][0], 1, &outarray[i][0], 1);
257 std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x;
261 (*x)->Apply(
m_fields, inarray, outarray, time);
279 int nvariables = inarray.num_elements();
288 for(i = 0; i < nvariables; ++i)
301 ASSERTL0(
false,
"Unknown projection scheme");
321 int nvariables = inarray.num_elements();
336 for (
int n = 1; n < nvariables; ++n)
345 for (
int i = 0; i < nvariables; ++i)
349 inarray[i], 1, F[i], 1);
362 for(
int i = 0; i <
m_fields.num_elements(); ++i)
364 m_fields[i]->ClearGlobalLinSysManager();
371 for (
int i = 0; i < nvariables; ++i)
392 const int nq =
m_fields[0]->GetNpoints();
394 for (
int i = 0; i < flux.num_elements(); ++i)
396 for (
int j = 0; j < flux[0].num_elements(); ++j)
420 for (
int k = 0; k < flux.num_elements(); ++k)
430 AdvectionSystem::v_GenerateSummary(s);
void GetFluxVectorDiff(const int i, const int j, const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &derivatives, Array< OneD, Array< OneD, NekDouble > > &flux)
Evaluate the flux at each solution point for the diffusion part.
#define ASSERTL0(condition, msg)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
std::vector< std::pair< std::string, std::string > > SummaryList
NekDouble m_sVVCutoffRatio
DiffusionFactory & GetDiffusionFactory()
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
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
Array< OneD, NekDouble > m_traceVn
std::map< ConstFactorType, NekDouble > ConstFactorMap
StdRegions::VarCoeffMap m_varCoeffLap
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise...
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class.
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 t...
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
virtual ~UnsteadyViscousBurgers()
Destructor.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
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)
UnsteadyViscousBurgers(const LibUtilities::SessionReaderSharedPtr &pSession)
Session reader.
Base class for unsteady solvers.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
RiemannSolverFactory & GetRiemannSolverFactory()
void GetFluxVectorAdv(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point for the advection part.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Perform the projection.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
void Neg(int n, T *x, const int incx)
Negate x = -x.
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print Summary.
EquationSystemFactory & GetEquationSystemFactory()
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
virtual void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
Array< OneD, NekDouble > & GetNormalVelocity(Array< OneD, Array< OneD, NekDouble > > &inarray)
Get the normal velocity.
SolverUtils::DiffusionSharedPtr m_diffusion
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
virtual void v_InitObject()
Initialise the object.
virtual void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
Solve implicitly the diffusion term.
void Zero(int n, T *x, const int incx)
Zero vector.
A base class for PDEs which include an advection component.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
bool m_useSpecVanViscVarDiff
static std::string className
Name of class.
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.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
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.
static FlagList NullFlagList
An empty flag list.
static VarCoeffMap NullVarCoeffMap
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.