MUQ  0.4.3
MultiIndexFactory.h
Go to the documentation of this file.
1 #ifndef MULTIINDEXFACTORY_H_
2 #define MULTIINDEXFACTORY_H_
3 
4 #include "MUQ/config.h"
5 
8 
9 namespace muq {
10  namespace Utilities{
11 
19 
20  public:
21 
30  static std::shared_ptr<MultiIndexSet> CreateTotalOrder(unsigned int const length,
31  unsigned int const maxOrder,
32  unsigned int const minOrder = 0,
33  std::shared_ptr<MultiIndexLimiter> limiter = std::make_shared<NoLimiter>());
34 
35 
36  static std::vector<std::shared_ptr<MultiIndexSet>> CreateTriTotalOrder(unsigned int const length,
37  unsigned int const maxOrder,
38  unsigned int const minOrder = 0,
39  std::shared_ptr<MultiIndexLimiter> limiter = std::make_shared<NoLimiter>());
40 
48  static std::shared_ptr<MultiIndexSet> CreateHyperbolic(unsigned int const length,
49  unsigned int const maxOrder,
50  const double q = 0.5,
51  std::shared_ptr<MultiIndexLimiter> limiter = std::make_shared<NoLimiter>());
52 
53 
54  static std::vector<std::shared_ptr<MultiIndexSet>> CreateTriHyperbolic(unsigned int const length,
55  unsigned int const maxOrder,
56  const double q = 0.5,
57  std::shared_ptr<MultiIndexLimiter> limiter = std::make_shared<NoLimiter>());
58 
65  static std::shared_ptr<MultiIndexSet> CreateFullTensor(unsigned int const length,
66  unsigned int const order,
67  std::shared_ptr<MultiIndexLimiter> limiter = std::make_shared<NoLimiter>());
68 
75  static std::shared_ptr<MultiIndexSet> CreateFullTensor(const Eigen::RowVectorXi& orders,
76  std::shared_ptr<MultiIndexLimiter> limiter = std::make_shared<NoLimiter>());
77 
87  static std::shared_ptr<MultiIndexSet> CreateAnisotropic(const Eigen::RowVectorXf& weights, const double epsilon);
88 
89 
99  static std::shared_ptr<MultiIndex> CreateSingleTerm(int totalDim, int nonzeroDim, int order);
100 
101 
102  private:
103 
104  static void RecursiveHyperbolicFill(const double maxOrderPow,
105  std::shared_ptr<MultiIndexSet> output,
106  unsigned int const currDim,
107  Eigen::RowVectorXi &base,
108  const double q,
109  std::shared_ptr<MultiIndexLimiter> limiter);
110 
111 
112  static void RecursiveTotalOrderFill(unsigned int const maxOrder,
113  unsigned int const minOrder,
114  std::shared_ptr<MultiIndexSet> output,
115  unsigned int const currDim,
116  Eigen::RowVectorXi &base,
117  std::shared_ptr<MultiIndexLimiter> limiter);
118 
119 
120  static void RecursiveTensor(const Eigen::RowVectorXi& orders,
121  std::shared_ptr<MultiIndexSet> output,
122  unsigned int const currDim,
123  Eigen::RowVectorXi &base,
124  std::shared_ptr<MultiIndexLimiter> limiter,
125  bool allInactive);
126 
127  };// class MultiIndexFactory
128  }// namespace Utilities
129 }// namespace muq
130 
131 
132 #endif
A factory class with static methods for generating MultiIndexSets.
static std::shared_ptr< MultiIndexSet > CreateTotalOrder(unsigned int const length, unsigned int const maxOrder, unsigned int const minOrder=0, std::shared_ptr< MultiIndexLimiter > limiter=std::make_shared< NoLimiter >())
Construct a total order limited MultiIndex.
static std::vector< std::shared_ptr< MultiIndexSet > > CreateTriHyperbolic(unsigned int const length, unsigned int const maxOrder, const double q=0.5, std::shared_ptr< MultiIndexLimiter > limiter=std::make_shared< NoLimiter >())
static void RecursiveTensor(const Eigen::RowVectorXi &orders, std::shared_ptr< MultiIndexSet > output, unsigned int const currDim, Eigen::RowVectorXi &base, std::shared_ptr< MultiIndexLimiter > limiter, bool allInactive)
static std::shared_ptr< MultiIndexSet > CreateFullTensor(unsigned int const length, unsigned int const order, std::shared_ptr< MultiIndexLimiter > limiter=std::make_shared< NoLimiter >())
Construct a full tensor product multiindex set.
static std::shared_ptr< MultiIndexSet > CreateHyperbolic(unsigned int const length, unsigned int const maxOrder, const double q=0.5, std::shared_ptr< MultiIndexLimiter > limiter=std::make_shared< NoLimiter >())
Construct a general hyperbolic MultiIndex.
static std::shared_ptr< MultiIndex > CreateSingleTerm(int totalDim, int nonzeroDim, int order)
Creates a single multiindex with one nonzero term.
static std::vector< std::shared_ptr< MultiIndexSet > > CreateTriTotalOrder(unsigned int const length, unsigned int const maxOrder, unsigned int const minOrder=0, std::shared_ptr< MultiIndexLimiter > limiter=std::make_shared< NoLimiter >())
static void RecursiveTotalOrderFill(unsigned int const maxOrder, unsigned int const minOrder, std::shared_ptr< MultiIndexSet > output, unsigned int const currDim, Eigen::RowVectorXi &base, std::shared_ptr< MultiIndexLimiter > limiter)
static void RecursiveHyperbolicFill(const double maxOrderPow, std::shared_ptr< MultiIndexSet > output, unsigned int const currDim, Eigen::RowVectorXi &base, const double q, std::shared_ptr< MultiIndexLimiter > limiter)
static std::shared_ptr< MultiIndexSet > CreateAnisotropic(const Eigen::RowVectorXf &weights, const double epsilon)
Construct an anisotropic multiindex set based on a priori information on the importance of each dimen...