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 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399
|
/* Copyright (c) 2020, Dyssol Development Team. All rights reserved. This file is part of Dyssol. See LICENSE file for license information. */
#pragma once
#include "TDArray.h"
#include "DenseMDMatrix.h"
#include "Matrix2D.h"
#include "TransformMatrix.h"
#include "MDMatrCacher.h"
#define DATA_SAVE_BLOCK 100
class CH5Handler;
struct sFraction
{
CTDArray tdArray;
sFraction *pNext;
};
/** This class is used to describe multidimensional distributed data. All data depends on time.
* One matrix is defined for all time points.*/
class CMDMatrix
{
private:
static const unsigned m_cnSaveVersion{ 2 };
std::vector<unsigned> m_vDimensions; ///< Types of the distributions
std::vector<unsigned> m_vClasses; ///< Number of classes of the distributions
std::vector<double> m_vTimePoints; ///< Vector of current time points
mutable sFraction *m_data{ nullptr }; ///< Current data itself
double m_dMinFraction{ DEFAULT_MIN_FRACTION }; ///< Minimal fraction. All smaller values are interpreted as 0
// ===== Variables for temporary use in recursive functions
mutable double m_dTempT1{ 0.0 };
double m_dTempT2{ 0.0 };
mutable std::vector<unsigned> m_vTempCoords;
mutable std::vector<unsigned> m_vTempDims;
mutable std::vector<double> m_vTempValues;
double m_dTempValue{ 0.0 };
CMDMatrix* m_pSortMatr{ nullptr };
std::wstring m_sCachePath{ L"" };
bool m_bCacheEnabled{ false };
mutable CMDMatrCacher m_cacheHandler;
mutable unsigned m_nCounter{ 0 };
mutable double m_dCurrWinStart{ 0.0 };
mutable double m_dCurrWinEnd{ 0.0 };
size_t m_nCacheWindow{ DEFAULT_CACHE_WINDOW };
mutable unsigned m_nNonCachedTPNum{ 0 };
mutable size_t m_nCurrOffset{ 0 };
mutable bool m_bCacheCoherent{ false };
public:
CMDMatrix() = default;
CMDMatrix(const CMDMatrix& _other);
~CMDMatrix();
/** Clears all data, time points and dimensions.*/
void Clear();
// ========== Functions to work with DIMENSIONS
/** Adds new dimension to a matrix. If dimension already exists, than nothing will be done.*/
void AddDimension( unsigned _nDim, unsigned _nClasses );
/** Removes specified dimension and erases ALL data.*/
void DeleteDimension( unsigned _nDim );
/** Removes specified dimensions and erases ALL data.*/
void DeleteDimensions( const std::vector<unsigned>& _vDims );
/** Returns vector with all current dimensions.*/
std::vector<unsigned> GetDimensions() const;
/** Returns vector with classes of all current dimensions.*/
std::vector<unsigned> GetClasses() const;
/** Sets dimension to a matrix and erases ALL data.*/
void SetDimension( unsigned _nDim, unsigned _nClasses );
/** Sets dimensions to a matrix and erases ALL data.*/
void SetDimensions( const std::vector<unsigned>& _vDims, const std::vector<unsigned>& _vClasses );
/** Returns current number of dimensions.*/
size_t GetDimensionsNumber() const;
/** Returns dimension by index. Returns -1 if there is no such dimension.*/
int GetDimensionTypeByIndex( unsigned _nIndex ) const;
/** Returns dimension size by index. Returns -1 if there is no such dimension.*/
int GetDimensionSizeByIndex( unsigned _nIndex ) const;
/** Returns dimension size by type. Returns 0 if there is no such dimension.*/
unsigned GetDimensionSizeByType( unsigned _nDim ) const;
/** Updates matrix with new set of dimensions.*/
void UpdateDimensions( const std::vector<unsigned>& _vDims, const std::vector<unsigned>& _vClasses );
/** Adds class to a dimension.*/
void AddClass( unsigned _nDim );
/** Removes specified class from dimension.*/
void RemoveClass(unsigned _nDim, size_t _nClassIndex);
/** Returns true if matrices have the same dimensions.*/
bool CompareDims( const CMDMatrix& _matr ) const;
// ========== Functions to work with TIME POINTS
/** Adds new time point. Data to a new time point will be copied from _dSourceTimePoint.
* If _dSourceTimePoint == -1, than data will be copied from the previous time point. If such time point already exists, than nothing will be done.*/
void AddTimePoint( double _dTime, double _dSrcTimePoint = -1 );
/** Changes specified time point. For UI purposes.*/
void ChangeTimePoint( unsigned _nTimePointIndex, double _dNewTime );
/** Removes specified time point.*/
void RemoveTimePoint( double _dTime );
/** Removes time points from interval, including or excluding boundaries.*/
void RemoveTimePoints( double _dStart, double _dEnd, bool _inclusive = true);
/** Removes all time points after the specified time.
* If _bIncludeTime == true, than _dTime will be includede in this time interval.*/
void RemoveTimePointsAfter( double _dTime, bool _bIncludeTime = false );
/** Removes all defined time points.*/
void RemoveAllTimePoints();
/** Returns vector of time points for interval.*/
std::vector<double> GetTimePoints( double _dStart, double _dEnd ) const;
/** Returns all defined time points.*/
std::vector<double> GetAllTimePoints() const;
/** Returns time point with specified index. If there is no such time point, than -1 will be returned.*/
double GetTimeForIndex( unsigned _nIndex ) const;
/** Returns number of defined time points.*/
size_t GetTimePointsNumber() const;
// ========== Functions to work with MINIMAL FRACTION
/** Returns current minimal fraction.*/
double GetMinimalFraction() const;
/** Sets minimal fraction for all matrices.s*/
void SetMinimalFraction( double _dValue );
// ========== Functions to GET and SET VALUES
/** Returns single value according to specified dimension and coordinate for time point with specified index.
* A dimensions set can be reduced. If dimension, coordinate or index doesn't exist, than 0 will be returned.
* Works with non-negative values only. If value is less then m_dMinFraction, than 0 will be returned.*/
double GetValue(size_t _nTimeIndex, unsigned _nDim, unsigned _nCoord) const;
/** Returns single value according to specified dimension and coordinate for existing time point _dTime
* or its approximation if _dTime wasn't defined. A dimensions set can be reduced.
* If dimension or coordinate doesn't exist, than 0 will be returned. Works with non-negative values only.
* If value is less then m_dMinFraction, than 0 will be returned.*/
double GetValue(double _dTime, unsigned _nDim, unsigned _nCoord) const;
/** Returns single value according to specified dimensions and coordinates for existing time point _dTime
* or its approximation if _dTime wasn't defined. A dimensions set can be reduced.
* If dimension or coordinate doesn't exist, than 0 will be returned. Works with non-negative values only.
* If value is less then m_dMinFraction, than 0 will be returned.*/
double GetValue(double _dTime, unsigned _nDim1, unsigned _nCoord1, unsigned _nDim2, unsigned _nCoord2) const;
/** Returns single value according to specified dimensions and coordinates for existing time point _dTime
* or its approximation if _dTime wasn't defined. A dimensions set can be reduced.
* If dimension or coordinate doesn't exist, than 0 will be returned. Works with non-negative values only.
* If value is less then m_dMinFraction, than 0 will be returned.*/
double GetValue(double _dTime, unsigned _nDim1, unsigned _nCoord1, unsigned _nDim2, unsigned _nCoord2, unsigned _nDim3, unsigned _nCoord3) const;
/** Returns single value according to all defined dimensions and specified coordinates for existing time point _dTime
* or its approximation if _dTime wasn't defined.
* If coordinate doesn't exist, than 0 will be returned. Works with non-negative values only.
* If value is less then m_dMinFraction, than 0 will be returned.*/
double GetValue(double _dTime, const std::vector<unsigned>& _vCoords) const;
/** Returns single value according to specified dimensions and coordinates for existing time point _dTime
* or its approximation if _dTime wasn't defined. A dimensions set can be reduced.
* If dimension or coordinate doesn't exist, than 0 will be returned. Works with non-negative values only.
* If value is less then m_dMinFraction, than 0 will be returned.*/
double GetValue(double _dTime, const std::vector<unsigned>& _vDims, const std::vector<unsigned>& _vCoords) const;
/** Sets value according to specified dimension and coordinate for the time point with specified index.
* Dimensions set can be reduced. All negative values will be transformed to a 0. Returns false on error.*/
bool SetValue( unsigned _nTimeIndex, unsigned _nDim, unsigned _nCoord, double _dValue, bool _bExternal = true );
/** Sets value according to specified dimension and coordinate for the existing time point. Dimensions set can be reduced.
* If time point wasn't defined, than nothing will be done. All negative values will be transformed to a 0. Returns false on error.*/
bool SetValue( double _dTime, unsigned _nDim, unsigned _nCoord, double _dValue, bool _bExternal = true );
/** Sets value according to specified dimensions and coordinates for the existing time point. Dimensions set can be reduced.
* If time point wasn't defined, than nothing will be done. All negative values will be transformed to a 0. Returns false on error.*/
bool SetValue( double _dTime, unsigned _nDim1, unsigned _nCoord1, unsigned _nDim2, unsigned _nCoord2, double _dValue, bool _bExternal = true );
/** Sets value according to specified dimensions and coordinates for the existing time point. Dimensions set can be reduced.
* If time point wasn't defined, than nothing will be done. All negative values will be transformed to a 0. Returns false on error.*/
bool SetValue( double _dTime, unsigned _nDim1, unsigned _nCoord1, unsigned _nDim2, unsigned _nCoord2, unsigned _nDim3, unsigned _nCoord3, double _dValue, bool _bExternal = true );
/** Sets value according to all defined dimensions and specified coordinates for the existing time point.
* If time point wasn't defined, than nothing will be done. All negative values will be transformed to a 0. Returns false on error.*/
bool SetValue( double _dTime, const std::vector<unsigned>& _vCoords, double _dValue, bool _bExternal = true );
/** Sets value according to specified dimensions and coordinates for the existing time point. Dimensions set can be reduced.
* If time point wasn't defined, than nothing will be done. All negative values will be transformed to a 0. Returns false on error.*/
bool SetValue( double _dTime, const std::vector<unsigned>& _vDims, const std::vector<unsigned>& _vCoords, double _dValue, bool _bExternal = true );
// ========== Functions to GET and SET VECTORS
/** Returns vector value according to specified dimension for time point with specified index.
* A dimensions set can be reduced. Returns false on error and on vector of zeros in _vResult.*/
std::vector<double> GetVectorValue(unsigned _nTimeIndex, unsigned _nDim) const;
/** Returns vector value according to specified dimension for time point with specified index.
* A dimensions set can be reduced. Returns false on error and on vector of zeros in _vResult.*/
bool GetVectorValue(unsigned _nTimeIndex, unsigned _nDim, std::vector<double>& _vResult) const;
/** Returns vector value according to specified dimension for existing time point _dTime
* or its approximation if _dTime wasn't defined. A dimensions set can be reduced. */
std::vector<double> GetVectorValue(double _dTime, unsigned _nDim) const;
/** Returns vector value according to specified dimension for existing time point _dTime
* or its approximation if _dTime wasn't defined. A dimensions set can be reduced. Returns false on error and on vector of zeros in _vResult.*/
bool GetVectorValue(double _dTime, unsigned _nDim, std::vector<double>& _vResult) const;
/** Returns vector value according to specified dimensions and coordinate for existing time point _dTime
* or its approximation if _dTime wasn't defined. A dimensions set can be reduced. */
std::vector<double> GetVectorValue(double _dTime, unsigned _nDim1, unsigned _nCoord1, unsigned _nDim2) const;
/** Returns vector value according to specified dimensions and coordinate for existing time point _dTime
* or its approximation if _dTime wasn't defined. A dimensions set can be reduced. Returns false on error and on vector of zeros in _vResult.*/
bool GetVectorValue(double _dTime, unsigned _nDim1, unsigned _nCoord1, unsigned _nDim2, std::vector<double>& _vResult) const;
/** Returns vector value according to specified dimensions and coordinate for existing time point _dTime
* or its approximation if _dTime wasn't defined. A dimensions set can be reduced. Number of coordinates must be one less then number of dimensions.*/
std::vector<double> GetVectorValue(double _dTime, const std::vector<unsigned>& _vDims, const std::vector<unsigned>& _vCoords) const;
/** Returns vector value according to specified dimensions and coordinate for existing time point _dTime
* or its approximation if _dTime wasn't defined. A dimensions set can be reduced. Number of coordinates must be one less then number of dimensions.
* Returns false on error and on vector of zeros in _vResult.*/
bool GetVectorValue(double _dTime, const std::vector<unsigned>& _vDims, const std::vector<unsigned>& _vCoords, std::vector<double>& _vResult) const;
std::vector<double> GetValues(unsigned _nDim, unsigned _nCoord) const;
/** Sets vector value according to specified dimension for the time point with specified index.
* Dimensions set can be reduced. Sets value only if the value of the previous level (dimension in structure) is already set.
* Returns false on error.*/
bool SetVectorValue( unsigned _nTimeIndex, unsigned _nDim, const std::vector<double>& _vValue, bool _bExternal = false );
/** Sets vector value according to specified dimension for the existing time point. Dimensions set can be reduced.
* If time point wasn't defined, than nothing will be done. Sets value only if the value of the previous level (dimension in structure)
* is already set. Returns false on error.*/
bool SetVectorValue( double _dTime, unsigned _nDim, const std::vector<double>& _vValue, bool _bExternal = false );
/** Sets vector value according to specified dimensions and coordinate for the existing time point. Dimensions set can be reduced.
* If time point wasn't defined, than nothing will be done. Sets value only if the value of the previous level (dimension in structure)
* is already set. Returns false on error.*/
bool SetVectorValue( double _dTime, unsigned _nDim1, unsigned _nCoord1, unsigned _nDim2, const std::vector<double>& _vValue, bool _bExternal = false );
/** Sets vector value according to specified dimensions and coordinates for the existing time point. Dimensions set can be reduced.
* If time point wasn't defined, than nothing will be done. Number of coordinates must be one less then number of dimensions.
* Sets value only if the value of the previous level (dimension in structure) is already set. Returns false on error.*/
bool SetVectorValue( double _dTime, const std::vector<unsigned>& _vDims, const std::vector<unsigned>& _vCoords, const std::vector<double>& _vValue, bool _bExternal = false );
// ========== Functions to GET and SET MATRICES
/** Gets 2-dimensional matrix according to last 2 dimensions in _vDimType. Returns false on error.
* Must be ( _vDimType.size() == _vCoord.size() + 2 ) and full size of dimensions vector. For UI purposes.*/
bool GetMatrixValue( double _dTime, const std::vector<unsigned>& _vDims, const std::vector<unsigned>& _vCoords, std::vector<std::vector<double>>& _vResult );
/** Gets 2-dimensional matrix according to last 2 dimensions in _vDimType. Returns false on error.
* Must be ( _vDimType.size() == _vCoord.size() + 2 ) and full size of dimensions vector. For UI purposes.*/
std::vector<std::vector<double>> GetMatrixValue(double _dTime, const std::vector<unsigned>& _vDims, const std::vector<unsigned>& _vCoords);
/** Sets 2-dimensional matrix according to last 2 dimensions in _vDimType. Returns false on error.
* Must be ( _vDimType.size() == _vCoord.size() + 2 ) and full size of dimensions vector. For UI purposes.*/
bool SetMatrixValue( double _dTime, const std::vector<unsigned>& _vDims, const std::vector<unsigned>& _vCoords, const std::vector<std::vector<double>>& _vValue );
// ========== Functions to GET and SET DISTRIBUTIONS
/** Return distribution by all defined dimensions and specified time. Uses GetVectorValue() functions.*/
CDenseMDMatrix GetDistribution(double _dTime) const;
/** Return distribution by all defined dimensions and specified time. Uses GetVectorValue() functions.*/
bool GetDistribution(double _dTime, CDenseMDMatrix& _Result) const;
/** Return distribution by specified dimension and time. Uses GetVectorValue() functions.*/
bool GetDistribution(double _dTime, unsigned _nDim, std::vector<double>& _vResult) const;
/** Return distribution by specified dimension and time. Uses GetVectorValue() functions.*/
std::vector<double> GetDistribution(double _dTime, unsigned _nDim) const;
/** Return distribution by specified dimensions and time. Uses GetVectorValue() functions.*/
bool GetDistribution(double _dTime, unsigned _nDim1, unsigned _nDim2, CMatrix2D& _Result) const;
/** Return distribution by specified dimensions and time. Uses GetVectorValue() functions.*/
CMatrix2D GetDistribution(double _dTime, unsigned _nDim1, unsigned _nDim2) const;
/** Return distribution by specified dimensions and time. Uses GetVectorValue() functions.*/
bool GetDistribution(double _dTime, unsigned _nDim1, unsigned _nDim2, unsigned _nDim3, CDenseMDMatrix& _Result) const;
/** Return distribution by specified dimensions and time. Uses GetVectorValue() functions.*/
CDenseMDMatrix GetDistribution(double _dTime, unsigned _nDim1, unsigned _nDim2, unsigned _nDim3) const;
/** Return distribution by specified dimensions and time. Uses GetVectorValue() functions.*/
bool GetDistribution(double _dTime, const std::vector<unsigned>& _vDims, CDenseMDMatrix& _Result) const;
/** Return distribution by specified dimensions and time. Uses GetVectorValue() functions.*/
CDenseMDMatrix GetDistribution(double _dTime, const std::vector<unsigned>& _vDims) const;
/** Sets distribution by specified dimension and time. Uses SetVectorValue() functions.*/
bool SetDistribution( double _dTime, unsigned _nDim, const std::vector<double>& _vDistr );
/** Sets distribution by specified dimensions and time. Uses SetVectorValue() functions.*/
bool SetDistribution( double _dTime, unsigned _nDim1, unsigned _nDim2, const CMatrix2D& _Distr );
/** Sets distribution by specified dimensions and time. Uses SetVectorValue() functions.*/
bool SetDistribution( double _dTime, const CDenseMDMatrix& _Distr );
// ========== Functions to work with matrix TRANSFORMATIONS
/** Transforms matrix according to matrix _Transformation. Returns false on error.*/
bool Transform( double _dTime, const CTransformMatrix& _TMatrix );
// ========== Functions for matrix NORMALIZATION
/** Normalizes data in matrix for specified time point. If time point wasn't defined, than nothing will be done.*/
void NormalizeMatrix( double _dTime );
/** Normalizes data in matrix for time interval (incl).*/
void NormalizeMatrix( double _dStart, double _dEnd );
/** Normalizes data in matrix for all time points.*/
void NormalizeMatrix();
/** Returns true if matrix is normalized in specified time point.*/
bool IsNormalized( double _dTime );
// ========== Functions to work with ANOTHER MATRICES
/** Copy data from another MDMatrix on specified time point. If time point doesn't exist, it will be created.*/
bool CopyFrom( const CMDMatrix& _Source, double _dTime );
/** Copy data from another MDMatrix on specified time interval.*/
bool CopyFrom( const CMDMatrix& _Source, double _dStart, double _dEnd );
/** Copy data from time point _dTimeSrc of another MDMatrix to time point _dTimeDest of this matrix.
* If time point doesn't exist, it will be created.*/
bool CopyFromTimePoint( const CMDMatrix& _Source, double _dTimeSrc, double _dTimeDest );
///** Adds _srcMatr to a matrix with specified factors for time point _dTime.*/
//void AddMatrix( CMDMatrix& _srcMatr, double _dFactorDst, double _dFactorSrc, double _dTime );
///** Adds _srcMatr to a matrix with specified factors for time intervals.*/
//void AddMatrix( CMDMatrix& _srcMatr, const std::vector<double>& _vFactorsDst, const std::vector<double>& _vFactorsSrc, const std::vector<double>& _vTimePoints );
// ========== Functions to SAVE / LOAD matrix
/** Save data to file.*/
void SaveToFile( CH5Handler& _h5File, const std::string& _sPath ) const;
void SaveMDBlockToFile(CH5Handler& _h5File, const std::string& _sPath, unsigned _iFirst, unsigned _iLast, std::vector<std::vector<double>>& _vvBuf) const;
/** Load data from file*/
void LoadFromFile(const CH5Handler& _h5File, const std::string& _sPath );
void LoadMDBlockFromFile(const CH5Handler& _h5File, const std::string& _sPath, unsigned _iFirst, unsigned _iLast, std::vector<std::vector<double>>& vvBuf);
void SetCachePath(const std::wstring& _sPath);
void SetCacheParams(bool _bEnabled, size_t _nWindow);
/** Removes all data, which can be approximated.*/
void CompressData( double _dStartTime, double _dEndTime, double _dATol, double _dRTol );
void ExtrapolateToPoint( double _dT1, double _dT2, double _dTExtra );
void ExtrapolateToPoint( double _dT0, double _dT1, double _dT2, double _dTExtra );
private:
/** Checks the duplicates in vector. Return true if there are no duplicates.*/
bool CheckDuplicates( const std::vector<unsigned>& _vVec ) const;
/** Returns index of time point. Strict search returns -1 if there is no such time, not strict search returns index to paste.*/
unsigned GetTimeIndex(double _dTime, bool _bIsStrict = true) const;
/** Initializes current fraction with specified size.*/
sFraction* IitialiseDimension( unsigned _nSize ) const;
/* Increments last coordinate for getting/setting vectors according to a dimensions set. Must be _vCoords.size()+1 == _vSizes.size().
* Returns false if the end is reached.*/
bool IncrementCoords( std::vector<unsigned>& _vCoords, const std::vector<unsigned>& _vSizes ) const;
/** Sorts matrix for time point _dSrcTime according to the order of dimensions in _dstMatrix and sets it to time-point _dDstTime.
* Returns false on error.*/
bool SortMatrix( double _dSrcTime, double _dDstTime, CMDMatrix& _dstMatrix );
/** Sorts matrix according to order in _dstMatrix for all time points.*/
bool SortMatrix( CMDMatrix& _dstMatrix );
/** Removes dimensions and puts new sorted matrix in _sortMatr. New dimensions and classes set will be in _vNewDims and _vNewClasses.*/
void DeleteDimsWithSort( const std::vector<unsigned>& _vDims, std::vector<unsigned>& _vNewDims, std::vector<unsigned>& _vNewClasses, CMDMatrix& _sortMatr );
// ========== Recursive functions
void NormalizeWrapper( const std::vector<double>& _vTimes );
/** Copies all data from _pSrc to a current matrix.*/
sFraction* CopyFractionsRecursive( sFraction *_pSrc, unsigned _nNesting = 0 );
/** Removes all fractions starting from _pFraction.*/
sFraction* RemoveFractionsRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
/** Adds specified time point m_dTempT1 to each fraction. Data to a new time point will be copied from m_dTempT2.*/
void AddTimePointRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
/** Changes time point m_dTempT1 to a m_dTempT2.*/
void ChangeTimePointRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
/** Removes time points from interval [m_dTempT1..m_dTempT2] or single time point m_dTempT1 (if m_dTempT2 == -1)*/
sFraction* RemoveTimePointsRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
/** Returns value for time m_dTempT1 according to specified dimensions m_vTempDims and coordinates m_vTempCoords. Dimensions set can be reduced.*/
double GetValueRecursive(sFraction *_pFraction, unsigned _nLevel = 1, unsigned _nNesting = 0) const;
/** Sets value m_dTempValue for time m_dTempT1 according to specified dimensions m_vTempDims and coordinates m_vTempCoords.
* Dimensions set can be reduced. If time point wasn't defined, nothing will be done. Returns false on error.*/
bool SetValueRecursive( sFraction *_pFraction, bool _bExternal = true, unsigned _nLevel = 1, unsigned _nNesting = 0 );
/** Returns vector value for time m_dTempT1 according to specified dimensions m_vTempDims and coordinates m_vTempCoords.
* Dimensions set can be reduced. Returns false on error.*/
bool GetVectorValueRecursive(sFraction *_pFraction, std::vector<double>& _vRes, unsigned _nNesting = 0) const;
/* Sets vector value m_vTempValues for time m_dTempT1 according to specified dimensions m_vTempDims and coordinates m_vTempCoords.
* Dimensions set can be reduced. Returns false on error.*/
bool SetVectorValueRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
/** Transforms matrix m_pSortMatr according to a transformation matrix _TMatr.
* Doesn't check the correspondence of dimensions between m_pSortMatr, _TMatr and this matrix.*/
void TransformRecurcive( const CTransformMatrix& _TMatr );
/** Normalizes matrix for time point m_vTempValues.*/
void NormalizeMatrixBySumRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
bool NormalizeMatrixRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
bool CheckNormalizationRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
/** Copies data from another matrix for time interval [m_dTempT1, m_dTempT2].*/
sFraction* CopyFromRecursive( sFraction *_pDest, sFraction *_pSource, unsigned _nNesting = 0 );
/** Copies data from another matrix from time point m_dTempT1 to m_dTempT2.*/
sFraction* CopyFromTimePointRecursive( sFraction *_pDest, sFraction *_pSource, unsigned _nNesting = 0 );
/** Adds new class to a dimension with index _nDimIndex.*/
sFraction* AddClassRecursive( sFraction *_pFraction, unsigned _nDimIndex, unsigned _nNesting = 0 );
/** Removes class _nClassIndex from a dimension _nDimIndex.*/
sFraction* RemoveClassRecursive(sFraction *_pFraction, unsigned _nDimIndex, size_t _nClassIndex, unsigned _nNesting = 0);
/** Removes all data, which can be approximated in time interval [m_dTempT1, m_dTempT2].*/
void CompressDataRecursive( sFraction *_pFraction, double _dATol, double _dRTol, unsigned _nNesting = 0 );
//void NormalizeToOneRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
//void NormalizeToZeroRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
//sFraction* CreateWithOneRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
sFraction* SetToZeroRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
// Recursively sets all values between selected time points to zero.
void SetToZero2Recursive(sFraction* _pFraction, unsigned _nNesting = 0);
sFraction* SetToOneRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
void Extrapolate2ToPointRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
void Extrapolate3ToPointRecursive( sFraction *_pFraction, unsigned _nNesting = 0 );
void UnCacheData(double _dTP) const;
void UnCacheData(double _dT1, double _dT2) const;
void FlushToCache() const;
void CheckCacheNeed() const;
void CacheData() const;
void CorrectWinBoundary() const;
void ClearCache() const;
sFraction* UnCacheDataRecursive( sFraction *_pFraction, std::vector<std::vector<double>>& _vData, unsigned _nNesting = 0 ) const;
void CacheDataRecursive( sFraction *_pFraction, std::vector<std::vector<double>>& _vData, unsigned _nNesting = 0 ) const;
void GetDataForSaveRecursive( sFraction *_pFraction, std::vector<std::vector<double>>& _vData, unsigned _nNesting = 0 ) const;
sFraction* SetDataForLoadRecursive( sFraction *_pFraction, std::vector<std::vector<double>>& _vData, unsigned _nNesting = 0 ) const;
};
|