35#include <boost/core/ignore_unused.hpp>
47 "Image warping system.");
97 m_session->LoadSolverInfo(
"AdvectionType", advName,
"WeakDG");
101 m_session->LoadSolverInfo(
"UpwindType", riemName,
"Upwind");
117 ASSERTL0(
false,
"Implicit unsteady Advection not set up.");
135 boost::ignore_unused(time);
138 int ncoeffs = inarray[0].size();
145 "CG not implemented yet.");
155 m_fields[0]->PhysDeriv(inarray[0], dIdx1, dIdx2);
163 for (
int i = 0; i < 2; ++i)
165 Vmath::Vmul(npoints, &alloc[i * npoints], 1, inarray[1].get(), 1,
166 m_fields[i + 2]->UpdatePhys().get(), 1);
168 m_fields[i + 2]->GetPhys().get(), 1,
169 m_fields[i + 2]->UpdatePhys().get(), 1);
178 for (
int i = 0; i < 2; ++i)
188 Vmath::Vvtvvtp(npoints, dIdx1.get(), 1, inarray[0].get(), 1, dIdx2.get(), 1,
189 inarray[0].get(), 1, dIdx1.get(), 1);
193 m_fields[0]->IProductWRTBase(dIdx1, tmp2);
196 m_fields[0]->MultiplyByElmtInvMass(tmp2, tmp2);
198 Vmath::Vadd(npoints, outarray[0], 1, tmp, 1, outarray[0], 1);
208 int nvariables = inarray.size();
216 if (inarray != outarray)
220 for (
int i = 0; i < nvariables; ++i)
228 ASSERTL0(
false,
"Unknown projection scheme");
266 "Dimension of flux array and velocity array do not match");
268 int nq = physfield[0].size();
270 for (
int i = 0; i < flux.size(); ++i)
272 for (
int j = 0; j < flux[0].size(); ++j)
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
ImageWarpingSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Array< OneD, Array< OneD, NekDouble > > m_velocity
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
static std::string className
Name of class.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
virtual void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
Array< OneD, NekDouble > m_traceVn
virtual ~ImageWarpingSystem()
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
A base class for PDEs which include an advection component.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT int GetNpoints()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
std::vector< int > m_intVariables
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ContField > ContFieldSharedPtr
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
std::vector< std::pair< std::string, std::string > > SummaryList
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::map< ConstFactorType, NekDouble > ConstFactorMap
StdRegions::ConstFactorMap factors
The above copyright notice and this permission notice shall be included.
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.
void Neg(int n, T *x, const int incx)
Negate x = -x.
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
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.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)