35#ifndef NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_IMPLICIT_TIME_INTEGRATION_SCHEME_SDC
36#define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_IMPLICIT_TIME_INTEGRATION_SCHEME_SDC
38#define LUE LIB_UTILITIES_EXPORT
50 std::vector<NekDouble> freeParams)
54 "Quadrature type that include the left end point (e.g. "
55 "GaussLobattoLegendre) should not be used for ImplicitSDC");
57 m_name =
"ImplicitSpectralDeferredCorrection";
62 std::string variant,
size_t order, std::vector<NekDouble> freeParams)
66 variant, order, freeParams);
99 for (
size_t i = 0; i <
m_nvars; ++i)
110 for (
size_t i = 0; i <
m_nvars; ++i)
121 [[maybe_unused]]
const NekDouble &delta_t, [[maybe_unused]]
const size_t n)
155 for (
size_t i = 0; i <
m_nvars; ++i)
183 for (
size_t i = 0; i <
m_nvars; ++i)
199 for (
size_t i = 0; i <
m_nvars; ++i)
#define ASSERTL0(condition, msg)
LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Array for temporary storage.
LUE void v_ComputeInitialGuess(const NekDouble &delta_t) override
Worker method to compute the initial SDC guess.
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
ImplicitTimeIntegrationSchemeSDC(std::string variant, size_t order, std::vector< NekDouble > freeParams)
LUE void v_SDCIterationLoop(const NekDouble &delta_t) override
Worker method to compute the SDC iteration.
LUE void v_ResidualEval(const NekDouble &delta_t, const size_t n) override
Worker method to compute the residual.
static std::string className
Binds a set of functions for use by time integration schemes.
void DoImplicitSolve(InArrayType &inarray, OutArrayType &outarray, const NekDouble time, const NekDouble lambda) const
void DoProjection(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
void DoOdeRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
Class for spectral deferred correction integration.
LUE void UpdateIntegratedResidualSFint(const NekDouble &delta_t)
Worker method that compute residual integral.
TripleArray m_Y
Array containing the last stage values.
TimeIntegrationSchemeOperators m_op
size_t m_nQuadPts
Order of the integration scheme.
LUE void AddFASCorrectionToSFint(void)
Worker method that add the FASCorrection.
bool m_first_quadrature
Number of points in the integration scheme.
size_t m_npoints
Number of variables in the integration scheme.
TripleArray m_SFint
Array containing the FAS correction term.
NekDouble m_theta
Array containing the interpolation coefficients.
TimeIntegrationSchemeType m_schemeType
SingleArray m_tau
Object containing quadrature data.
TripleArray m_F
Array containing the stage values.
LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Worker method to initialize the integration scheme.
size_t m_nvars
Number of quadrature points.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
AT< AT< NekDouble > > DoubleArray
AT< NekDouble > SingleArray
@ eImplicit
Fully implicit scheme.
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
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.