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
|
// Copyright (C) 2016, 2017 EDF
// All Rights Reserved
// This code is published under the GNU Lesser General Public License (GNU LGPL)
#ifndef LOCALSAMESIZEREGRESSION_H
#define LOCALSAMESIZEREGRESSION_H
#include <vector>
#include <memory>
#include <array>
#include <Eigen/Dense>
#include "StOpt/regression/BaseRegression.h"
#include "StOpt/core/grids/InterpolatorSpectral.h"
/** \file LocalSameSizeRegression.h
* \brief Base class for local regressions with size of constant size
* Useful mostly to check theoretical convergence rate obtained.
* \author Xavier Warin
*/
namespace StOpt
{
/**
* \defgroup LocalSameSize Piecewise regression with mesh of same size
* \brief It implements local local regression
*@{
*/
/// \class LocalSameSizeRegression LocalSameSizeRegression.h
/// To be used in Monte Carlo methods regression on each cell with cells with same size
class LocalSameSizeRegression : public BaseRegression
{
protected :
Eigen::ArrayXd m_lowValues ; ///< minimal value of the mesh in each direction
Eigen::ArrayXd m_step; ///< Step in each direction
Eigen::ArrayXi m_nbStep ; ///< Number of steps in each dimension
int m_nbMeshTotal ; ///< Total number of meshes
std::vector< Eigen::Array2i > m_simAndCell; ///< For element if gives the active simulation and the cell it belongs to
Eigen::ArrayXi m_simToCell; ///< if each simulation give it cell number of -1 if the simulation doesn't belong to any cell
/// \brief fill in arrays m_simAndCell, m_simToCell
/// \param p_particles particles used for the meshes.
/// First dimension : dimension of the problem,
void fillInSimCell(const Eigen::ArrayXXd &p_particles);
/// \brief find point location
/// \param p_coordinates point coordinates
/// \return cell where the point below to or -1
inline int pointLocation(const Eigen::ArrayXd &p_coordinates) const
{
int iret = 0;
int idec = 1 ;
for (int id = 0; id < p_coordinates.size(); ++id)
{
double position = (p_coordinates(id) - m_lowValues(id)) / m_step(id);
if ((position < m_nbStep(id)) && (position >= 0.))
iret += static_cast<int>(position) * idec;
else
return -1;
idec *= m_nbStep(id);
}
return iret;
}
public :
/// \brief Default ructor
LocalSameSizeRegression() {}
/// \brief Constructor
/// \param p_lowValues in each dimension minimal value of the grid
/// \param p_step in each dimension the step size
/// \param p_nbStep in each dimension the number of steps
LocalSameSizeRegression(const Eigen::ArrayXd &p_lowValues, const Eigen::ArrayXd &p_step, const Eigen::ArrayXi &p_nbStep);
/// \brief Constructor for object constructed at each time step
/// \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_lowValues in each dimension minimal value of the grid
/// \param p_step in each dimension the step size
/// \param p_nbStep in each dimension the number of steps
LocalSameSizeRegression(const bool &p_bZeroDate,
const Eigen::ArrayXXd &p_particles,
const Eigen::ArrayXd &p_lowValues,
const Eigen::ArrayXd &p_step,
const Eigen::ArrayXi &p_nbStep);
/// \brief Constructor only used for serialization for simulation part
/// \param p_bZeroDate first date is 0?
/// \param p_lowValues in each dimension minimal value of the grid
/// \param p_step in each dimension the step size
/// \param p_nbStep in each dimension the number of steps
/// \param p_meanX scaled factor in each direction (average of particles values in each direction)
LocalSameSizeRegression(const bool &p_bZeroDate,
const Eigen::ArrayXd &p_lowValues,
const Eigen::ArrayXd &p_step,
const Eigen::ArrayXi &p_nbStep);
/// \brief Copy constructor
/// \param p_object object to copy
LocalSameSizeRegression(const LocalSameSizeRegression &p_object) ;
/// \name get back value
///@{
const Eigen::ArrayXd &getLowValues() const
{
return m_lowValues ;
}
const Eigen::ArrayXd &getStep() const
{
return m_step ;
}
const Eigen::ArrayXi &getNbStep() const
{
return m_nbStep;
}
///@}
/// \brief get number of steps in one direction
/// \param p_id dimension concerned
inline int getNbStep(const int &p_id) const
{
return m_nbStep(p_id);
}
///@}
/// \brief get back the simulations belong to the cells and the cell numbers
inline const std::vector< Eigen::Array2i > &getSimAndCell() const
{
return m_simAndCell ;
};
/// \brief get back for each simulation its cell number
inline const Eigen::ArrayXi &getSimToCell() const
{
return m_simToCell ;
};
/// \brief get the total number of meshes
inline int getNbMeshTotal() const
{
return m_nbMeshTotal;
}
/// \brief get the total number of cells
inline int getNumberOfCells() const
{
if (m_nbStep.size() > 0)
return m_nbMeshTotal;
return 1;
};
};
/**@}*/
}
#endif /*LOCALSAMESIZEREGRESSION_H*/
|