File: localLinearDiscrLastDimMatrixOperation.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 (110 lines) | stat: -rw-r--r-- 6,877 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
// Copyright (C) 2021 EDF
// All Rights Reserved
// This code is published under the GNU Lesser General Public License (GNU LGPL)
#ifndef LOCALLINEARDISCRLASTDIMMATRIXOPERATION_H
#define LOCALLINEARDISCRLASTDIMMATRIXOPERATION_H
/** \file localLinearDiscrLastDimMatrixOperation.h
 *  \brief  Develop the normal operators for regression, permits inversion by Choleski of \f$ A^tA  x = A^t b \f$  and reconstruct the solution
 *          Here approximation is linear for first dimensions and last dimension has discrete values gathered in constant approximation
 *  \author Xavier Warin
 */
#include <vector>
#include <array>
#include <memory>
#include <Eigen/Dense>
#include "StOpt/core/grids/InterpolatorSpectral.h"

namespace StOpt
{
/**
 * \addtogroup LocalLinearDiscrLastDim
 * @{
 */
/// \brief Regression matrix calculation for local linear functions (except last dimension)
/// \param p_particles   Particles used for regression
/// \param p_simToCell   To each particle associates the cell it belongs to.
/// \param p_mesh        Description of the mesh coordinates for each mesh
/// \return             Regression matrix to calculate
Eigen::ArrayXXd localLinearDiscrLastDimMatrixCalculation(const Eigen::ArrayXXd &p_particles,
        const Eigen::ArrayXi &p_simToCell,
        const Eigen::Array< std::array< double, 2>, Eigen::Dynamic, Eigen::Dynamic > &p_mesh);

/// \brief Calculate the second member of the normal problem \f$ A^t b \f$
/// \param p_particles   Particles used for regression
/// \param p_simToCell   To each particle associates the cell it belongs to.
/// \param p_mesh        Description of the mesh coordinates for each mesh
/// \param p_fToRegress  simulations to regress
/// \return parameters of the regressed function  (
Eigen::ArrayXXd localLinearDiscrLastDimSecondMemberCalculation(const Eigen::ArrayXXd &p_particles,
        const Eigen::ArrayXi &p_simToCell,
        const Eigen::Array< std::array< double, 2>, Eigen::Dynamic, Eigen::Dynamic > &p_mesh,
        const Eigen::ArrayXXd &p_fToRegress);

/// \brief Calculate the second member of the normal problem \f$ A^t b \f$
/// \param p_particles   Particles used for regression
/// \param p_SimulBelongingToCell  simulations belonging to the cell
/// \param p_mesh        Description of coordinates of the cell
/// \param p_fToRegress  simulations to regress corresponding to particles belonging to the cell
/// \return parameters of the regressed function  (
Eigen::ArrayXXd localLinearDiscrLastDimSecondMemberCalculationOneCell(const Eigen::ArrayXXd &p_particles,
        const std::vector<int>   &p_SimulBelongingToCell,
        const Eigen::Ref<const Eigen::Array< std::array<double, 2 >, Eigen::Dynamic, 1 > >  &p_mesh,
        const Eigen::ArrayXXd &p_fToRegress);


/// \brief Reconstruct the conditional expectation for the simulations used in the optimization part (multiple functions possible)
/// \param p_particles   Particles used for regression
/// \param p_simToCell   To each particle associates the cell it belongs to.
/// \param p_mesh        Description of the mesh coordinates for each mesh
/// \param p_foncBasisCoef Coefficients associated to each function basis (first dimension is the function number)
/// \return  Reconstructed conditional expectation for each simulation (first dimension is the function number)
Eigen::ArrayXXd  localLinearDiscrLastDimReconstruction(const Eigen::ArrayXXd &p_particles, const Eigen::ArrayXi &p_simToCell,
        const Eigen::Array< std::array< double, 2>, Eigen::Dynamic, Eigen::Dynamic > &p_mesh,
        const Eigen::ArrayXXd   &p_foncBasisCoef);

/// \brief Reconstruct the conditional expectation for a given simulation  used in the optimization part
/// \param p_isim        Particle number used
/// \param p_particles   Particles used for regression
/// \param p_simToCell   To each particle associates the cell it belongs to.
/// \param p_mesh        Description of the mesh coordinates for each mesh
/// \param p_foncBasisCoef Coefficients associated to the  function basis
/// \return  Reconstructed conditional expectation for the current simulation
double  localLinearDiscrLastDimReconstructionASim(const int &p_isim, const Eigen::ArrayXXd &p_particles, const Eigen::ArrayXi &p_simToCell,
        const Eigen::Array< std::array< double, 2>, Eigen::Dynamic, Eigen::Dynamic > &p_mesh,
        const Eigen::ArrayXd   &p_foncBasisCoef);



/// \brief  Calculate conditional expectation for one particular particle
/// \param p_oneParticle  One point where to calculate the conditional expectation (multiple functions possible)
/// \param p_mesh1D        Discretization of the domain in each dimension
/// \param p_foncBasisCoef   Coefficients associated to each function basis (first dimension is the function number)
/// \return  Solution calculated for each function
Eigen::ArrayXd  localLinearDiscrLastDimReconstructionOnePoint(const Eigen::ArrayXd &p_oneParticle,
        const std::vector< std::shared_ptr< Eigen::ArrayXd >  > &p_mesh1D,
        const Eigen::ArrayXXd   &p_foncBasisCoef);

/// \brief  Permit to reconstruct a value  for a state (pt_stock, uncertainty) : the uncertainty  and pt_stock are
///         given and interpolator are given to interpolate function basis at pt_stock
/// \param p_oneParticle      One uncertainty point where to calculate the function (using regression)
/// \param p_stock            One stock point where to calculate the function
/// \param p_interpBaseFunc   spectral interpolator to interpolate the basis functions used in regression at p_stock
/// \param p_mesh1D           Discretization of the domain in each dimension
/// \return  Solution calculated  at  (pt_stock, uncertainty)
double  localLinearDiscrLastDimReconsOnePointSimStock(const Eigen::ArrayXd &p_oneParticle, const Eigen::ArrayXd &p_stock,
        const std::vector< std::shared_ptr<InterpolatorSpectral> > &p_interpBaseFunc,
        const std::vector< std::shared_ptr< Eigen::ArrayXd >  > &p_mesh1D);

/// \brief Given a particle and the coordinates of the mesh it belongs to, reconstruct the regressed value
/// \param p_oneParticle  The particle generated
/// \param p_mesh  the mesh it belongs to (each coordinate of the particle should be inside the low and max coordinate of the cell in each dimension)
///                In each dimension  gives the low and high coordinates values
/// \param p_foncBasisCoef  function basis on the current cell (the row is the number of reconstruction to achieve, the columns the number of function basis)
/// \return send back the array of reconstructed values
Eigen::ArrayXd  localLinearDiscrLastDimReconstructionOnePointOneCell(const Eigen::ArrayXd &p_oneParticle,
        const Eigen::Array< std::array< double, 2>, Eigen::Dynamic, 1 > &p_mesh,
        const Eigen::ArrayXXd   &p_foncBasisCoef);

/**@}*/
}
#endif /*LOCALLINEARDISCRLASTDIMMATRIXOPERATION_H*/