STK++ 0.9.13
STK::LocalVariance< Array > Class Template Reference

A LocalVariance is an implementation of the abstract ILinearReduct class. More...

#include <STK_LocalVariance.h>

Inheritance diagram for STK::LocalVariance< Array >:
Inheritance graph

Public Types

typedef hidden::Traits< Array >::Row RowVector
 
typedef ILinearReduct< Array, VectorXBase
 
- Public Types inherited from STK::ILinearReduct< Array, VectorX >
typedef IReducer< Array, VectorXBase
 

Public Member Functions

 LocalVariance (Array const *p_data, Reduct::TypeGraph type=Reduct::distance_, int nbNeighbor=1)
 Constructor.
 
 LocalVariance (Array const &data, Reduct::TypeGraph type=Reduct::distance_, int nbNeighbor=1)
 Constructor.
 
 LocalVariance (LocalVariance const &reducer)
 copy Constructor.
 
virtual LocalVarianceclone () const
 clone pattern
 
virtual ~LocalVariance ()
 Destructor.
 
int nbNeighbor () const
 
ArrayXXi constpred () const
 
ArraySquareX constcovariance () const
 
ArraySquareX constlocalCovariance () const
 
- Public Member Functions inherited from STK::ILinearReduct< Array, VectorX >
 ILinearReduct ()
 Default Constructor.
 
 ILinearReduct (Array const *p_data)
 Constructor.
 
 ILinearReduct (Array const &data)
 Constructor.
 
 ILinearReduct (ILinearReduct const &reducer)
 copy Constructor.
 
virtual ~ILinearReduct ()
 virtual destructor

 
VectorX constcriteriaValues () const
 
Array constaxis () const
 
virtual bool run ()
 Compute the projection matrix axis_ by maximizing the criteria and project the data set in order to obtain p_projected_.
 
virtual bool run (VectorX const &weights)
 Compute the projection matrix set by maximizing the weighted criteria and project the data set in order to obtain p_projected_.
 
- Public Member Functions inherited from STK::IReducer< Array, Weights >
virtual ~IReducer ()
 virtual destructor.
 
int dim () const
 get the number of dimension.
 
Array * p_reduced () const
 get a pointer on the reduced data set
 
void setDimension (const int &dim)
 set the number of dimension.
 
void clear ()
 clear allocated memory
 
- Public Member Functions inherited from STK::IRunnerUnsupervised< Array, Weights >
Array constp_data () const
 get the data set
 
virtual void setData (Array const *p_data)
 Set the data set.
 
virtual void setData (Array const &data)
 Set the data set.
 
- Public Member Functions inherited from STK::IRunnerBase
String consterror () const
 get the last error message.
 

Protected Member Functions

virtual void update ()
 Compute the proximity graph if the data set is modified.
 
virtual void maximizeStep ()
 Compute the axis by maximizing the ratio of the local variance on the total variance of the data set.
 
virtual void maximizeStep (VectorX const &weights)
 Compute the axis by maximizing the ratio of the weighted local variance on the weighted total variance of the data set.
 
void prim ()
 compute the minimal spanning tree
 
void minimalDistance ()
 compute the minimal distance graph
 
void computeCovarianceMatrices ()
 compute the covariances matrices of the data set
 
void computeCovarianceMatrices (VectorX const &weights)
 compute the weighted covariances matrices of the data set
 
- Protected Member Functions inherited from STK::IReducer< Array, Weights >
 IReducer ()
 Default constructor.
 
 IReducer (Array const *p_data)
 Constructor with a pointer on the constant data set.
 
 IReducer (Array const &data)
 Constructor with a constant reference on the data set.
 
 IReducer (IReducer const &reducer)
 Copy constructor.
 
- Protected Member Functions inherited from STK::IRunnerUnsupervised< Array, Weights >
 IRunnerUnsupervised ()
 default constructor.
 
 IRunnerUnsupervised (Array const *const p_data)
 constructor with a pointer on the constant data set
 
 IRunnerUnsupervised (Array const &data)
 constructor with a constant reference on the data set
 
 IRunnerUnsupervised (IRunnerUnsupervised const &runner)
 copy constructor
 
 ~IRunnerUnsupervised ()
 destructor
 
