54 std::string variant,
size_t order, std::vector<NekDouble> freeParams)
56 m_name(
"FractionalInTime")
64 "FractionalInTime Time integration scheme bad order: " +
65 std::to_string(order));
68 freeParams.size() == 1 ||
69 freeParams.size() == 2 ||
70 freeParams.size() == 6,
71 "FractionalInTime Time integration scheme invalid number "
72 "of free parameters, expected zero, one <alpha>, "
73 "two <alpha, base>, or "
74 "six <alpha, base, nQuadPts, sigma, mu0, nu> received " +
75 std::to_string(freeParams.size()));
77 if (freeParams.size() >= 1)
82 if (freeParams.size() >= 2)
87 if (freeParams.size() == 6)
91 m_mu0 = freeParams[4];
133 for (
size_t m = 0; m <=
m_order; ++m)
137 for (
size_t i = 0; i <
m_nvars; ++i)
164 for (
size_t i = 0; i <
m_nvars; ++i)
176 for (
size_t m = 1, m_1 = 0; m <
m_order; ++m, ++m_1)
185 for (
size_t m = 1; m <=
m_order; ++m)
189 for (
size_t n = 0; n < m; ++n)
245 for (
size_t j = 2; j <= m; ++j)
247 for (
size_t i = 0; i < m; ++i)
249 m_Ahats[m][j - 1][i] = pow((1 - j), i);
253 ASSERTL1(
false,
"No matrix inverse.");
264 for (
size_t i = 0; i <
m_order; ++i)
266 for (
size_t j = 0; j <
m_order; ++j)
274 for (
size_t l = 0; l <
m_Lmax; ++l)
284 const size_t timestep, [[maybe_unused]]
const NekDouble delta_t)
287 "Delta T has changed which is not permitted.");
290 size_t timeStep = timestep + 1;
294 for (
size_t l = 0; l <
m_Lmax; ++l)
306 for (
size_t l = 0; l <
L; ++l)
312 for (
size_t i = 0; i <
m_nvars; ++i)
322 for (
size_t m =
m_order; m > 0; --m)
324 for (
size_t i = 0; i <
m_nvars; ++i)
328 m_u[m][i][j] =
m_u[m - 1][i][j];
334 for (
size_t i = 0; i <
m_nvars; ++i)
346 for (
size_t i = 0; i <
m_Lmax; ++i)
359 const size_t base)
const
361 return (counter + 1) % base;
369 const size_t l)
const
371 size_t L = ceil(
log(l / 2.0) /
log(base));
373 if (l % (
size_t)(2 * pow(base,
L)) == 0)
396 for (
size_t i = 0; i <
L - 1; ++i)
398 m_qml[i] = floor(m / pow(base, i + 1)) - 1;
429 for (
size_t i = 1; i <
L; ++i)
459 for (
size_t q = 0;
q < nQuadPts; ++
q)
464 lamb[
q] = sigma + mu * th * std::complex<NekDouble>(1. / tan(th), nu);
466 w[
q] = std::complex<NekDouble>(0, -1. /
NekDouble(nQuadPts)) * mu *
467 std::complex<NekDouble>(1. / tan(th) - th / (sin(th) * sin(th)),
473 if (nQuadPts % 2 == 1)
475 size_t q = (nQuadPts + 1) / 2;
477 lamb[
q] = std::complex<NekDouble>(sigma + mu, 0);
479 w[
q] = std::complex<NekDouble>(nu * mu / nQuadPts, 0);
487 const size_t index,
Instance &instance)
const
522 instance.
index = index;
554 std::pair<size_t, size_t>(0, 0);
573 std::pair<size_t, size_t>(0, 0);
580 instance.
csandbox_ind = std::pair<size_t, size_t>(0, 0);
584 instance.
fstash_ind = std::pair<size_t, size_t>(0, 0);
590 instance.
fsandbox_ind = std::pair<size_t, size_t>(0, 0);
622 for (
size_t m = 1; m <=
m_order + 1; ++m)
626 for (
size_t n = 0; n < m; ++n)
648 As[m][1][0] = 1. / 2.;
650 As[m][1][2] = -1. / 2.;
651 As[m][2][0] = 1. / 2.;
653 As[m][2][2] = 1. / 2.;
662 As[m][1][0] = 1. / 3.;
663 As[m][1][1] = 1. / 2.;
665 As[m][1][3] = 1. / 6.;
667 As[m][2][0] = 1. / 2.;
669 As[m][2][2] = 1. / 2.;
672 As[m][3][0] = 1. / 6.;
673 As[m][3][1] = -1. / 2.;
674 As[m][3][2] = 1. / 2.;
675 As[m][3][3] = -1. / 6.;
685 As[m][1][0] = 1. / 4.;
686 As[m][1][1] = 5. / 6.;
687 As[m][1][2] = -3. / 2.;
688 As[m][1][3] = 1. / 2.;
689 As[m][1][4] = -1. / 12.;
691 As[m][2][0] = 11. / 24.;
692 As[m][2][1] = -5. / 6.;
693 As[m][2][2] = 1. / 4.;
694 As[m][2][3] = 1. / 6.;
695 As[m][2][4] = -1. / 24.;
697 As[m][3][0] = 1. / 4.;
698 As[m][3][1] = -5. / 6.;
700 As[m][3][3] = -1. / 2.;
701 As[m][3][4] = 1. / 12.;
703 As[m][4][0] = 1. / 24.;
704 As[m][4][1] = -1. / 6.;
705 As[m][4][2] = 1. / 4.;
706 As[m][4][3] = -1. / 6.;
707 As[m][4][4] = 1. / 24.;
715 ASSERTL1(
false,
"No matrix inverse.");
730 for (
size_t m = 0; m <
m_order + 1; ++m)
739 1. / instance.
z[
q] * (exp(instance.
z[
q] *
m_deltaT) - 1.0);
743 instance.
Eh[m][
q] = -1. / instance.
z[
q] +
745 instance.
Eh[m - 1][
q];
755 for (
size_t m = 0; m <=
m_order; ++m)
761 for (
size_t i = 0; i <=
m_order; ++i)
763 instance.
AtEh[m][
q] +=
805 instance.
fsandbox_ind = std::pair<size_t, size_t>(0, 0);
809 for (
size_t i = 0; i <
m_nvars; ++i)
845 for (
size_t m = 0; m <
m_order; ++m)
849 for (
size_t i = 0; i <
m_nvars; ++i)
864 const size_t timeStep,
const size_t tauml,
const Instance &instance)
874 for (
size_t i = 0; i <
m_nvars; ++i)
885 if (
m_uInt[i][j].real() < 1e8)
912 for (
size_t m = 0; m <= order; ++m)
916 instance.
AtEh[m][
q] = 0;
918 for (
size_t i = 0; i <= order; ++i)
920 instance.
AtEh[m][
q] +=
921 instance.
As[order + 1][m][i] * instance.
Eh[i][
q];
932 for (
size_t m = 0; m <= order; ++m)
936 for (
size_t i = 0; i <
m_nvars; ++i)
945 y[i][j][
q] *= instance.
E[
q];
949 y[i][j][
q] +=
m_F[i][j] * instance.
AtEh[m][
q];
987 for (
size_t i = 0; i <
m_nvars; ++i)
1017 for (
size_t i = 0; i <
m_nvars; ++i)
1037 std::pair<size_t, size_t>(timeStep, timeStep);
1047 os <<
"Time Integration Scheme: " <<
GetFullName() << std::endl
1048 <<
" Alpha " <<
m_alpha << std::endl
1049 <<
" Base " <<
m_base << std::endl
1050 <<
" Number of instances " <<
m_Lmax << std::endl
1051 <<
" Number of quadature points " <<
m_nQuadPts << std::endl
1052 <<
" Talbot Parameter: sigma " <<
m_sigma << std::endl
1053 <<
" Talbot Parameter: mu0 " <<
m_mu0 << std::endl
1054 <<
" Talbot Parameter: nu " <<
m_nu << std::endl;
1059 os <<
"Time Integration Scheme: " <<
GetFullName() << std::endl
1060 <<
" Alpha " <<
m_alpha << std::endl
1061 <<
" Base " <<
m_base << std::endl
1062 <<
" Number of instances " <<
m_Lmax << std::endl
1063 <<
" Number of quadature points " <<
m_nQuadPts << std::endl
1064 <<
" Talbot Parameter: sigma " <<
m_sigma << std::endl
1065 <<
" Talbot Parameter: mu0 " <<
m_mu0 << std::endl
1066 <<
" Talbot Parameter: nu " <<
m_nu << std::endl;
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Class for fractional-in-time integration.
size_t computeTaus(const size_t base, const size_t m)
Method to compute the demarcation interval marker tau_{m, ell}.
Array< OneD, size_t > m_taus
LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Worker method to initialize the integration scheme.
void timeAdvance(const size_t timeStep, Instance &instance, ComplexTripleArray &y)
Method to get the solution to y' = z*y + f(u), using an exponential integrator with implicit order (m...
TimeIntegrationSchemeOperators m_op
Array< OneD, size_t > m_qml
void updateStage(const size_t timeStep, Instance &instance)
Method to rearrange of staging/stashing for current time.
void advanceSandbox(const size_t timeStep, Instance &instance)
Method to update sandboxes to the current time.
void integralClassInitialize(const size_t index, Instance &instance) const
Method to initialize the integral class.
LUE void v_print(std::ostream &os) const override
Worker method to print details on the integration scheme.
size_t modIncrement(const size_t counter, const size_t base) const
Method that increments the counter then performs mod calculation.
ComplexSingleArray m_expFactor
void finalIncrement(const size_t timeStep)
Method to approximate the integral.
ComplexDoubleArray m_uInt
void integralContribution(const size_t timeStep, const size_t tauml, const Instance &instance)
Method to get the integral contribution over [taus(i+1) taus(i)]. Stored in m_uInt.
void talbotQuadrature(const size_t nQuadPts, const NekDouble mu, const NekDouble nu, const NekDouble sigma, ComplexSingleArray &lamb, ComplexSingleArray &w) const
Method to compute the quadrature rule over Tablot contour.
size_t computeL(const size_t base, const size_t m) const
Method to compute the smallest integer L such that base < 2.
Array< OneD, Instance > m_integral_classes
LUE void v_printFull(std::ostream &os) const override
LUE ConstDoubleArray & v_TimeIntegrate(const size_t timestep, const NekDouble delta_t) override
Worker method that performs the time integration.
FractionalInTimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
Constructor.
size_t computeQML(const size_t base, const size_t m)
Method to compute the demarcation integers q_{m, ell}.
std::vector< NekDouble > m_freeParams
Base class for time integration schemes.
LUE void print(std::ostream &os) const
LUE std::string GetFullName() const
Binds a set of functions for use by time integration schemes.
void DoOdeRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
AT< AT< NekDouble > > DoubleArray
AT< AT< std::complex< NekDouble > > > ComplexDoubleArray
AT< NekDouble > SingleArray
AT< std::complex< NekDouble > > ComplexSingleArray
AT< AT< AT< NekDouble > > > TripleArray
std::shared_ptr< FractionalInTimeIntegrationScheme > FractionalInTimeIntegrationSchemeSharedPtr
AT< AT< AT< std::complex< NekDouble > > > > ComplexTripleArray
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
std::vector< double > w(NPUPPER)
std::vector< double > q(NPUPPER *NPUPPER)
scalarT< T > log(scalarT< T > in)
ComplexTripleArray fstash_y
std::pair< size_t, size_t > fsandbox_ind
ComplexTripleArray csandbox_y
ComplexTripleArray cstash_y
ComplexTripleArray stage_y
std::pair< size_t, size_t > cstash_ind
size_t fsandbox_activebase
std::pair< size_t, size_t > stage_ind
std::pair< size_t, size_t > fstash_ind
std::pair< size_t, size_t > csandbox_ind
size_t fsandbox_stashincrement
ComplexTripleArray fsandbox_y