File: gridKernelRegression.h

package info (click to toggle)
stopt 6.1%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 9,264 kB
  • sloc: cpp: 75,778; python: 6,012; makefile: 72; sh: 57
file content (163 lines) | stat: -rw-r--r-- 9,442 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
// Copyright (C) 2017 EDF
// All Rights Reserved
// This code is published under the GNU Lesser General Public License (GNU LGPL)
#ifndef GRIDKERNELREGRESSION_H
#define GRIDKERNELREGRESSION_H
#include <vector>
#include <memory>
#include <Eigen/Dense>
#include <boost/multi_array.hpp>
#include "StOpt/core/grids/InterpolatorSpectral.h"



/// \file gridKernelConstruction.h
namespace StOpt
{


/// \brief  calculate the regressed values  on the grid p_z
/// \param  p_newX    renormalized particle
/// \param  p_y       second member to regress
/// \param  p_z       grid points for regressions (in each dimension)
/// \param  p_h       bandwith for each point in p_z (in each dimension)
/// \param  p_g       storing all point in p_z + /- p_h (in each dimension)
/// \param  p_xG      permits to affect each particle to a slice defined by p_g in each direction
/// \param  p_zl      for each dimension, indicate in which point in p_g  correspond a  point  in  p_z-p_h
/// \param  p_zr      for each dimension, indicate in which point in p_g   correspond a  point  in  p_z+p_h
/// \param  p_tick    tick to help regresssion matrix
/// \param  p_bLinear true if we use linear regressions
/// \return An array with all the regressed value at p_z points
Eigen::ArrayXd  gridKernelRegressedValuesOnGrid(const Eigen::ArrayXXd   &p_newX,
        const Eigen::ArrayXd &p_y,
        const std::vector< std::shared_ptr<Eigen::ArrayXd> >   &p_z,
        const std::vector< std::shared_ptr<Eigen::ArrayXd> >   &p_h,
        const std::vector<  std::shared_ptr<Eigen::ArrayXd> > &p_g,
        const Eigen::ArrayXXi &p_xG,
        const std::vector< std::shared_ptr<Eigen::ArrayXi> > &p_zl,
        const std::vector< std::shared_ptr<Eigen::ArrayXi> > &p_zr,
        const bool &p_bLinear,
        double p_tick = 1e-7)  ;

/// \brief  calculate he regressed values  on the grid p_z
/// \param  p_newX    renormalized particle
/// \param  p_y       second member to regress
/// \param  p_z       grid points for regressions (in each dimension)
/// \param  p_h       bandwith for each point in p_z (in each dimension)
/// \param  p_g       storing all point in p_z + /- p_h (in each dimension)
/// \param  p_xG      permits to affect each particle to a slice defined by p_g in each direction
/// \param  p_zl      for each dimension, indicate in which point in p_g  correspond a  point  in  p_z-p_h
/// \param  p_zr      for each dimension, indicate in which point in p_g   correspond a  point  in  p_z+p_h
/// \param  p_bLinear true if we use linear regressions
/// \param  p_tick    tick to help regresssion matrix
/// \return An array with all the regressed value at p_z points
Eigen::ArrayXd  gridKernelRegressedValuesOnGridStable(const Eigen::ArrayXXd   &p_newX,
        const Eigen::ArrayXd &p_y,
        const std::vector< std::shared_ptr<Eigen::ArrayXd> >   &p_z,
        const std::vector< std::shared_ptr<Eigen::ArrayXd> >   &p_h,
        const std::vector<  std::shared_ptr<Eigen::ArrayXd> > &p_g,
        const Eigen::ArrayXXi &p_xG,
        const std::vector< std::shared_ptr<Eigen::ArrayXi> > &p_zl,
        const std::vector< std::shared_ptr<Eigen::ArrayXi> > &p_zr,
        const bool &p_bLinear,
        double p_tick = 1e-7)    ;

/// \brief  calculate the regressed values  associated to (p_newX,p_y)
/// \param  p_newX    renormalized particle
/// \param  p_newY    renormalized second member to regress
/// \param  p_z       grid points for regressions (in each dimension)
/// \param  p_h       bandwith for each point in p_z (in each dimension)
/// \param  p_g       storing all point in p_z + /- p_h (in each dimension)
/// \param  p_xG      permits to affect each particle to a slice defined by p_g in each direction
/// \param  p_zl      for each dimension, indicate in which point in p_g  correspond a  point  in  p_z-p_h
/// \param  p_zr      for each dimension, indicate in which point in p_g   correspond a  point  in  p_z+p_h
/// \param  p_iSort   in each direction sorted particles by their number
/// \param  p_tick    tick to help regresssion matrix
/// \param  p_bLinear true if linear regression
/// \return An array with all the regressed value at p_z points
Eigen::ArrayXd  gridKernelRegressedValues(const Eigen::ArrayXXd   &p_newX,
        const Eigen::ArrayXd &p_y,
        const std::vector< std::shared_ptr<Eigen::ArrayXd> >   &p_z,
        const std::vector< std::shared_ptr<Eigen::ArrayXd> >   &p_h,
        const std::vector<  std::shared_ptr<Eigen::ArrayXd> > &p_g,
        const Eigen::ArrayXXi &p_xG,
        const std::vector< std::shared_ptr<Eigen::ArrayXi> > &p_zl,
        const std::vector< std::shared_ptr<Eigen::ArrayXi> > &p_zr,
        const Eigen::ArrayXXi &p_iSort,
        const bool &p_bLinear,
        double p_tick = 1e-7)    ;


/// \brief General Subroutine : naive implementation for all points
/// \param p_X        sample point
/// \param p_Y        sample point to regress with respect to p_X
/// \param p_prop     Proportion of point taken by each local regression (bandwith)
/// \param p_bLinear do we use linera regressions
/// \param p_q        p_q by number of simulations correspond to the number of points in the grid
Eigen::ArrayXd  locAdapRegNaive(const Eigen::ArrayXXd &p_X, const Eigen::ArrayXd &p_Y, const double &p_prop, const bool &p_bLinear, double p_q = 1., double p_tick = 1e-7);

/// \brief General Subroutine : fast sum implementation for all points
/// \param p_X        sample point
/// \param p_Y        sample point to regress with respect to p_X
/// \param p_prop     Proportion of point taken by each local regression (bandwith)
/// \param p_bLinear do we use linera regressions
/// \param p_q        p_q by number of simulations correspond to the number of points in the grid
Eigen::ArrayXd  locAdapReg(const Eigen::ArrayXXd &p_X, const Eigen::ArrayXd &p_Y, const double &p_prop, const bool &p_bLinear, double p_q = 1., double p_tick = 1e-7);


/// \brief General Subroutine : naive implementation for all points
/// \param p_X        sample point
/// \param p_Y        sample point to regress with respect to p_X
/// \param p_prop     Proportion of point taken by each local regression (bandwith)
/// \param p_bLinear  do we use linera regressions
/// \param p_q        p_q by number of simulations correspond to the number of points in the grid
Eigen::ArrayXd  locAdapRegNaiveOnGrid(const Eigen::ArrayXXd &p_X, const Eigen::ArrayXd &p_Y, const double &p_prop, const bool &p_bLinear, double p_q = 1., double p_tick = 1e-7);

/// \brief General Subroutine : naive implementation for all points (stabilized version)
/// \param p_X        sample point
/// \param p_Y        sample point to regress with respect to p_X
/// \param p_prop     Proportion of point taken by each local regression (bandwith)
/// \param p_bLinear  do we use linera regressions
/// \param p_q        p_q by number of simulations correspond to the number of points in the grid
Eigen::ArrayXd  locAdapRegNaiveOnGridStab(const Eigen::ArrayXXd &p_X, const Eigen::ArrayXd &p_Y, const double &p_prop, const bool &p_bLinear, double p_q = 1., double p_tick = 1e-7);

/// \brief  Calculate for all points used by the kernel method the regressed values from the interplated values on the grid
/// \param  p_newX    renormalized particle
/// \param  p_yHatZ   regressed values on the grid  defined by p_z plus scaling (mean, std)
/// \param  p_z       grid points for regressions (in each dimension)
/// \param  p_iSort   in each direction sorted particles by their number
/// \return An array with all the regressed value for the particle p_newX
Eigen::ArrayXd  fromGridValuesGetRegAllSim(const Eigen::ArrayXXd   &p_newX,
        const Eigen::ArrayXd   &p_yHatZ,
        const std::vector< std::shared_ptr<Eigen::ArrayXd> >   &p_z,
        const Eigen::ArrayXXi &p_iSort);



/// \brief  Calculate for one point  used by the kernel method the regressed value from the interpolated values on the grid
///         by a naive method
/// \param  p_newX    renormalized particle
/// \param  p_yHatZ   regressed values on the grid  defined by p_z plus scaling (mean, std)
/// \param  p_z       grid points for regressions (in each dimension)
/// \return  the interpolated value
double   fromGridValuesGetRegASim(const Eigen::ArrayXd    &p_newX,
                                  const Eigen::ArrayXd   &p_yHatZ,
                                  const std::vector< std::shared_ptr<Eigen::ArrayXd> >   &p_z);


/// \brief  Calculate for one point  used by the kernel method the regressed value from the interpolated values on the grid
///         Used to avoid costly interpolation on all function basis
/// \param  p_newX    renormalized particle
/// \param  p_yHatZ   regressed values on the grid  defined by p_z plus scaling (mean, std)
/// \param  p_z       grid points for regressions (in each dimension)
/// \param  p_ptOfStock            grid point
/// \param  p_interpFuncBasis      spectral interpolator to interpolate on the basis the function value (regressed functions on all grid with
///                                spectral representation)
/// \return  the interpolated value
double   fromGridValuesGetRegASimOnBasis(const Eigen::ArrayXd    &p_newX,
        const std::vector< std::shared_ptr<Eigen::ArrayXd> >   &p_z,
        const Eigen::ArrayXd &p_ptOfStock,
        const std::vector< std::shared_ptr<InterpolatorSpectral> > &p_interpFuncBasis);
}

#endif /* GRIDKERNELREGRESSION_H */