- Protected Member Functions inherited from STK::IRunnerBase
 IRunnerBase ()
 default constructor
 
 IRunnerBase (IRunnerBase const &runner)
 copy constructor
 
virtual ~IRunnerBase ()
 destructor
 

Protected Attributes

Reduct::TypeGraph type_
 number of neighbors
 
int nbNeighbor_
 number of neighbors
 
ArrayXXi neighbors_
 Predecessor array : store the spanning tree or the minimal distance to the neighbors.
 
ArrayXX dist_
 distance matrix : store the minimal distance to the neighbors
 
ArraySquareX localCovariance_
 the local covariance Array
 
ArraySquareX covariance_
 the covariance Array
 
- Protected Attributes inherited from STK::ILinearReduct< Array, VectorX >
VectorX idx_values_
 The values of the index for each axis.
 
Array axis_
 The computed axis.
 
- Protected Attributes inherited from STK::IReducer< Array, Weights >
int dim_
 dimension of the reduced data set
 
Array * p_reduced_
 The reduced data set.
 
- Protected Attributes inherited from STK::IRunnerUnsupervised< Array, Weights >
Array constp_data_
 A pointer on the original data set.
 
- Protected Attributes inherited from STK::IRunnerBase
String msg_error_
 String with the last error message.
 
bool hasRun_
 true if run has been used, false otherwise
 

Private Member Functions

void computeAxis ()
 compute the axis using the first eigenvectors of the matrix of covariance and local covariance
 
 LocalVariance (Reduct::TypeGraph type=Reduct::distance_, int nbNeighbor=1)
 Default constructor.
 

Additional Inherited Members

- Public Attributes inherited from STK::ILinearReduct< Array, VectorX >
Array * p_reduced_
 The reduced data set.
 
- Protected Types inherited from STK::IReducer< Array, Weights >
typedef IRunnerUnsupervised< Array, WeightsRunner
 

Detailed Description

template<class Array>
class STK::LocalVariance< Array >

A LocalVariance is an implementation of the abstract ILinearReduct class.

A LocalVariance is an Index which maximize the projected local variance of the data set. The class can use the algorithm of Prim or the minimal distance in order to compute the proximity graph defining the local variance.

This class derive from ILinearReduct which derive itself from IRunnerUnsupervised. The run() and run(weights) methods have been implemented in the ILinearReduct interface using the pure virtual methods maximizeStep() and MaximizeCriteria(weights).

Definition at line 62 of file STK_LocalVariance.h.

Member Typedef Documentation

◆ Base

template<class Array >
typedef ILinearReduct<Array, VectorX> STK::LocalVariance< Array >::Base

Definition at line 66 of file STK_LocalVariance.h.

◆ RowVector

template<class Array >
typedef hidden::Traits<Array>::Row STK::LocalVariance< Array >::RowVector

Definition at line 65 of file STK_LocalVariance.h.

Constructor & Destructor Documentation

◆ LocalVariance() [1/4]

template<class Array >
STK::LocalVariance< Array >::LocalVariance ( Array const p_data,
Reduct::TypeGraph  type = Reduct::distance_,
int  nbNeighbor = 1 
)

Constructor.

the TypeGraph and the number of neighbors are given by the user and cannot modified.

Parameters
p_datapointyer on the data set to process
typetype of proximity graph to build
nbNeighbornumber of neighbors to use in the proximity graph

Definition at line 165 of file STK_LocalVariance.h.

