File: LocalGridKernelRegression.h

package info (click to toggle)
stopt 5.12%2Bdfsg-3
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 8,860 kB
  • sloc: cpp: 70,456; python: 5,950; makefile: 72; sh: 57
file content (195 lines) | stat: -rw-r--r-- 9,419 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
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
// Copyright (C) 2017 EDF
// All Rights Reserved
// This code is published under the GNU Lesser General Public License (GNU LGPL)
#ifndef LOCALGRIDKERNELREGRESSION_H
#define LOCALGRIDKERNELREGRESSION_H
#include <vector>
#include <memory>
#include <Eigen/Dense>
#include <Eigen/SVD>
#include "StOpt/regression/BaseRegression.h"
#include <iostream>


/** \file LocalGridKernelRegression.h
 * \brief Computation conditional expectation using additive kernel
 *  - in  1D  using real Epanechnikov kernel
 *  - in nD using the grid kernel method
 *  .
 * When the grid is created, the number number of grid points in  a direction id is proportional
 * to the SVD value  S(id) in this  direction obtained with
 * the simulation particles renormalized to a Unit Centered Gaussian.
 * The grid is chosen such than  it constains roughly m_coefNbGridPoint nbSimul  points
 * In each direction the grid is adapted such that in each slice each mesh contains roughly the
 * same number of particles.
 * m_coeffBandWith permits to adapt the bandwidth.
 * \author  Xavier Warin
 */
namespace StOpt
{
class LocalGridKernelRegression : public BaseRegression
{
private:

    bool m_bZeroDate ;                    ///< Is the regression date zero ?
    Eigen::ArrayXXi m_iSort ; ///<  Permits to get for a particle sorted in a direction its position in the initial particles array
    std::vector< std::shared_ptr<Eigen::ArrayXd> > m_z ; ///< In each direction, define the grid points coordinate
    std::vector< std::shared_ptr<Eigen::ArrayXd> > m_h ;  ///< In each direction, for each point define the bandwidth
    std::vector< std::shared_ptr<Eigen::ArrayXd> > m_g ; ///< in each direction store the union of \f$ m_z-m_h\f$  and \f$ m_z+m_h \$ : the grid obtained is modified sugc that each space contains some particles
    std::vector< std::shared_ptr<Eigen::ArrayXi> > m_zl ; ///< permits to locat the position of \f$ m_z-m_h\f$ in \f$ m_g \f$.
    std::vector< std::shared_ptr<Eigen::ArrayXi> > m_zr ; ///< permits to locat the position of \f$ m_z+m_h\f$ in \f$ m_g \f$.
    Eigen::ArrayXXi m_xG ; ///< In each dimension, affect a particle to its slide in the grid defined by $m_g$
    double m_coeffBandWidth ; ///<Multiplicative coefficient to define bandwidth
    double   m_coefNbGridPoint ; ///< Multiplicative coefficient for the number of grid points
    bool  m_bLinear ; ///< true if use linear regression

    /// \brief Precompute grid points, bandwidth
    void  createGrid();


public :

    ///\brief Default constructor
    LocalGridKernelRegression() {}


    /// \brief Constructor for grid kernel regression
    /// \param  p_bZeroDate          first date is 0?
    /// \param  p_particles          particles used for the meshes.
    ///                              First dimension  : dimension of the problem,
    ///                              second dimension : the  number of particles
    /// \param p_coeffBandWidth      Multiplicative coefficient to define bandwidth
    /// \param p_coefNbGridPoint     Multiplicative coefficient for the number of grid points
    /// \param p_bLinear             True if linear regression
    LocalGridKernelRegression(const bool &p_bZeroDate,
                              const Eigen::ArrayXXd  &p_particles,
                              const double &p_coeffBandWidth,
                              const double   &p_coefNbGridPoint,
                              const bool &p_bLinear);


    ///\brief second  constructor
    /// \param p_coeffBandWidth      Multiplicative coefficient to define bandwidth
    /// \param p_coefNbGridPoint     Multiplicative coefficient for the number of grid points
    /// \param p_bLinear             True if linear regression
    LocalGridKernelRegression(const double &p_coeffBandWidth,
                              const double   &p_coefNbGridPoint,
                              const bool &p_bLinear): BaseRegression(true), m_coeffBandWidth(p_coeffBandWidth), m_coefNbGridPoint(p_coefNbGridPoint), m_bLinear(p_bLinear) {}



    /// \brief Last  constructor only used for out of sample simulations
    LocalGridKernelRegression(const bool &p_bZeroDate,
                              const Eigen::ArrayXd &p_meanX,
                              const Eigen::ArrayXd &p_etypX,
                              const Eigen::MatrixXd &p_svdMatrix,
                              const std::vector< std::shared_ptr<Eigen::ArrayXd> > &p_z,
                              const bool &p_bLinear): BaseRegression(p_bZeroDate, p_meanX, p_etypX, p_svdMatrix, true), m_z(p_z), m_bLinear(p_bLinear) {};


