35#ifndef NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_TIME_INTEGRATION_SCHEME_SDC
36#define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_TIME_INTEGRATION_SCHEME_SDC
38#define LUE LIB_UTILITIES_EXPORT
60 std::vector<NekDouble> freeParams)
62 m_name(
"SpectralDeferredCorrection")
65 "SDC Time integration scheme invalid number "
66 "of free parameters, expected two "
67 "<theta, number of quadrature>, received " +
68 std::to_string(freeParams.size()));
70 ASSERTL0(0.0 <= freeParams[0] && freeParams[0] <= 1.0,
71 "Spectral Deferred Correction Time integration "
72 "scheme bad parameter numbers (0.0 - 1.0): " +
73 std::to_string(freeParams[0]));
86 " quadrature require quadrature "
88 "): " + std::to_string(freeParams[1]));
97 else if (
m_variant ==
"GaussLobattoLegendre")
101 " quadrature require quadrature "
103 "): " + std::to_string(freeParams[1]));
112 else if (
m_variant ==
"GaussRadauLegendre")
116 " quadrature require quadrature "
118 "): " + std::to_string(freeParams[1]));
127 else if (
m_variant ==
"GaussGaussLegendre")
131 " quadrature require quadrature "
133 "): " + std::to_string(freeParams[1]));
144 ASSERTL0(
false,
"unknow variant (quadrature) type");
148 "Spectral Deferred Correction Time integration "
149 "scheme bad order numbers (>=" +
151 "): " + std::to_string(
m_order));
154 "Spectral Deferred Correction Time integration "
155 "scheme bad order numbers (<=" +
157 "): " + std::to_string(
m_order));
171 m_tau[i] = (tau + 1.0) / 2.0;
181 std::string variant,
size_t order, std::vector<NekDouble> freeParams)
185 variant, order, freeParams);
341 const size_t timestep,
const NekDouble delta_t)
override;
345 ASSERTL0(
false,
"Specific version of spectral deferred correction "
347 boost::ignore_unused(delta_t, n);
352 ASSERTL0(
false,
"Specific version of spectral deferred correction "
354 boost::ignore_unused(delta_t);
359 ASSERTL0(
false,
"Specific version of spectral deferred correction "
361 boost::ignore_unused(delta_t);
366 ASSERTL0(
false,
"Specific version of spectral deferred correction "
368 boost::ignore_unused(delta_t);
371 LUE virtual void v_print(std::ostream &os)
const override;
#define ASSERTL0(condition, msg)
Defines a specification for a set of points.
Base class for time integration schemes.
Binds a set of functions for use by time integration schemes.
Class for spectral deferred correction integration.
std::vector< NekDouble > m_freeParams
LUE TripleArray & UpdateIntegratedResidualQFintVector()
size_t m_ordermax
Minimum order of the integration scheme.
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
LUE void UpdateIntegratedResidualSFint(const NekDouble &delta_t)
Worker method that compute residual integral.
virtual LUE size_t v_GetNumIntegrationPhases() const override
LUE const DoubleArray & GetLastQuadratureSolutionVector() const
SingleArray m_interp
Array containing the integration matrix.
LUE void SetPFASST(bool pfasst)
virtual LUE void v_ComputeInitialGuess(const NekDouble &delta_t)
virtual LUE size_t v_GetOrder() const override
LUE const TripleArray & GetFAScorrectionVector() const
virtual LUE TimeIntegrationSchemeType v_GetIntegrationSchemeType() const override
virtual ~TimeIntegrationSchemeSDC()
Destructor.
TripleArray m_Y
Array containing the last stage values.
virtual LUE void v_print(std::ostream &os) const override
Worker method to print details on the integration scheme.
LUE size_t GetNpoints() const
TimeIntegrationSchemeOperators m_op
virtual LUE std::string v_GetName() const override
TripleArray m_QFint
Array containing the integrated residual term.
LUE void UpdateIntegratedResidualQFint(const NekDouble &delta_t)
size_t m_nQuadPts
Order of the integration scheme.
LUE size_t GetMaxOrder() const
LUE void AddFASCorrectionToSFint(void)
Worker method that add the FASCorrection.
virtual LUE TripleArray & v_UpdateSolutionVector() override
bool m_first_quadrature
Number of points in the integration scheme.
size_t m_ordermin
SDC parameter.
LUE DoubleArray & UpdateFirstQuadratureSolutionVector()
LUE const TripleArray & GetIntegratedResidualQFintVector() const
size_t m_npoints
Number of variables in the integration scheme.
SingleArray m_QMat
Array containing the integrated residual term.
LUE void ComputeInitialGuess(const NekDouble &delta_t)
virtual LUE void v_printFull(std::ostream &os) const override
LUE void SetTime(double time)
LUE bool HasFirstQuadrature() const
TripleArray m_FAScorr
Array containing the stage derivatives.
TripleArray m_SFint
Array containing the FAS correction term.
virtual LUE void v_SDCIterationLoop(const NekDouble &delta_t)
DoubleArray m_Y_f
Array containing the quadrature points.
LUE DoubleArray & UpdateLastQuadratureSolutionVector()
LUE TripleArray & UpdateFAScorrectionVector()
NekDouble m_theta
Array containing the interpolation coefficients.
virtual LUE std::vector< NekDouble > v_GetFreeParams() const override
TimeIntegrationSchemeType m_schemeType
LUE const TripleArray & GetResidualVector() const
LUE void ResidualEval(const NekDouble &delta_t)
size_t m_order
Maximum order of the integration scheme.
SingleArray m_tau
Object containing quadrature data.
LUE void UpdateFirstQuadrature(void)
Worker method that update the first quadrature.
TripleArray m_F
Array containing the stage values.
LUE const DoubleArray & GetFirstQuadratureSolutionVector() const
LUE size_t GetNvars() const
virtual LUE void v_ResidualEval(const NekDouble &delta_t, const size_t n)
virtual LUE void v_SetSolutionVector(const size_t Offset, const DoubleArray &y) override
Sets the solution vector of the ODE.
virtual LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Worker method to initialize the integration scheme.
LUE void UpdateLastQuadrature(void)
Worker method that update the last quadrature.
virtual LUE ConstDoubleArray & v_TimeIntegrate(const size_t timestep, const NekDouble delta_t) override
Worker method that performs the time integration.
LUE const PointsKey & GetPointsKey() const
virtual LUE std::string v_GetVariant() const override
virtual LUE void v_ResidualEval(const NekDouble &delta_t)
virtual LUE NekDouble v_GetTimeStability() const override
LUE TripleArray & UpdateResidualVector()
LUE void ResidualEval(const NekDouble &delta_t, const size_t n)
LUE size_t GetQuadPtsNumber() const
LUE void SDCIterationLoop(const NekDouble &delta_t)
TimeIntegrationSchemeSDC(std::string variant, size_t order, std::vector< NekDouble > freeParams)
LUE bool HasLastQuadrature() const
virtual LUE const TripleArray & v_GetSolutionVector() const override
Gets the solution vector of the ODE.
size_t m_nvars
Number of quadrature points.
static std::string className
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
PointsManagerT & PointsManager(void)
AT< NekDouble > SingleArray
TimeIntegrationSchemeType
@ eNoTimeIntegrationSchemeType
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
@ eGaussRadauPLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
The above copyright notice and this permission notice shall be included.