167 : Base(p_data)
168 , type_(type)
170 , neighbors_()
171 , dist_()
172
173{
174 if (!p_data) return;
176 dist_.resize(p_data->rows(), Range(1,nbNeighbor_));
177 dist_ = Arithmetic<Real>::max();
178 // compute minimal proximity graph of the data set
179 switch (type_)
180 {
181 case Reduct::prim_:
182 prim();
183 break;
186 break;
188 STKRUNTIME_ERROR_NO_ARG(LocalVariance::LocalVariance(data, type, nbNeighbor),unknown proximity graph);
189 break;
190 };
191}
#define STKRUNTIME_ERROR_NO_ARG(Where, Error)
Definition STK_Macros.h:138
Derived & resize(Range const &I, Range const &J)
resize the array.
Array const * p_data() const
get the data set
ILinearReduct< Array, VectorX > Base
void minimalDistance()
compute the minimal distance graph
Reduct::TypeGraph type_
number of neighbors
LocalVariance(Array const *p_data, Reduct::TypeGraph type=Reduct::distance_, int nbNeighbor=1)
Constructor.
ArrayXXi neighbors_
Predecessor array : store the spanning tree or the minimal distance to the neighbors.
void prim()
compute the minimal spanning tree
int nbNeighbor_
number of neighbors
ArrayXX dist_
distance matrix : store the minimal distance to the neighbors
TRange< UnknownSize > Range
Definition STK_Range.h:59

References STK::LocalVariance< Array >::dist_, STK::Reduct::distance_, STK::LocalVariance< Array >::LocalVariance(), STK::LocalVariance< Array >::minimalDistance(), STK::LocalVariance< Array >::nbNeighbor(), STK::LocalVariance< Array >::nbNeighbor_, STK::LocalVariance< Array >::neighbors_, STK::IRunnerUnsupervised< Array, Weights >::p_data(), STK::LocalVariance< Array >::prim(), STK::Reduct::prim_, STK::IArray2D< Derived >::resize(), STKRUNTIME_ERROR_NO_ARG, STK::LocalVariance< Array >::type_, and STK::Reduct::unknown_graph_.

Referenced by STK::LocalVariance< Array >::LocalVariance(), and STK::LocalVariance< Array >::LocalVariance().

◆ LocalVariance() [2/4]

template<class Array >
STK::LocalVariance< Array >::LocalVariance ( Array const data,
Reduct::TypeGraph  type = Reduct::distance_,
int  nbNeighbor = 1 
)

Constructor.

the TypeGraph and the number of neighbors are given by the user and cannot modified.

Parameters
datadata set to process
typetype of proximity graph to build
nbNeighbornumber of neighbors to use in the proximity graph

Definition at line 198 of file STK_LocalVariance.h.

200 : Base(data)
201 , type_(type)
203 , neighbors_(data.rows(), Range(1,nbNeighbor_))
204 , dist_(data.rows(), Range(1,nbNeighbor_), Arithmetic<Real>::max())
205{
206 // compute minimal proximity graph of the data set
207 switch (type_)
208 {
209 case Reduct::prim_:
210 prim();
211 break;
214 break;
216 STKRUNTIME_ERROR_NO_ARG(LocalVariance::LocalVariance(data, type, nbNeighbor),unknown proximity graph);
217 break;
218 };
219}

References STK::Reduct::distance_, STK::LocalVariance< Array >::LocalVariance(), STK::LocalVariance< Array >::minimalDistance(), STK::LocalVariance< Array >::nbNeighbor(), STK::LocalVariance< Array >::prim(), STK::Reduct::prim_, STKRUNTIME_ERROR_NO_ARG, STK::LocalVariance< Array >::type_, and STK::Reduct::unknown_graph_.

◆ LocalVariance() [3/4]

template<class Array >
STK::LocalVariance< Array >::LocalVariance ( LocalVariance< Array > const reducer)

copy Constructor.

Parameters
reducerthe reducer to copy

Definition at line 225 of file STK_LocalVariance.h.

226 : Base(reducer)
227 , type_(reducer.type_)
228 , nbNeighbor_(reducer.nbNeighbor_)
229 , neighbors_(reducer.neighbors_)
230 , dist_(reducer.dist_)
231{}

◆ ~LocalVariance()

template<class Array >
STK::LocalVariance< Array >::~LocalVariance ( )
virtual

Destructor.

Definition at line 237 of file STK_LocalVariance.h.

