Lodestar
An integrated real-time control package in C++
ButcherTableau.hpp
1 //
2 // Created by Hamza El-Kebir on 6/18/21.
3 //
4 
5 #ifndef LODESTAR_BUTCHERTABLEAU_HPP
6 #define LODESTAR_BUTCHERTABLEAU_HPP
7 
8 #include <type_traits>
9 #include <tuple>
10 #include <array>
11 #include <string>
12 #include <string_view>
13 #include <functional>
14 #include "Lodestar/aux/Indices.hpp"
15 #include "Lodestar/aux/Conjunction.hpp"
16 #include <iostream>
17 #include <Eigen/Dense>
18 
19 namespace ls {
20  namespace primitives {
21  namespace detail {
34  template<typename TScalarType, size_t TStages, size_t TStage, bool TIsWeights>
35  struct ButcherRowImpl {
36  };
37 
47  template<typename TScalarType, size_t TStages, size_t TStage>
48  struct ButcherRowImpl<TScalarType, TStages, TStage, false> {
49  static_assert(TStage >= 0, "Butcher tableau row number must be non-negative.");
50  static_assert(TStage < (TStages - 1),
51  "Butcher tableau row number must be smaller than the number of stages minus one.");
52 
53  TScalarType node;
54  std::array<TScalarType, TStage + 1> rkCoefficients;
55  };
56 
66  template<typename TScalarType, size_t TStages, size_t TStage>
67  struct ButcherRowImpl<TScalarType, TStages, TStage, true> {
68  std::array<TScalarType, TStages> weights;
69  };
70  }
71 
88  template<size_t TStages, bool TExtended = true, typename TScalarType = double>
90  static_assert(TStages > 1, "Butcher tableau must have more than one stage.");
91  };
92 
93  template<size_t TStages, typename TScalarType>
94  class ButcherTableau<TStages, false, TScalarType> {
95  public:
96  static_assert(TStages > 1, "Butcher tableau must have more than one stage.");
97 
98  static const size_t stages = TStages;
99  using type = TScalarType;
100 
107  template<size_t TStage, bool TIsWeights = (TStage >= (TStages - 1))>
109 
118  template<size_t TRow, size_t TCoeffIdx>
119  typename std::enable_if<(TRow < (TStages - 1)) && (TCoeffIdx < (TRow + 1)), TScalarType>::type
120  inline getCoefficient()
121  {
122  return std::get<TCoeffIdx>(std::get<TRow>(rows_).rkCoefficients);
123  }
124 
133  template<size_t TRow, size_t TCoeffIdx>
134  typename std::enable_if<!((TRow < (TStages - 1)) && (TCoeffIdx < (TRow + 1))), TScalarType>::type
135  inline getCoefficient()
136  {
137  static_assert(TRow < (TStages - 1),
138  "Row index must be smaller than the number of stages minus one to access coefficients.");
139  static_assert(TCoeffIdx < (TRow + 1), "Coefficient index must be less than or equal to the row index.");
140 
141  return TScalarType{};
142  }
143 
152  template<size_t TRow, size_t TCoeffIdx>
153  typename std::enable_if<(TRow < (TStages - 1)) && (TCoeffIdx < (TRow + 1)), void>::type
154  inline setCoefficient(TScalarType coeff)
155  {
156  std::get<TCoeffIdx>(std::get<TRow>(rows_).rkCoefficients) = coeff;
157  }
158 
167  template<size_t TRow, size_t TCoeffIdx>
168  typename std::enable_if<!((TRow < (TStages - 1)) && (TCoeffIdx < (TRow + 1))), void>::type
169  inline setCoefficient(TScalarType coeff)
170  {
171  static_assert(TRow < (TStages - 1),
172  "Row index must be smaller than the number of stages minus one to access coefficients.");
173  static_assert(TCoeffIdx < (TRow + 1), "Coefficient index must be less than or equal to the row index.");
174  }
175 
183  template<size_t TRow>
184  typename std::enable_if<TRow < (TStages - 1), TScalarType>::type
185  inline getNode()
186  {
187  return std::get<TRow>(rows_).node;
188  }
189 
197  template<size_t TRow>
198  typename std::enable_if<TRow >= (TStages - 1), TScalarType>::type
199  inline getNode()
200  {
201  static_assert(TRow < (TStages - 1),
202  "Row index must be smaller than the number of stages minus one to access coefficients.");
203 
204  return TScalarType{};
205  }
206 
214  template<size_t TRow>
215  typename std::enable_if<TRow < (TStages - 1), void>::type
216  inline setNode(TScalarType node)
217  {
218  std::get<TRow>(rows_).node = node;
219  }
220 
228  template<size_t TRow>
229  typename std::enable_if<TRow >= (TStages - 1), void>::type
230  inline setNode(TScalarType node)
231  {
232  static_assert(TRow < (TStages - 1),
233  "Row index must be smaller than the number of stages minus one to access coefficients.");
234  }
235 
243  template<size_t TIdx>
244  typename std::enable_if<TIdx < TStages, TScalarType>::type
245  inline getWeight()
246  {
247  return std::get<TIdx>(std::get<TStages - 1>(rows_).weights);
248  }
249 
257  template<size_t TIdx>
258  typename std::enable_if<TIdx >= TStages, TScalarType>::type
259  inline getWeight()
260  {
261  static_assert(TIdx < TStages, "Weight index must be smaller than the number of stages.");
262 
263  return TScalarType{};
264  }
265 
273  template<size_t TIdx>
274  typename std::enable_if<TIdx < TStages, void>::type
275  inline setWeight(TScalarType weight)
276  {
277  std::get<TIdx>(std::get<TStages - 1>(rows_).weights) = weight;
278  }
279 
287  template<size_t TIdx>
288  typename std::enable_if<TIdx >= TStages, void>::type
289  inline setWeight(TScalarType weight)
290  {
291  static_assert(TIdx < TStages, "Weight index must be smaller than the number of stages.");
292  }
293 
312  template<typename TType, size_t TStage = 0>
313  typename std::enable_if<TStage == 0, TType>::type
314  inline execute(const std::function<TType(TScalarType, TType)> &f, const TType &y, const TScalarType t,
315  const TScalarType h)
316  {
317 // TType kCurr;
318 // std::cout << "f(t,y) in Butcher tableau: " << f(t, y) << std::endl;
319 // kCurr = f(t, y);
320  auto kCurr = f(t,y);
321 
322  return execute<TType, TStage + 1>(f, y, t, h, kCurr);
323  }
324 
342  template<typename TType, size_t TStage = 0, typename... TArgs>
343  typename std::enable_if<(TStage > 0) && (TStage < TStages) &&
344  (Conjunction<std::is_convertible<TArgs, TType>...>::value), TType>::type
345  inline execute(const std::function<TType(TScalarType, TType)> &f, const TType &y, const TScalarType t,
346  const TScalarType h,
347  TArgs... vars)
348  {
349  TType kCurr{}, yCurr = y;
350 
351  yCurr += h * sumCoefficients<TType, TStage>(vars...);
352  TScalarType tCurr = t + getNode<TStage - 1>() * h;
353 
354  kCurr = f(tCurr, yCurr);
355 
356  return execute<TType, TStage + 1>(f, y, t, h, vars..., kCurr);
357 // return execute<TType, TStage + 1>(f, y, t, h, vars..., f(t + getNode<TStage - 1>() * h, y + h * sumCoefficients<TType, TStage>(vars...)));
358  }
359 
377  template<typename TType, size_t TStage = 0, typename... TArgs>
378  typename std::enable_if<(TStage > 0) && (TStage == TStages) &&
379  (Conjunction<std::is_convertible<TArgs, TType>...>::value), TType>::type
380  inline execute(const std::function<TType(TScalarType, TType)> &f, const TType &y, const TScalarType t,
381  const TScalarType h,
382  TArgs... vars)
383  {
384  TType yFinal = y;
385  yFinal += h * sumWeights<TType>(vars...);
386 
387  return yFinal;
388  }
389 
390  protected:
400  template<int TTimes, typename TIndices = typename Indices<TTimes>::type>
401  struct Rows;
402 
409  template<int TTimes, int... TIndices>
410  struct Rows<TTimes, IndexSequence<TIndices...>> {
411  using type = std::tuple<ButcherRow<TIndices>...>;
412  };
413 
414  typename Rows<TStages>::type rows_;
415 
430  template<typename TType, size_t TStage = 0, typename... TArgs, typename TIndex = typename Indices<sizeof...(TArgs)>::type>
431  typename std::enable_if<
432  (TStage > 0) && (TStage < TStages) &&
433  (Conjunction<std::is_convertible<TArgs, TType>...>::value), TType>::type
434  inline sumCoefficients(TArgs... vars)
435  {
436  return sumCoefficientsImpl<TType, TStage - 1>(TIndex{}, vars...);
437  }
438 
453  template<typename TType, size_t TStage = 0, typename... TArgs, int... TIndices>
454  TType
456  {
457  return sum<TType>((getCoefficient<TStage, (size_t) TIndices>() * vars)...);
458  }
459 
471  template<typename TType, typename... TArgs, typename TIndex = typename Indices<sizeof...(TArgs)>::type>
472  typename std::enable_if<(Conjunction<std::is_convertible<TArgs, TType>...>::value), TType>::type
473  inline sumWeights(TArgs... vars)
474  {
475  return sumWeightsImpl<TType>(TIndex{}, vars...);
476  }
477 
491  template<typename TType, typename... TArgs, int... TIndices>
492  TType
494  {
495  return sum<TType>((getWeight<(size_t) TIndices>() * vars)...);
496  }
497 
510  template<typename TType, typename TArg>
511  typename std::enable_if<std::is_convertible<TArg, TType>::value, TType>::type
512  inline sum(TArg var)
513  {
514  return var;
515  }
516 
527  template<typename TType, typename TArg>
528  typename std::enable_if<!std::is_convertible<TArg, TType>::value, TType>::type
529  inline sum(TArg var)
530  {
531  static_assert(std::is_convertible<TArg, TType>::value, "Summed values must be convertible.");
532  return var;
533  }
534 
547  template<typename TType, typename TArg, typename... TArgs>
548  typename std::enable_if<std::is_convertible<TArg, TType>::value &&
549  (Conjunction<std::is_convertible<TArgs, TType>...>::value), TType>::type
550  inline sum(TArg var, TArgs... vars)
551  {
552  return var + sum<TType>(vars...);
553  }
554 
565  template<typename TType, typename TArg, typename... TArgs>
566  typename std::enable_if<!(std::is_convertible<TArg, TType>::value &&
567  (Conjunction<std::is_convertible<TArgs, TType>...>::value)), TType>::type
568  inline sum(TArg var, TArgs... vars)
569  {
570  static_assert(std::is_convertible<TArg, TType>::value &&
571  (Conjunction<std::is_convertible<TArgs, TType>...>::value),
572  "Summed values must be convertible.");
573 
574  return var + sum<TType>(vars...);
575  }
576  };
577 
578 // --------------- Extended Butcher tableau -----------------------
579 
595  template<size_t TStages, typename TScalarType>
596  class ButcherTableau<TStages, true, TScalarType> {
597  public:
598  static_assert(TStages > 1, "Butcher tableau must have more than one stage.");
599 
600  static const size_t stages = TStages;
601  using type = TScalarType;
602 
609  template<size_t TStage, bool TIsWeights = (TStage >= (TStages - 1))>
611 
620  template<size_t TRow, size_t TCoeffIdx>
621  typename std::enable_if<(TRow < (TStages - 1)) && (TCoeffIdx < (TRow + 1)), TScalarType>::type
622  inline getCoefficient()
623  {
624  return std::get<TCoeffIdx>(std::get<TRow>(rows_).rkCoefficients);
625  }
626 
635  template<size_t TRow, size_t TCoeffIdx>
636  typename std::enable_if<!((TRow < (TStages - 1)) && (TCoeffIdx < (TRow + 1))), TScalarType>::type
637  inline getCoefficient()
638  {
639  static_assert(TRow < (TStages - 1),
640  "Row index must be smaller than the number of stages minus one to access coefficients.");
641  static_assert(TCoeffIdx < (TRow + 1), "Coefficient index must be less than or equal to the row index.");
642 
643  return TScalarType{};
644  }
645 
654  template<size_t TRow, size_t TCoeffIdx>
655  typename std::enable_if<(TRow < (TStages - 1)) && (TCoeffIdx < (TRow + 1)), void>::type
656  inline setCoefficient(TScalarType coeff)
657  {
658  std::get<TCoeffIdx>(std::get<TRow>(rows_).rkCoefficients) = coeff;
659  }
660 
669  template<size_t TRow, size_t TCoeffIdx>
670  typename std::enable_if<!((TRow < (TStages - 1)) && (TCoeffIdx < (TRow + 1))), void>::type
671  inline setCoefficient(TScalarType coeff)
672  {
673  static_assert(TRow < (TStages - 1),
674  "Row index must be smaller than the number of stages minus one to access coefficients.");
675  static_assert(TCoeffIdx < (TRow + 1), "Coefficient index must be less than or equal to the row index.");
676  }
677 
685  template<size_t TRow>
686  typename std::enable_if<TRow < (TStages - 1), TScalarType>::type
687  inline getNode()
688  {
689  return std::get<TRow>(rows_).node;
690  }
691 
699  template<size_t TRow>
700  typename std::enable_if<TRow >= (TStages - 1), TScalarType>::type
701  inline getNode()
702  {
703  static_assert(TRow < (TStages - 1),
704  "Row index must be smaller than the number of stages minus one to access coefficients.");
705 
706  return TScalarType{};
707  }
708 
716  template<size_t TRow>
717  typename std::enable_if<TRow < (TStages - 1), void>::type
718  inline setNode(TScalarType node)
719  {
720  std::get<TRow>(rows_).node = node;
721  }
722 
730  template<size_t TRow>
731  typename std::enable_if<TRow >= (TStages - 1), void>::type
732  inline setNode(TScalarType node)
733  {
734  static_assert(TRow < (TStages - 1),
735  "Row index must be smaller than the number of stages minus one to access coefficients.");
736  }
737 
748  template<size_t TIdx, bool THigherOrder = true>
749  typename std::enable_if<(TIdx < TStages) && THigherOrder, TScalarType>::type
750  inline getWeight()
751  {
752  return std::get<TIdx>(std::get<TStages - 1>(rows_).weights);
753  }
754 
765  template<size_t TIdx, bool THigherOrder = true>
766  typename std::enable_if<(TIdx < TStages) && !THigherOrder, TScalarType>::type
767  inline getWeight()
768  {
769  return std::get<TIdx>(std::get<TStages>(rows_).weights);
770  }
771 
780  template<size_t TIdx, bool THigherOrder = true>
781  typename std::enable_if<TIdx >= TStages, TScalarType>::type
782  inline getWeight()
783  {
784  static_assert(TIdx < TStages, "Weight index must be smaller than the number of stages.");
785 
786  return TScalarType{};
787  }
788 
799  template<size_t TIdx, bool THigherOrder = true>
800  typename std::enable_if<(TIdx < TStages) && THigherOrder, void>::type
801  inline setWeight(TScalarType weight)
802  {
803  std::get<TIdx>(std::get<TStages - 1>(rows_).weights) = weight;
804  }
805 
816  template<size_t TIdx, bool THigherOrder = true>
817  typename std::enable_if<(TIdx < TStages) && !THigherOrder, void>::type
818  inline setWeight(TScalarType weight)
819  {
820  std::get<TIdx>(std::get<TStages>(rows_).weights) = weight;
821  }
822 
831  template<size_t TIdx, bool THigherOrder = true>
832  typename std::enable_if<TIdx >= TStages, void>::type
833  inline setWeight(TScalarType weight)
834  {
835  static_assert(TIdx < TStages, "Weight index must be smaller than the number of stages.");
836  }
837 
857  template<typename TType, size_t TStage = 0, bool TIsEmbedded = false>
858  typename std::enable_if<(TStage == 0) && !TIsEmbedded, TType>::type
859  inline execute(const std::function<TType(TScalarType, TType)> &f, const TType &y, const TScalarType t,
860  const TScalarType h)
861  {
862  TType kCurr = f(t, y);
863 
864  return execute<TType, TStage + 1, TIsEmbedded>(f, y, t, h, kCurr);
865  }
866 
886  template<typename TType, size_t TStage = 0, bool TIsEmbedded = false>
887  typename std::enable_if<(TStage == 0) && TIsEmbedded, std::pair<TType, TType>>::type
888  inline execute(const std::function<TType(TScalarType, TType)> &f, const TType &y, const TScalarType t,
889  const TScalarType h)
890  {
891  TType kCurr = f(t, y);
892 
893  return execute<TType, TStage + 1, TIsEmbedded>(f, y, t, h, kCurr);
894  }
895 
914  template<typename TType, size_t TStage = 0, bool TIsEmbedded = false, typename... TArgs>
915  typename std::enable_if<(TStage > 0) && (TStage < TStages) && !TIsEmbedded &&
916  (Conjunction<std::is_convertible<TArgs, TType>...>::value), TType>::type
917  inline execute(const std::function<TType(TScalarType, TType)> &f, const TType &y, const TScalarType t,
918  const TScalarType h,
919  TArgs... vars)
920  {
921  TType kCurr{}, yCurr = y;
922 
923  yCurr += h * sumCoefficients<TType, TStage>(vars...);
924  TScalarType tCurr = t + getNode<TStage - 1>() * h;
925 
926  kCurr = f(tCurr, yCurr);
927 
928  return execute<TType, TStage + 1, TIsEmbedded>(f, y, t, h, vars..., kCurr);
929  }
930 
949  template<typename TType, size_t TStage = 0, bool TIsEmbedded = false, typename... TArgs>
950  typename std::enable_if<(TStage > 0) && (TStage < TStages) && TIsEmbedded &&
951  (Conjunction<std::is_convertible<TArgs, TType>...>::value), std::pair<TType, TType>>::type
952  inline execute(const std::function<TType(TScalarType, TType)> &f, const TType &y, const TScalarType t,
953  const TScalarType h,
954  TArgs... vars)
955  {
956  TType kCurr{}, yCurr = y;
957 
958  yCurr += h * sumCoefficients<TType, TStage>(vars...);
959  TScalarType tCurr = t + getNode<TStage - 1>() * h;
960 
961  kCurr = f(tCurr, yCurr);
962 
963  return execute<TType, TStage + 1, TIsEmbedded>(f, y, t, h, vars..., kCurr);
964  }
965 
984  template<typename TType, size_t TStage = 0, bool TIsEmbedded = false, typename... TArgs>
985  typename std::enable_if<(TStage > 0) && (TStage == TStages) && !TIsEmbedded &&
986  (Conjunction<std::is_convertible<TArgs, TType>...>::value), TType>::type
987  inline execute(const std::function<TType(TScalarType, TType)> &f, const TType &y, const TScalarType t,
988  const TScalarType h,
989  TArgs... vars)
990  {
991  TType yFinal = y;
992  yFinal += h * sumWeights<TType>(vars...);
993 
994  return yFinal;
995  }
996 
1015  template<typename TType, size_t TStage = 0, bool TIsEmbedded = false, typename... TArgs>
1016  typename std::enable_if<(TStage > 0) && (TStage == TStages) && TIsEmbedded &&
1017  (Conjunction<std::is_convertible<TArgs, TType>...>::value), std::pair<TType, TType>>::type
1018  inline execute(const std::function<TType(TScalarType, TType)> &f, const TType &y, const TScalarType t,
1019  const TScalarType h,
1020  TArgs... vars)
1021  {
1022  TType yFinal = y;
1023  TType summedWeightsHigher = sumWeights<TType, true>(vars...);
1024  TType summedWeightsLower = sumWeights<TType, false>(vars...);
1025  yFinal += h * summedWeightsHigher;
1026 
1027  std::pair<TType, TType> pair{yFinal, h * (summedWeightsHigher - summedWeightsLower)};
1028 
1029  return pair;
1030  }
1031 
1032  protected:
1042  template<int TTimes, typename TIndices = typename Indices<TTimes>::type>
1043  struct Rows;
1044 
1051  template<int TTimes, int... TIndices>
1052  struct Rows<TTimes, IndexSequence<TIndices...>> {
1053  using type = std::tuple<ButcherRow<TIndices>...>;
1054  };
1055 
1056  typename Rows<TStages + 1>::type rows_;
1057 
1072  template<typename TType, size_t TStage = 0, typename... TArgs, typename TIndex = typename Indices<sizeof...(TArgs)>::type>
1073  typename std::enable_if<
1074  (TStage > 0) && (TStage < TStages) &&
1075  (Conjunction<std::is_convertible<TArgs, TType>...>::value), TType>::type
1076  inline sumCoefficients(TArgs... vars)
1077  {
1078  return sumCoefficientsImpl<TType, TStage - 1>(TIndex{}, vars...);
1079  }
1080 
1095  template<typename TType, size_t TStage = 0, typename... TArgs, int... TIndices>
1096  TType
1098  {
1099  return sum<TType>((getCoefficient<TStage, (size_t) TIndices>() * vars)...);
1100  }
1101 
1114  template<typename TType, bool THigherOrder = true, typename... TArgs, typename TIndex = typename Indices<sizeof...(TArgs)>::type>
1115  typename std::enable_if<(Conjunction<std::is_convertible<TArgs, TType>...>::value), TType>::type
1116  inline sumWeights(TArgs... vars)
1117  {
1118  return sumWeightsImpl<TType, THigherOrder>(TIndex{}, vars...);
1119  }
1120 
1135  template<typename TType, bool THigherOrder, typename... TArgs, int... TIndices>
1136  typename std::enable_if<THigherOrder, TType>::type
1138  {
1139  return sum<TType>((getWeight<(size_t) TIndices, true>() * vars)...);
1140  }
1141 
1156  template<typename TType, bool THigherOrder, typename... TArgs, int... TIndices>
1157  typename std::enable_if<!THigherOrder, TType>::type
1159  {
1160  return sum<TType>((getWeight<(size_t) TIndices, false>() * vars)...);
1161  }
1162 
1175  template<typename TType, typename TArg>
1176  typename std::enable_if<std::is_convertible<TArg, TType>::value, TType>::type
1177  inline sum(TArg var)
1178  {
1179  return var;
1180  }
1181 
1192  template<typename TType, typename TArg>
1193  typename std::enable_if<!std::is_convertible<TArg, TType>::value, TType>::type
1194  inline sum(TArg var)
1195  {
1196  static_assert(std::is_convertible<TArg, TType>::value, "Summed values must be convertible.");
1197  return var;
1198  }
1199 
1212  template<typename TType, typename TArg, typename... TArgs>
1213  typename std::enable_if<std::is_convertible<TArg, TType>::value &&
1214  (Conjunction<std::is_convertible<TArgs, TType>...>::value), TType>::type
1215  inline sum(TArg var, TArgs... vars)
1216  {
1217  return var + sum<TType>(vars...);
1218  }
1219 
1230  template<typename TType, typename TArg, typename... TArgs>
1231  typename std::enable_if<!(std::is_convertible<TArg, TType>::value &&
1232  (Conjunction<std::is_convertible<TArgs, TType>...>::value)), TType>::type
1233  inline sum(TArg var, TArgs... vars)
1234  {
1235  static_assert(std::is_convertible<TArg, TType>::value &&
1236  (Conjunction<std::is_convertible<TArgs, TType>...>::value),
1237  "Summed values must be convertible.");
1238 
1239  return var + sum<TType>(vars...);
1240  }
1241  };
1242  }
1243 }
1244 
1245 #endif //LODESTAR_BUTCHERTABLEAU_HPP
ls::primitives::ButcherTableau< TStages, true, TScalarType >::getWeight
std::enable_if<(TIdx< TStages) &&!THigherOrder, TScalarType >::type getWeight()
Returns the weight value at the given index.
Definition: ButcherTableau.hpp:767
ls::primitives::ButcherTableau< TStages, true, TScalarType >::sumWeights
std::enable_if<(Conjunction< std::is_convertible< TArgs, TType >... >::value), TType >::type sumWeights(TArgs... vars)
Sums weight coefficients multiplied with results from previous stages.
Definition: ButcherTableau.hpp:1116
Indices
Definition: Indices.hpp:15
ls::primitives::ButcherTableau< TStages, false, TScalarType >::sum
std::enable_if<!std::is_convertible< TArg, TType >::value, TType >::type sum(TArg var)
Degenerate case of helper function that sums arguments.
Definition: ButcherTableau.hpp:529
ls::primitives::ButcherTableau< TStages, true, TScalarType >::getWeight
std::enable_if<(TIdx< TStages) &&THigherOrder, TScalarType >::type getWeight()
Returns the weight value at the given index.
Definition: ButcherTableau.hpp:750
ls::primitives::ButcherTableau< TStages, false, TScalarType >::type
TScalarType type
Number of stages.
Definition: ButcherTableau.hpp:99
ls::primitives::ButcherTableau< TStages, false, TScalarType >::getCoefficient
std::enable_if<(TRow<(TStages - 1)) &&(TCoeffIdx<(TRow+1)), TScalarType >::type getCoefficient()
Returns the Runge-Kutta coefficient for the given row number and coefficient index.
Definition: ButcherTableau.hpp:120
ls::primitives::ButcherTableau< TStages, false, TScalarType >::setNode
std::enable_if< TRow<(TStages - 1), void >::type inline setNode(TScalarType node) { std::get< TRow >rows_).node=node;} template< size_t TRow > typename std::enable_if< TRow >=(TStages - 1), void >::type setNode(TScalarType node)
Sets the node value for the given row number.
Definition: ButcherTableau.hpp:230
ls::primitives::ButcherTableau< TStages, true, TScalarType >::sumWeightsImpl
std::enable_if<!THigherOrder, TType >::type sumWeightsImpl(IndexSequence< TIndices... >, TArgs... vars)
Implementation of sumWeights.
Definition: ButcherTableau.hpp:1158
ls::primitives::ButcherTableau< TStages, true, TScalarType >::execute
std::enable_if<(TStage > 0) &&(TStage==TStages) &&TIsEmbedded &&(Conjunction< std::is_convertible< TArgs, TType >... >::value), std::pair< TType, TType > >::type execute(const std::function< TType(TScalarType, TType)> &f, const TType &y, const TScalarType t, const TScalarType h, TArgs... vars)
Executes the explicit Runge-Kutta scheme defined in the Butcher tableau.
Definition: ButcherTableau.hpp:1018
ls::primitives::ButcherTableau
Implements a compile-time structure for simple/extended Butcher tableaus.
Definition: ButcherTableau.hpp:89
ls::primitives::ButcherTableau< TStages, true, TScalarType >::getCoefficient
std::enable_if<(TRow<(TStages - 1)) &&(TCoeffIdx<(TRow+1)), TScalarType >::type getCoefficient()
Returns the Runge-Kutta coefficient for the given row number and coefficient index.
Definition: ButcherTableau.hpp:622
ls::primitives::ButcherTableau< TStages, true, TScalarType >::getCoefficient
std::enable_if<!((TRow<(TStages - 1)) &&(TCoeffIdx<(TRow+1))), TScalarType >::type getCoefficient()
Degenerate Runge-Kutta coefficient getter.
Definition: ButcherTableau.hpp:637
Conjunction
Definition: Conjunction.hpp:12
ls::primitives::ButcherTableau< TStages, true, TScalarType >::sum
std::enable_if< std::is_convertible< TArg, TType >::value, TType >::type sum(TArg var)
Helper function that sums arguments.
Definition: ButcherTableau.hpp:1177
ls::primitives::ButcherTableau< TStages, false, TScalarType >::sumWeights
std::enable_if<(Conjunction< std::is_convertible< TArgs, TType >... >::value), TType >::type sumWeights(TArgs... vars)
Sums weight coefficients multiplied with results from previous stages.
Definition: ButcherTableau.hpp:473
ls::primitives::ButcherTableau< TStages, true, TScalarType >::setWeight
std::enable_if<(TIdx< TStages) &&!THigherOrder, void >::type setWeight(TScalarType weight)
Sets the weight value for the given index.
Definition: ButcherTableau.hpp:818
ls::primitives::ButcherTableau< TStages, true, TScalarType >::setWeight
std::enable_if<(TIdx< TStages) &&THigherOrder, void >::type setWeight(TScalarType weight)
Sets the weight value for the given index.
Definition: ButcherTableau.hpp:801
ls::primitives::ButcherTableau< TStages, true, TScalarType >::setNode
std::enable_if< TRow<(TStages - 1), void >::type inline setNode(TScalarType node) { std::get< TRow >rows_).node=node;} template< size_t TRow > typename std::enable_if< TRow >=(TStages - 1), void >::type setNode(TScalarType node)
Sets the node value for the given row number.
Definition: ButcherTableau.hpp:732
ls
Main Lodestar code.
Definition: BilinearTransformation.hpp:12
ls::primitives::ButcherTableau< TStages, false, TScalarType >::sum
std::enable_if<!(std::is_convertible< TArg, TType >::value &&(Conjunction< std::is_convertible< TArgs, TType >... >::value)), TType >::type sum(TArg var, TArgs... vars)
Degenerate case of helper function that sums arguments.
Definition: ButcherTableau.hpp:568
ls::primitives::ButcherTableau< TStages, true, TScalarType >::execute
std::enable_if<(TStage > 0) &&(TStage< TStages) &&!TIsEmbedded &&(Conjunction< std::is_convertible< TArgs, TType >... >::value), TType >::type execute(const std::function< TType(TScalarType, TType)> &f, const TType &y, const TScalarType t, const TScalarType h, TArgs... vars)
Executes the explicit Runge-Kutta scheme defined in the Butcher tableau.
Definition: ButcherTableau.hpp:917
ls::primitives::ButcherTableau< TStages, false, TScalarType >::setCoefficient
std::enable_if<!((TRow<(TStages - 1)) &&(TCoeffIdx<(TRow+1))), void >::type setCoefficient(TScalarType coeff)
Degenerate Runge-Kutta coefficient setter.
Definition: ButcherTableau.hpp:169
ls::primitives::ButcherTableau< TStages, true, TScalarType >::execute
std::enable_if<(TStage > 0) &&(TStage==TStages) &&!TIsEmbedded &&(Conjunction< std::is_convertible< TArgs, TType >... >::value), TType >::type execute(const std::function< TType(TScalarType, TType)> &f, const TType &y, const TScalarType t, const TScalarType h, TArgs... vars)
Executes the explicit Runge-Kutta scheme defined in the Butcher tableau.
Definition: ButcherTableau.hpp:987
ls::primitives::ButcherTableau< TStages, false, TScalarType >::sumCoefficientsImpl
TType sumCoefficientsImpl(IndexSequence< TIndices... >, TArgs... vars)
Implementation of sumCoefficients.
Definition: ButcherTableau.hpp:455
ls::primitives::ButcherTableau< TStages, false, TScalarType >::setCoefficient
std::enable_if<(TRow<(TStages - 1)) &&(TCoeffIdx<(TRow+1)), void >::type setCoefficient(TScalarType coeff)
Sets the Runge-Kutta coefficient at the given row number and coefficient index.
Definition: ButcherTableau.hpp:154
ls::primitives::ButcherTableau< TStages, true, TScalarType >::sumWeightsImpl
std::enable_if< THigherOrder, TType >::type sumWeightsImpl(IndexSequence< TIndices... >, TArgs... vars)
Implementation of sumWeights.
Definition: ButcherTableau.hpp:1137
ls::primitives::ButcherTableau< TStages, false, TScalarType >::setWeight
std::enable_if< TIdx< TStages, void >::type inline setWeight(TScalarType weight) { std::get< TIdx >std::get< TStages - 1 >rows_).weights)=weight;} template< size_t TIdx > typename std::enable_if< TIdx >=TStages, void >::type setWeight(TScalarType weight)
Sets the weight value for the given index.
Definition: ButcherTableau.hpp:289
ls::primitives::ButcherTableau< TStages, true, TScalarType >::getNode
std::enable_if< TRow<(TStages - 1), TScalarType >::type inline getNode() { return std::get< TRow >rows_).node;} template< size_t TRow > typename std::enable_if< TRow >=(TStages - 1), TScalarType >::type getNode()
Returns the node value for the given row number.
Definition: ButcherTableau.hpp:701
ls::primitives::detail::ButcherRowImpl
A Butcher tableau row.
Definition: ButcherTableau.hpp:35
ls::primitives::ButcherTableau< TStages, true, TScalarType >::sumCoefficients
std::enable_if<(TStage > 0) &&(TStage< TStages) &&(Conjunction< std::is_convertible< TArgs, TType >... >::value), TType >::type sumCoefficients(TArgs... vars)
Rows of the Butcher tableau.
Definition: ButcherTableau.hpp:1076
ls::primitives::ButcherTableau< TStages, false, TScalarType >::sum
std::enable_if< std::is_convertible< TArg, TType >::value, TType >::type sum(TArg var)
Helper function that sums arguments.
Definition: ButcherTableau.hpp:512
ls::primitives::ButcherTableau< TStages, true, TScalarType >::sum
std::enable_if< std::is_convertible< TArg, TType >::value &&(Conjunction< std::is_convertible< TArgs, TType >... >::value), TType >::type sum(TArg var, TArgs... vars)
Helper function that sums arguments.
Definition: ButcherTableau.hpp:1215
ls::primitives::ButcherTableau< TStages, false, TScalarType >::sumCoefficients
std::enable_if<(TStage > 0) &&(TStage< TStages) &&(Conjunction< std::is_convertible< TArgs, TType >... >::value), TType >::type sumCoefficients(TArgs... vars)
Rows of the Butcher tableau.
Definition: ButcherTableau.hpp:434
ls::primitives::ButcherTableau< TStages, true, TScalarType >::sum
std::enable_if<!std::is_convertible< TArg, TType >::value, TType >::type sum(TArg var)
Degenerate case of helper function that sums arguments.
Definition: ButcherTableau.hpp:1194
ls::primitives::ButcherTableau< TStages, false, TScalarType >::execute
std::enable_if<(TStage > 0) &&(TStage==TStages) &&(Conjunction< std::is_convertible< TArgs, TType >... >::value), TType >::type execute(const std::function< TType(TScalarType, TType)> &f, const TType &y, const TScalarType t, const TScalarType h, TArgs... vars)
Executes the explicit Runge-Kutta scheme defined in the Butcher tableau.
Definition: ButcherTableau.hpp:380
ls::primitives::ButcherTableau< TStages, false, TScalarType >::getCoefficient
std::enable_if<!((TRow<(TStages - 1)) &&(TCoeffIdx<(TRow+1))), TScalarType >::type getCoefficient()
Degenerate Runge-Kutta coefficient getter.
Definition: ButcherTableau.hpp:135
ls::primitives::ButcherTableau< TStages, false, TScalarType >::execute
std::enable_if< TStage==0, TType >::type execute(const std::function< TType(TScalarType, TType)> &f, const TType &y, const TScalarType t, const TScalarType h)
Executes the explicit Runge-Kutta scheme defined in the Butcher tableau.
Definition: ButcherTableau.hpp:314
ls::primitives::ButcherTableau< TStages, true, TScalarType >::execute
std::enable_if<(TStage==0) &&!TIsEmbedded, TType >::type execute(const std::function< TType(TScalarType, TType)> &f, const TType &y, const TScalarType t, const TScalarType h)
Executes the explicit Runge-Kutta scheme defined in the Butcher tableau.
Definition: ButcherTableau.hpp:859
ls::primitives::ButcherTableau< TStages, true, TScalarType >::execute
std::enable_if<(TStage > 0) &&(TStage< TStages) &&TIsEmbedded &&(Conjunction< std::is_convertible< TArgs, TType >... >::value), std::pair< TType, TType > >::type execute(const std::function< TType(TScalarType, TType)> &f, const TType &y, const TScalarType t, const TScalarType h, TArgs... vars)
Executes the explicit Runge-Kutta scheme defined in the Butcher tableau.
Definition: ButcherTableau.hpp:952
ls::primitives::ButcherTableau< TStages, false, TScalarType >::getWeight
std::enable_if< TIdx< TStages, TScalarType >::type inline getWeight() { return std::get< TIdx >std::get< TStages - 1 >rows_).weights);} template< size_t TIdx > typename std::enable_if< TIdx >=TStages, TScalarType >::type getWeight()
Returns the weight value at the given index.
Definition: ButcherTableau.hpp:259
ls::primitives::ButcherTableau< TStages, false, TScalarType >::sumWeightsImpl
TType sumWeightsImpl(IndexSequence< TIndices... >, TArgs... vars)
Implementation of sumWeights.
Definition: ButcherTableau.hpp:493
ls::primitives::ButcherTableau< TStages, false, TScalarType >::execute
std::enable_if<(TStage > 0) &&(TStage< TStages) &&(Conjunction< std::is_convertible< TArgs, TType >... >::value), TType >::type execute(const std::function< TType(TScalarType, TType)> &f, const TType &y, const TScalarType t, const TScalarType h, TArgs... vars)
Executes the explicit Runge-Kutta scheme defined in the Butcher tableau.
Definition: ButcherTableau.hpp:345
ls::primitives::ButcherTableau< TStages, true, TScalarType >::setCoefficient
std::enable_if<(TRow<(TStages - 1)) &&(TCoeffIdx<(TRow+1)), void >::type setCoefficient(TScalarType coeff)
Sets the Runge-Kutta coefficient at the given row number and coefficient index.
Definition: ButcherTableau.hpp:656
ls::primitives::ButcherTableau< TStages, true, TScalarType >::execute
std::enable_if<(TStage==0) &&TIsEmbedded, std::pair< TType, TType > >::type execute(const std::function< TType(TScalarType, TType)> &f, const TType &y, const TScalarType t, const TScalarType h)
Executes the explicit Runge-Kutta scheme defined in the Butcher tableau.
Definition: ButcherTableau.hpp:888
ls::primitives::ButcherTableau< TStages, false, TScalarType >::getNode
std::enable_if< TRow<(TStages - 1), TScalarType >::type inline getNode() { return std::get< TRow >rows_).node;} template< size_t TRow > typename std::enable_if< TRow >=(TStages - 1), TScalarType >::type getNode()
Returns the node value for the given row number.
Definition: ButcherTableau.hpp:199
ls::primitives::ButcherTableau< TStages, false, TScalarType >::sum
std::enable_if< std::is_convertible< TArg, TType >::value &&(Conjunction< std::is_convertible< TArgs, TType >... >::value), TType >::type sum(TArg var, TArgs... vars)
Helper function that sums arguments.
Definition: ButcherTableau.hpp:550
ls::primitives::ButcherTableau< TStages, true, TScalarType >::sum
std::enable_if<!(std::is_convertible< TArg, TType >::value &&(Conjunction< std::is_convertible< TArgs, TType >... >::value)), TType >::type sum(TArg var, TArgs... vars)
Degenerate case of helper function that sums arguments.
Definition: ButcherTableau.hpp:1233
ls::primitives::ButcherTableau< TStages, true, TScalarType >::setCoefficient
std::enable_if<!((TRow<(TStages - 1)) &&(TCoeffIdx<(TRow+1))), void >::type setCoefficient(TScalarType coeff)
Degenerate Runge-Kutta coefficient setter.
Definition: ButcherTableau.hpp:671
ls::primitives::ButcherTableau< TStages, true, TScalarType >::sumCoefficientsImpl
TType sumCoefficientsImpl(IndexSequence< TIndices... >, TArgs... vars)
Implementation of sumCoefficients.
Definition: ButcherTableau.hpp:1097
ls::primitives::ButcherTableau< TStages, true, TScalarType >::type
TScalarType type
Number of stages.
Definition: ButcherTableau.hpp:601
IndexSequence
Definition: Indices.hpp:11