    /// \brief update the particles used in regression  and construct the matrices
    /// \param  p_bZeroDate    first date is 0?
    /// \param  p_particles    particles used for the meshes.
    ///                        First dimension  : dimension of the problem,
    ///                        second dimension : the  number of particles
    void updateSimulations(const bool &p_bZeroDate, const Eigen::ArrayXXd &p_particles);

    /// \brief  For this kernel method get back the regressed values on a deterministic grid with coordinates given by m_z
    /// \param  p_fToRegress  function to regress associated to each simulation used in optimization
    /// \return regressed values on the grid
    /// @{
    Eigen::ArrayXd getCoordBasisFunction(const Eigen::ArrayXd &p_fToRegress) const;
    Eigen::ArrayXd getCoordBasisFunctionStable(const Eigen::ArrayXd &p_fToRegress) const;
    ///@}
    /// \brief  For this kernel method get back the regressed values on a deterministic grid with coordinates given by m_z
    /// \param  p_fToRegress  function to regress associated to each simulation used in optimization (size : number of functions to regress \times the number of Monte Carlo simulations)
    /// \return regressed values on the grid  (size :  number of function to regress  \times number of grids points  )
    /// @{
    Eigen::ArrayXXd getCoordBasisFunctionMultiple(const Eigen::ArrayXXd &p_fToRegress) const ;
    ///@}

    /// \brief conditional expectation calculation
    /// \param  p_fToRegress  simulations  to regress used in optimization
    /// \return regressed value function
    /// @{
    Eigen::ArrayXd getAllSimulations(const Eigen::ArrayXd &p_fToRegress) const ;
    Eigen::ArrayXXd getAllSimulationsMultiple(const Eigen::ArrayXXd &p_fToRegress) const;
    ///@}


    /// \brief Use basis functions to reconstruct the solution
    /// \param p_basisCoefficients basis coefficients
    ///@{
    Eigen::ArrayXd reconstruction(const Eigen::ArrayXd   &p_basisCoefficients) const;
    Eigen::ArrayXXd reconstructionMultiple(const Eigen::ArrayXXd   &p_basisCoefficients) const;
    /// @}
    /// \brief use basis function to reconstruct a given simulation
    /// \param p_isim               simulation number
    /// \param p_basisCoefficients  basis coefficients to reconstruct a given conditional expectation
    double reconstructionASim(const int &p_isim, const Eigen::ArrayXd   &p_basisCoefficients) const ;

    /// \brief conditional expectation reconstruction
    /// \param  p_coordinates        coordinates to interpolate (uncertainty sample)
    /// \param  p_coordBasisFunction regression coordinates on the basis  (size: number of meshes multiplied by the dimension plus one)
    /// \return regressed value function reconstructed for each simulation
    double getValue(const Eigen::ArrayXd   &p_coordinates,
                    const Eigen::ArrayXd   &p_coordBasisFunction) const ;


    /// \brief permits to reconstruct a function with basis functions coefficients values given on a grid
    /// \param  p_coordinates          coordinates  (uncertainty sample)
    /// \param  p_ptOfStock            grid point
    /// \param  p_interpFuncBasis      spectral interpolator to interpolate the basis functions  coefficients used in regression on the grid (given for each basis function)
    double getAValue(const Eigen::ArrayXd &p_coordinates,  const Eigen::ArrayXd &p_ptOfStock,
                     const std::vector< std::shared_ptr<InterpolatorSpectral> > &p_interpFuncBasis) const;


    /// \brief get the number of basis functions
    inline int getNumberOfFunction() const
    {
        if (m_bZeroDate)
            return 1;
        else
        {
            int nPtGrid = m_z.size() + 1;
            for (int id = 0 ; id < m_z.size(); ++id)
                nPtGrid *= m_z[id]->size();
            return nPtGrid;
        }
    }

    /// \brief Clone the regressor
    virtual std::shared_ptr<BaseRegression> clone() const
    {
        return std::static_pointer_cast<BaseRegression>(std::make_shared<LocalGridKernelRegression>(*this));
    }

    /// \breif get back meshes
    inline std::vector< std::shared_ptr<Eigen::ArrayXd> > getZ() const
    {
        return m_z;
    }
    /// \brief get back scale factor for particles
    inline Eigen::ArrayXd getMeanX() const
    {
        return m_meanX;
    }
    /// \brief get back scale factor for particles
    inline Eigen::ArrayXd getEtypX() const
    {
        return m_etypX;
    }
    /// \brief get back rotation matrix for particles
    inline Eigen::MatrixXd getSvdMatrix() const
    {
        return m_svdMatrix;
    }
    /// \brief get back the king of regression (true if linear)
    inline bool getBLinear() const
    {
        return  m_bLinear;
    }
};
}

#endif /*   LOCALGRIDKERNELREGRESSION_H  */