237{}

◆ LocalVariance() [4/4]

template<class Array >
STK::LocalVariance< Array >::LocalVariance ( Reduct::TypeGraph  type = Reduct::distance_,
int  nbNeighbor = 1 
)
private

Default constructor.

The TypeGraph and the number of neighbors are given by the user and are not modified.

Parameters
typetype of proximity graph to build
nbNeighbornumber of neighbors to use in the proximity graph

Member Function Documentation

◆ clone()

template<class Array >
virtual LocalVariance * STK::LocalVariance< Array >::clone ( ) const
inlinevirtual

clone pattern

Definition at line 95 of file STK_LocalVariance.h.

95{ return new LocalVariance(*this);}

◆ computeAxis()

template<class Array >
void STK::LocalVariance< Array >::computeAxis ( )
private

compute the axis using the first eigenvectors of the matrix of covariance and local covariance

Definition at line 370 of file STK_LocalVariance.h.

371{
372 // compute the number of axis, cannot be greater than the dimension of the data
373 Range range(p_data_->beginCols(), std::min(this->dim_, p_data_->sizeCols()));
374 // compute the eigenvalues decomposition of the local covariance
375 SymEigen<ArraySquareX>* decomp = new SymEigen<ArraySquareX>(localCovariance_);
376 decomp->run();
377 // compute the generalized inverse of the local covariance
378 ArraySquareX inv_local;
379 decomp->ginv(inv_local);
380 // compute the product
381 ArraySquareX prod = inv_local * covariance_;
382 // compute the eigenvalues decomposition of the product
383 decomp->setData(prod);
384 decomp->run();
385 // save axis and index values
386 axis_.resize(p_data_->cols(), range);
387 idx_values_.resize(range);
388 axis_ = decomp->rotation().col(range);
389 idx_values_ = decomp->eigenValues()[range];
390 // we can safely remove the decomposition
391 delete decomp;
392}
VectorX idx_values_
The values of the index for each axis.
Array const * p_data_
A pointer on the original data set.
ArraySquareX localCovariance_
the local covariance Array
ArraySquareX covariance_
the covariance Array
Array2DSquare< Real > ArraySquareX

References STK::IRegression< YArray, XArray, Weights >::run(), and STK::IRunnerSupervised< YArray_, XArray_, Weights_ >::setData().

◆ computeCovarianceMatrices() [1/2]

template<class Array >
void STK::LocalVariance< Array >::computeCovarianceMatrices ( )
protected

compute the covariances matrices of the data set

Definition at line 304 of file STK_LocalVariance.h.

305{
306 Stat::Multivariate<Array, Real> stats_;
307 // compute the covariance matrix
308 stats_.setData(p_data_);
309 stats_.run();
310 covariance_.move(stats_.covariance());
311 // constants
312 const Real pond = 2* nbNeighbor_ * p_data_->sizeRows();
313
314 // compute local covariance matrix
316 for (int j=p_data_->beginCols(); j<p_data_->endCols(); j++)
317 {
318 for (int k=p_data_->beginCols(); k<p_data_->endCols(); k++)
319 {
320 Real sum = 0.0;
321 for (int i=p_data_->beginRows(); i<p_data_->endRows(); i++)
322 {
323 for (int l = 1; l <= nbNeighbor_; ++l)
324 {
325 sum += ((*p_data_)(i, j) - (*p_data_)(neighbors_(i, l), j))
326 * ((*p_data_)(i, k) - (*p_data_)(neighbors_(i, l), k));
327 }
328 }
329 localCovariance_(j, k) = sum/pond;
330 }
331 }
332}
Derived & move(Derived const &T)
move T to this.
Arrays::SumOp< Lhs, Rhs >::result_type sum(Lhs const &lhs, Rhs const &rhs)
convenience function for summing two arrays
double Real
STK fundamental type of Real values.

References STK::IRegression< YArray, XArray, Weights >::run(), STK::IRunnerSupervised< YArray_, XArray_, Weights_ >::setData(), and STK::sum().

