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
|
// Copyright (C) 2019 EDF
// All Rights Reserved
// This code is published under the GNU Lesser General Public License (GNU LGPL)
#ifndef LOCALADAPTCELLREGRESSION_H
#define LOCALADAPTCELLREGRESSION_H
#include <Eigen/Dense>
#include <array>
#include "StOpt/regression/BaseRegression.h"
/** \file LocalAdaptCelleRegression.h
* \brief Base class for all regressors on cells with adaptation of the size
* \author Xavier Warin
*
*/
namespace StOpt
{
/// \class LocalAdaptCellRegression LocalAdaptCellRegression.h
/// General class for local basis function with adaptation in each dimension
class LocalAdaptCellRegression : public BaseRegression
{
protected :
Eigen::ArrayXi m_nbMesh ; ///< Number of discretization meshes in each direction: copy because of python binding
int m_nbMeshTotal ; ///< Total number of meshes
Eigen::Array< std::array< double, 2>, Eigen::Dynamic, Eigen::Dynamic > m_mesh ; ///< In each dimension, for each mesh, give the coordinates min and max of the mesh.
Eigen::ArrayXi m_simToCell; ///< For each simulation gives its global mesh number
std::vector< std::shared_ptr< std::vector< int> > > m_simulBelongingToCell; ///< Utility : to each cell defines the particles lying inside
/// \brief To a particle affect to cell number
/// \param p_oneParticle One point
/// \return cell number
virtual int particleToMesh(const Eigen::ArrayXd &p_oneParticle) const = 0;
public :
/// \brief Default ructor
LocalAdaptCellRegression() {}
/// \brief Constructor
/// \param p_nbMesh discretization in each direction
/// \param p_bRotationAndRecale do we use SVD
LocalAdaptCellRegression(const Eigen::ArrayXi &p_nbMesh, const bool &p_bRotationAndRecale);
/// \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_nbMesh discretization in each direction
/// \param p_bRotationAndRecale do we use SVD
LocalAdaptCellRegression(const bool &p_bZeroDate,
const Eigen::ArrayXXd &p_particles,
const Eigen::ArrayXi &p_nbMesh,
const bool &p_bRotationAndRecale);
/// \brief Second constructor , only to be used in simulation
LocalAdaptCellRegression(const bool &p_bZeroDate,
const Eigen::ArrayXi &p_nbMesh,
const Eigen::Array< std::array< double, 2>, Eigen::Dynamic, Eigen::Dynamic > &p_mesh,
const Eigen::ArrayXd &p_meanX,
const Eigen::ArrayXd &p_etypX,
const Eigen::MatrixXd &p_svdMatrix,
const bool &p_bRotationAndRecale) ;
/// \brief Copy constructor
/// \param p_object object to copy
LocalAdaptCellRegression(const LocalAdaptCellRegression &p_object);
/// \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
virtual void updateSimulations(const bool &p_bZeroDate, const Eigen::ArrayXXd &p_particles) = 0;
/// \brief Get some local accessors
///@{
inline const Eigen::ArrayXi &getNbMesh() const
{
return m_nbMesh;
}
inline const Eigen::Array< std::array< double, 2>, Eigen::Dynamic, Eigen::Dynamic > &getMesh() const
{
return m_mesh;
}
inline const Eigen::ArrayXi &getSimToCell() const
{
return m_simToCell;
}
///@}
/// \brief Get some local accessors with copy (useful for python)
///@{
inline Eigen::ArrayXi getNbMeshCopy() const
{
return m_nbMesh;
}
inline Eigen::Array< std::array< double, 2>, Eigen::Dynamic, Eigen::Dynamic > getMeshCopy() const
{
return m_mesh;
}
inline Eigen::ArrayXi getSimToCellCopy() const
{
return m_simToCell;
}
///@}
/// \brief Get dimension of the problem
inline int getDimension() const
{
return m_nbMesh.size();
}
/// \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_nbMesh.size() > 0)
return m_nbMesh.prod();
return 1;
};
/// \brief To a particle get back the cell it belongs to
/// \param p_isim particle number
/// \return cell number it belongs to
inline int getCellAssociatedToSim(const int &p_isim) const
{
if (m_simToCell.size() > 0)
return m_simToCell(p_isim);
else
return 0;
}
/// \brief calculate a vector of vector of points giving for each cell the vector of points belonging to this cell
void evaluateSimulBelongingToCell();
/// \brief get particles belonging all mesh
inline const std::vector< std::shared_ptr< std::vector< int> > > &getSimulBelongingToCell() const
{
return m_simulBelongingToCell;
}
/// \brief to a particle gives its cell number
/// \param p_oneParticle One point
inline int getMeshNumberAssociatedTo(const Eigen::ArrayXd &p_oneParticle) const
{
if (((m_nbMesh.size() > 0) && (!m_bZeroDate)) && (p_oneParticle.size() != 0))
{
// rotation
Eigen::VectorXd x = p_oneParticle.matrix();
x = ((x.array() - m_meanX) / m_etypX).matrix();
x = m_svdMatrix * x;
return particleToMesh(x.array());
}
else
return 0;
}
/// \brief conditional expectation basis function coefficient calculation for a special cell
/// \param p_iCell cell number
/// \param p_fToRegress function to regress associated to each simulation used in optimization and corresponding to the cell
/// \return regression coordinates on the basis (size : the dimension of the problem plus one)
/// @{
virtual Eigen::ArrayXd getCoordBasisFunctionOneCell(const int &p_iCell, const Eigen::ArrayXd &p_fToRegress) const = 0;
virtual Eigen::ArrayXXd getCoordBasisFunctionMultipleOneCell(const int &p_iCell, const Eigen::ArrayXXd &p_fToRegress) const = 0 ;
///@}
/// \brief Given a particle and the coordinates of the mesh it belongs to, get back the conditional expectation
/// \param p_oneParticle The particle generated
/// \param p_cell the cell it belongs to
/// \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 conditional expectations
virtual Eigen::ArrayXd getValuesOneCell(const Eigen::ArrayXd &p_oneParticle, const int &p_cell, const Eigen::ArrayXXd &p_foncBasisCoef) const = 0;
};
}
#endif
|