◆ computeCovarianceMatrices() [2/2]

template<class Array >
void STK::LocalVariance< Array >::computeCovarianceMatrices ( VectorX const weights)
protected

compute the weighted covariances matrices of the data set

Definition at line 336 of file STK_LocalVariance.h.

337{
338 Stat::Multivariate<Array, Real> stats_;
339 // compute the covariance matrix
340 stats_.setData(p_data_);
341 stats_.run(weights);
342 covariance_.move(stats_.covariance());
343
344 // get dimensions
345 const Real pond = 2* nbNeighbor_ * p_data_->sizeRows() ;
346 // compute weighted local covariance matrix
348 for (int i=p_data_->beginCols(); i<p_data_->endCols(); i++)
349 {
350 for (int j=p_data_->beginCols(); j<p_data_->endCols(); j++)
351 {
352 Real sum = 0.0;
353 for (int k=p_data_->beginRows(); k<p_data_->endRows(); k++)
354 {
355 for (int l = 1; l <= nbNeighbor_; ++l)
356 {
357 sum += (weights[k] * weights[neighbors_(k, l)])
358 * ((*p_data_)(k, i) - (*p_data_)(neighbors_(k, l), i))
359 * ((*p_data_)(k, j) - (*p_data_)(neighbors_(k, l), j));
360 }
361 }
362 localCovariance_(i, j) = sum / pond;
363 }
364 }
365}

References STK::IRegression< YArray, XArray, Weights >::run(), STK::IRunnerSupervised< YArray_, XArray_, Weights_ >::setData(), and STK::sum().

◆ covariance()

template<class Array >
ArraySquareX const & STK::LocalVariance< Array >::covariance ( ) const
inline
Returns
the covariance matrix of the data set

Definition at line 103 of file STK_LocalVariance.h.

103{ return covariance_;}

References STK::LocalVariance< Array >::covariance_.

◆ localCovariance()

template<class Array >
ArraySquareX const & STK::LocalVariance< Array >::localCovariance ( ) const
inline
Returns
the local covariance matrix computed using the proximity graph.

Definition at line 105 of file STK_LocalVariance.h.

105{ return localCovariance_;}

References STK::LocalVariance< Array >::localCovariance_.

◆ maximizeStep() [1/2]

template<class Array >
void STK::LocalVariance< Array >::maximizeStep ( )
protectedvirtual

Compute the axis by maximizing the ratio of the local variance on the total variance of the data set.

Implements STK::ILinearReduct< Array, VectorX >.

Definition at line 271 of file STK_LocalVariance.h.

272{
273#ifdef STK_REDUCT_DEBUG
274 if (!p_data_)
276#endif
277 // compute covariance matrices
279 // compute the axis
280 computeAxis();
281}
void computeAxis()
compute the axis using the first eigenvectors of the matrix of covariance and local covariance
void computeCovarianceMatrices()
compute the covariances matrices of the data set
virtual void maximizeStep()
Compute the axis by maximizing the ratio of the local variance on the total variance of the data set.

References STK::LocalVariance< Array >::maximizeStep(), and STKRUNTIME_ERROR_NO_ARG.

Referenced by STK::LocalVariance< Array >::maximizeStep(), and STK::LocalVariance< Array >::maximizeStep().

◆ maximizeStep() [2/2]

template<class Array >
void STK::LocalVariance< Array >::maximizeStep ( VectorX const weights)
protectedvirtual

Compute the axis by maximizing the ratio of the weighted local variance on the weighted total variance of the data set.

Parameters
weightsthe weights to use

Implements STK::ILinearReduct< Array, VectorX >.

Definition at line 288 of file STK_LocalVariance.h.

289{
290#ifdef STK_REDUCT_DEBUG
291 if (!p_data_)
293#endif
294 // compute covariance matrices
296 // compute the axis
297 computeAxis();
298}

References STK::LocalVariance< Array >::maximizeStep(), and STKRUNTIME_ERROR_NO_ARG.

◆ minimalDistance()

template<class Array >
void STK::LocalVariance< Array >::minimalDistance ( )
protected

compute the minimal distance graph

Definition at line 442 of file STK_LocalVariance.h.

443{
444 dist_.resize(p_data_->rows(), Range(1,nbNeighbor_));
445 dist_ = Arithmetic<Real>::max();
446 // get dimensions
447 const int begin_ind = p_data_->beginRows();
448 const int last_ind = p_data_->lastIdxRows();
449 // start minimal distance algorithm
450 for (int j = begin_ind; j<last_ind; j++)
451 {
452 // reference on the current point
453 RowVector curPoint(p_data_->row(j), true);
454 // update distance of the neighbors point
455 for (int i=j+1; i<=last_ind; i++)
456 {
457 // compute distance between point i and point j
458 Real d=dist(curPoint, p_data_->row(i));
459 // check if we get a better distance for the point j
460 if (dist_(i, nbNeighbor_) > d )
461 {
462 // check if we get a better distance for the point i
463 int pos = nbNeighbor_;
464 while (dist_(i, pos) > d && pos-- > 1) {}
465 pos++;
466 // shift values
467 for (int k= nbNeighbor_ -1; k>= pos; k--)
468 {
469 dist_(i, k+1) = dist_(i, k);
470 neighbors_(i, k+1) = neighbors_(i, k);
471 }
472 // set minimal distance in place
473 dist_(i, pos) = d;
474 neighbors_(i, pos) = j;
475 }
476 // check if we get a better distance for the point j
477 if (dist_(j, nbNeighbor_) > d )
478 {
479 // insertion sorting algorihtm
480 int pos = nbNeighbor_;
481 while (dist_(j, pos) > d && pos-- > 1) {}
482 pos++;
483 // shift valuesconst
484 for (int k= nbNeighbor_ -1; k>= pos; k--)
485 {
486 dist_(j, k+1) = dist_(j, k);
487 neighbors_(j, k+1) = neighbors_(j, k);
488 }
489 // set minimal distance in place
490 dist_(j, pos) = d;
491 neighbors_(j, pos) = i;
492 }
493 }
494 }
495}
hidden::Traits< Array >::Row RowVector
Real dist(ExprBase< Container1D1 > const &x, ExprBase< Container1D2 > const &y)
Compute the distance between two vectors.

References STK::dist().

Referenced by STK::LocalVariance< Array >::LocalVariance(), and STK::LocalVariance< Array >::LocalVariance().

◆ nbNeighbor()

template<class Array >
int STK::LocalVariance< Array >::nbNeighbor ( ) const
inline
Returns
the number of neighbors used in the local covariance.

Definition at line 99 of file STK_LocalVariance.h.

99{ return nbNeighbor_;}

References STK::LocalVariance< Array >::nbNeighbor_.

Referenced by STK::LocalVariance< Array >::LocalVariance(), and STK::LocalVariance< Array >::LocalVariance().

◆ pred()

template<class Array >
ArrayXXi const & STK::LocalVariance< Array >::pred ( ) const
inline
Returns
an array with the index of the neighbors of an individual

Definition at line 101 of file STK_LocalVariance.h.

101{ return neighbors_;}

References STK::LocalVariance< Array >::neighbors_.

◆ prim()

template<class Array >
void STK::LocalVariance< Array >::prim ( )
protected

compute the minimal spanning tree

Definition at line 395 of file STK_LocalVariance.h.

396{
397 // get dimensions
398 const int begin_ind = p_data_->beginRows();
399 const int last_ind = p_data_->lastIdxRows();
400 /* value vector : store the minimal value. */
401 VectorX value(p_data_->rows(), Arithmetic<Real>::max());
402 /* position of the points */
403 Array1D<int> ipos(p_data_->rows());
404 // Initialization the position array
405 for (int i=begin_ind; i<=last_ind; i++) ipos[i] = i;
406
407 // start Prim algorithm
408 //Initialization of the root
409 value[begin_ind] = 0.0; // the root have value 0.0
410 neighbors_(begin_ind, 1) = begin_ind; // and have itself as predecessor
411 int imin = begin_ind; // the index of the current minimal value
412 Real kmin = 0.0; // the current minimal value
413 // begin iterations
414 for (int iter = last_ind; iter>=begin_ind; iter--)
415 {
416 // put the minimal key at the end of the array key_
417 value.swap(imin, iter); // put the minimal value to the end
418 ipos.swap(imin, iter); // update the position of current minimal point
419 // Update the value for the neighbors points and find minimal value
420 imin = begin_ind;
421 kmin = value[begin_ind];
422 // reference on the current point
423 int icur = ipos[iter];
424 RowVector P(p_data_->row(icur), true);
425 // update distance of the neighbors point
426 for (int i=begin_ind; i<iter; i++)
427 {
428 // check if we have a better distance for the neighbors
429 Real d=dist(P, p_data_->row(ipos[i]));
430 if (d < value[i])
431 {
432 value[i] = d;
433 neighbors_(ipos[i], 1) = icur;
434 }
435 // minimal key
436 if (kmin>value[i]) { imin=i; kmin = value[i];}
437 }
438 }
439}
Array2DVector< Real > VectorX

References STK::dist(), and STK::IArray2D< Derived >::swap().

Referenced by STK::LocalVariance< Array >::LocalVariance(), and STK::LocalVariance< Array >::LocalVariance().

◆ update()

template<class Array >
void STK::LocalVariance< Array >::update ( )
protectedvirtual

Compute the proximity graph if the data set is modified.

Reimplemented from STK::IRunnerBase.

Definition at line 243 of file STK_LocalVariance.h.

244{
245#ifdef STK_REDUCT_DEBUG
246 if (!p_data_)
248#endif
249 // update dimensions of the containers for the proximity graph
251 // compute minimal proximity graph of the data set
252 switch (type_)
253 {
254 case Reduct::prim_:
255 prim();
256 break;
259 break;
261 STKRUNTIME_ERROR_NO_ARG(LocalVariance::update,unknown proximity graph);
262 break;
263 };
264}
virtual void update()
Compute the proximity graph if the data set is modified.

References STK::Reduct::distance_, STK::Reduct::prim_, STKRUNTIME_ERROR_NO_ARG, STK::Reduct::unknown_graph_, and STK::LocalVariance< Array >::update().

Referenced by STK::LocalVariance< Array >::update().

Member Data Documentation

◆ covariance_

template<class Array >
ArraySquareX STK::LocalVariance< Array >::covariance_
protected

the covariance Array

Definition at line 126 of file STK_LocalVariance.h.

Referenced by STK::LocalVariance< Array >::covariance().

◆ dist_

template<class Array >
ArrayXX STK::LocalVariance< Array >::dist_
protected

distance matrix : store the minimal distance to the neighbors

Definition at line 121 of file STK_LocalVariance.h.

Referenced by STK::LocalVariance< Array >::LocalVariance().

◆ localCovariance_

template<class Array >
ArraySquareX STK::LocalVariance< Array >::localCovariance_
protected

the local covariance Array

Definition at line 124 of file STK_LocalVariance.h.

Referenced by STK::LocalVariance< Array >::localCovariance().

◆ nbNeighbor_

template<class Array >
int STK::LocalVariance< Array >::nbNeighbor_
protected

number of neighbors

Definition at line 116 of file STK_LocalVariance.h.

Referenced by STK::LocalVariance< Array >::LocalVariance(), and STK::LocalVariance< Array >::nbNeighbor().

◆ neighbors_

template<class Array >
ArrayXXi STK::LocalVariance< Array >::neighbors_
protected

Predecessor array : store the spanning tree or the minimal distance to the neighbors.

Definition at line 119 of file STK_LocalVariance.h.

Referenced by STK::LocalVariance< Array >::LocalVariance(), and STK::LocalVariance< Array >::pred().

◆ type_

template<class Array >
Reduct::TypeGraph STK::LocalVariance< Array >::type_
protected

The documentation for this class was generated from the following file: