STK++ 0.9.13
STK::IMixtureComposer Class Referenceabstract

Base class for Mixture (composed) model. More...

#include <STK_IMixtureComposer.h>

Inheritance diagram for STK::IMixtureComposer:
Inheritance graph

Public Member Functions

virtual ~IMixtureComposer ()
 destructor
 
Clust::modelState state () const
 
void setState (Clust::modelState state)
 set the state of the model : should be used by any strategy
 
virtual IMixtureComposercreate () const =0
 create pattern
 
virtual IMixtureComposerclone () const =0
 clone pattern
 
virtual void randomInit ()=0
 initialize randomly the parameters of the components of the model
 
virtual void paramUpdateStep ()=0
 Compute the proportions and the model parameters given the current tik mixture parameters.
 
virtual int cStep ()
 Replace tik by zik.
 
virtual int sStep ()
 Simulate zi accordingly to tik and replace tik by zik by calling cStep().
 
virtual Real eStep ()
 compute the zi, the lnLikelihood of the current estimates and the next value of the tik.
 
virtual void mapStep ()
 Compute zi using the Map estimate.
 
virtual void finalizeStep ()
 Finalize the estimation of the model.
 
virtual void pStep ()
 Compute proportions using the ML estimates, default implementation.
 
virtual void initializeStep ()
 Initialize the model before its first use.
 
void randomClassInit ()
 Initialize randomly the labels zi of the model.
 
void randomFuzzyInit ()
 Initialize randomly the posterior probabilities tik of the model, then compute the zi values with mapStep, compute the initial parameter values using paramUpdateStep, and compute the tik using eStep().
 
- Public Member Functions inherited from STK::IMixtureStatModel
virtual ~IMixtureStatModel ()
 destructor
 
int nbCluster () const
 
CPointX constpk () const
 
CArrayXX consttik () const
 
CPointX consttk () const
 
CVectorXi constzi () const
 
std::vector< IMixture * > constv_mixtures () const
 
Real computeLnLikelihood (int i) const
 
Real computeLikelihood (int i) const
 
Real computeLnLikelihood () const
 
Real computeICL () const
 
IMixturegetMixture (String const &idData) const
 Utility lookup function allowing to find a Mixture from its idData.
 
void registerMixture (IMixture *p_mixture)
 register a mixture to the composer.
 
void releaseMixture (String const &idData)
 release a mixture from the composer.
 
int computeNbFreeParameters () const
 compute the number of free parameters of the model.
 
int computeNbMissingValues () const
 compute the missing values of the model.
 
virtual Real lnComponentProbability (int i, int k) const =0
 
virtual void imputationStep ()
 Impute the missing values.
 
virtual void samplingStep ()
 Simulation of all the latent variables and/or missing data excluding class labels.
 
virtual void setParametersStep ()
 Utility method allowing to signal to a mixture to set its parameters.
 
virtual void storeIntermediateResults (int iteration)
 This step can be used to signal to the mixtures that they must store results.
 
virtual void releaseIntermediateResults ()
 This step can be used to signal to the mixtures that they must release the stored results.
 
virtual void writeParameters (ostream &os) const
 write the parameters of the model in the stream os.
 
template<class Array >
void setMixtureParameters (Array const &tik)
 set the mixture parameters using an array of posterior probabilities.
 
template<class Array , class RowVector >
void setMixtureParameters (Array const &tik, RowVector const &pk)
 set the mixture parameters giving the posterior probabilities and the proportions.
 
template<class RowVector >
void setProportions (RowVector const &pk)
 Set proportions of each classes.
 
template<class Manager , class Parameters >
void setParameters (IMixtureManager< Manager > const &manager, String const &idData, Parameters const &param)
 Utility method allowing to set the parameters to a specific mixture.
 
template<class Manager >
void createMixture (IMixtureManager< Manager > &manager)
 Utility method allowing to create all the mixtures registered in the data handler of a mixture manager and to register them.
 
template<class Manager >
IMixturecreateMixture (IMixtureManager< Manager > &manager, String const &idData)
 Utility method allowing to create a mixture with a given data set and register it.
 
template<class Manager >
void removeMixture (IMixtureManager< Manager > &manager, String const &idData)
 Utility method allowing to release completely a mixture with its data set.
 
template<class Manager , class Parameters >
void getParameters (IMixtureManager< Manager > const &manager, String const &idData, Parameters &param) const
 Utility method allowing to get the parameters of a specific mixture.
 
template<class Manager , class MissingValues >
void getMissingValues (IMixtureManager< Manager > const &manager, String const &idData, MissingValues &missing) const
 Utility method allowing to get the missing values of a specific mixture.
 
template<class DataHandler >
void createMixture (IMixtureManager< DataHandler > &manager)
 
- Public Member Functions inherited from STK::IStatModelBase
int nbSample () const
 
Real lnNbSample () const
 
int nbVariable () const
 
Real lnLikelihood () const
 
Real likelihood () const
 
int nbFreeParameter () const
 
Real computeBIC () const
 
Real computeAIC () const
 
Real computeML () const
 

Protected Member Functions

 IMixtureComposer (int nbSample, int nbCluster)
 Constructor.
 
 IMixtureComposer (IMixtureComposer const &model)
 copy constructor.
 
void sStep (int i)
 Simulate zi accordingly to tik.
 
void cStep (int i)
 Replace tik by zik.
 
Real eStep (int i)
 compute one zi and the next value of the tik for i fixed
 
void mapStep (int i)
 Compute zi using the Map estimate for i fixed.
 
virtual void initializeMixtureParameters ()
 Create the mixture model parameters pk_ and tik_.
 
virtual int randomTik ()
 generate random tik_
 
virtual int randomZi ()
 generate random zi_
 
- Protected Member Functions inherited from STK::IMixtureStatModel
 IMixtureStatModel (int nbSample, int nbCluster)
 Constructor.
 
 IMixtureStatModel (IMixtureStatModel const &model)
 copy constructor.
 
void setNbCluster (int nbCluster)
 set the number of cluster of the model
 
- Protected Member Functions inherited from STK::IStatModelBase
 IStatModelBase ()
 Default constructor.
 
 IStatModelBase (int nbSample)
 Constructor with specified dimension.
 
 IStatModelBase (int nbSample, int nbVariable)
 Constructor with specified dimension.
 
 IStatModelBase (IStatModelBase const &model)
 Copy constructor.
 
 ~IStatModelBase ()
 destructor
 
void setNbFreeParameter (int const &nbFreeParameter)
 set the number of free parameters of the model
 
void setNbSample (int const &nbSample)
 set the number of samples of the model
 
void setNbVariable (int const &nbVariable)
 set the number of variables of the model
 
void setLnLikelihood (Real const &lnLikelihood)
 set the log-likelihood of the model
 
void initialize (int nbSample, int nbVariable)
 set the dimensions of the parameters of the model
 

Private Attributes

Clust::modelState state_
 state of the model
 
CPointX lnComp_
 Auxiliary array used in the eStep.
 

Additional Inherited Members

- Public Types inherited from STK::IMixtureStatModel
typedef std::vector< IMixture * >::const_iterator ConstMixtIterator
 
typedef std::vector< IMixture * >::iterator MixtIterator
 
- Protected Attributes inherited from STK::IMixtureStatModel
int nbCluster_
 number of cluster.
 
CPointX pk_
 The proportions of each mixtures.
 
CArrayXX tik_
 The tik probabilities.
 
CPointX tk_
 The sum of the columns of tik_.
 
CVectorXi zi_
 The zi class label.
 
std::vector< IMixture * > v_mixtures_
 vector of pointers to the mixtures components
 

Detailed Description

Base class for Mixture (composed) model.

In this interface we assume there is an underline generative model that will be estimated using either an EM, SEM or CEM algorithm. All mixture parameters are created using the method

The MultidimRegression class allows to regress a multidimensional output variable among a multivariat...

in the constructor. They can be accessed from the mixtures using constant accessors.

The pure virtual function to implement in derived class are

virtual IMixtureComposer* create() const = 0;
virtual IMixtureComposer* clone() const = 0;
virtual void randomInit() = 0;
virtual void paramUpdateStep() = 0;
virtual int computeNbFreeParameters() const = 0;
Base class for Mixture (composed) model.
virtual IMixtureComposer * clone() const =0
clone pattern
virtual IMixtureComposer * create() const =0
create pattern
virtual void randomInit()=0
initialize randomly the parameters of the components of the model
virtual void paramUpdateStep()=0
Compute the proportions and the model parameters given the current tik mixture parameters.
int computeNbFreeParameters() const
compute the number of free parameters of the model.

The public virtual function that can be re-implemented in derived class for a specific behavior are:

virtual void initializeStep();
virtual void randomClassInit();
virtual void randomFuzzyInit();
virtual int cStep();
virtual int sStep();
virtual Real eStep();
virtual void mapStep();
virtual void finalizeStep();
virtual void pStep();
virtual int cStep()
Replace tik by zik.
virtual void initializeStep()
Initialize the model before its first use.
void randomFuzzyInit()
Initialize randomly the posterior probabilities tik of the model, then compute the zi values with map...
virtual void pStep()
Compute proportions using the ML estimates, default implementation.
virtual Real eStep()
compute the zi, the lnLikelihood of the current estimates and the next value of the tik.
virtual void mapStep()
Compute zi using the Map estimate.
void randomClassInit()
Initialize randomly the labels zi of the model.
virtual int sStep()
Simulate zi accordingly to tik and replace tik by zik by calling cStep().
virtual void finalizeStep()
Finalize the estimation of the model.
double Real
STK fundamental type of Real values.
See also
IMixture
Note
the virtual method IMixtureComposer::initializeStep is called in all the initialization method. Don't forget to called it in the randomInit implementation.

Definition at line 85 of file STK_IMixtureComposer.h.

Constructor & Destructor Documentation

◆ IMixtureComposer() [1/2]

STK::IMixtureComposer::IMixtureComposer ( int  nbSample,
int  nbCluster 
)
protected

Constructor.

Parameters
nbCluster,nbSamplenumber of clusters and samples

Definition at line 49 of file STK_IMixtureComposer.cpp.

52#ifndef _OPENMP
54#endif
55{}
CPointX lnComp_
Auxiliary array used in the eStep.
Clust::modelState state_
state of the model
IMixtureStatModel(int nbSample, int nbCluster)
Constructor.
@ modelCreated_
the model has been created but is not initialized

◆ IMixtureComposer() [2/2]

STK::IMixtureComposer::IMixtureComposer ( IMixtureComposer const model)
protected

copy constructor.

Parameters
modelthe model to clone

Definition at line 58 of file STK_IMixtureComposer.cpp.

59 : IMixtureStatModel(model)
60 , state_(model.state_)
61#ifndef _OPENMP
62 , lnComp_(model.lnComp_)
63#endif
64{}

◆ ~IMixtureComposer()

STK::IMixtureComposer::~IMixtureComposer ( )
virtual

destructor

Definition at line 66 of file STK_IMixtureComposer.cpp.

66{}

Member Function Documentation

◆ clone()

◆ create()

◆ cStep() [1/2]

int STK::IMixtureComposer::cStep ( )
virtual

Replace tik by zik.

Returns
the minimal value of individuals in a class

Reimplemented in STK::MixtureSemiLearner.

Definition at line 149 of file STK_IMixtureComposer.cpp.

150{
151#ifdef STK_MIXTURE_VERY_VERBOSE
152 stk_cout << _T("Entering IMixtureComposer::cStep()\n");
153#endif
154 for (int i=tik_.beginRows(); i < tik_.endRows(); i++)
155 { cStep(i);}
156 // count the number of individuals in each class
158#ifdef STK_MIXTURE_VERY_VERBOSE
159 stk_cout << _T("IMixtureComposer::cStep() tk_ = ") << tk_;
160 stk_cout << _T("IMixtureComposer::cStep() done\n");
161#endif
162 return tk_.minElt();
163}
#define stk_cout
Standard stk output stream.
#define _T(x)
Let x unmodified.
Type const minElt(int &row, int &col) const
CArrayXX tik_
The tik probabilities.
CPointX tk_
The sum of the columns of tik_.
hidden::FunctorTraits< Derived, SumOp >::Row sum(Derived const &A)
Compute the sum of A.

References _T, cStep(), STK::ExprBase< Derived >::minElt(), stk_cout, STK::Stat::sum(), STK::IMixtureStatModel::tik_, and STK::IMixtureStatModel::tk_.

Referenced by cStep(), STK::MixtureSemiLearner::cStep(), randomZi(), STK::CEMAlgo::run(), and sStep().

◆ cStep() [2/2]

void STK::IMixtureComposer::cStep ( int  i)
protected

Replace tik by zik.

Parameters
iindex of the the individual

Definition at line 219 of file STK_IMixtureComposer.cpp.

220{ tik_.row(i) = 0.; tik_.elt(i, zi_[i]) = 1.;}
hidden::CSlice< Derived, 1, sizeCols_ >::Result row(int i) const
implement the row operator using a reference on the row of the allocator
CVectorXi zi_
The zi class label.

References STK::ICArray< Derived >::row(), STK::IMixtureStatModel::tik_, and STK::IMixtureStatModel::zi_.

◆ eStep() [1/2]

Real STK::IMixtureComposer::eStep ( )
virtual

compute the zi, the lnLikelihood of the current estimates and the next value of the tik.

Returns
the minimal value of tk

Reimplemented in STK::MixtureSemiLearner.

Definition at line 180 of file STK_IMixtureComposer.cpp.

181{
182#ifdef STK_MIXTURE_VERY_VERBOSE
183 stk_cout << _T("Entering IMixtureComposer::eStep()\n");
184#endif
185 Real sum = 0.;
186 int i;
187#ifdef _OPENMP
188#pragma omp parallel for reduction (+:sum)
189#endif
190 for (i = tik_.beginRows(); i < tik_.endRows(); ++i) { sum += eStep(i);}
191 // update ln-likelihood
193 // compute proportions
195#ifdef STK_MIXTURE_VERY_VERBOSE
196 stk_cout << _T("IMixtureComposer::eStep() lnLikelihood =") << sum << _T("\n");
197 stk_cout << _T("IMixtureComposer::eStep() done\n");
198 if (!isFinite(sum))
199 { stk_cout << _T("tik_ =\n") << tik_.transpose();}
200#endif
201 return tk_.minElt();
202}
TransposeOperator< Derived > const transpose() const
void setLnLikelihood(Real const &lnLikelihood)
set the log-likelihood of the model
bool isFinite(Type const &x)
utility method allowing to know if a value is a finite value
Arrays::SumOp< Lhs, Rhs >::result_type sum(Lhs const &lhs, Rhs const &rhs)
convenience function for summing two arrays
hidden::FunctorTraits< Derived, SumOp >::Row sumByCol(Derived const &A)

References _T, eStep(), STK::isFinite(), STK::ExprBase< Derived >::minElt(), STK::IStatModelBase::setLnLikelihood(), stk_cout, STK::sum(), STK::Stat::sumByCol(), STK::IMixtureStatModel::tik_, STK::IMixtureStatModel::tk_, and STK::ArrayBase< Derived >::transpose().

Referenced by STK::IMixtureAlgoPredict::burnStep(), eStep(), STK::MixtureSemiLearner::eStep(), STK::IMixtureAlgoPredict::predictBayesClassifier(), randomClassInit(), randomFuzzyInit(), STK::MixtureComposer::randomInit(), STK::EMAlgo::run(), STK::CEMAlgo::run(), STK::SEMAlgo::run(), STK::SemiSEMAlgo::run(), STK::EMPredict::run(), and STK::SemiSEMPredict::run().

◆ eStep() [2/2]

Real STK::IMixtureComposer::eStep ( int  i)
protected

compute one zi and the next value of the tik for i fixed

Parameters
ithe individual
Returns
the contribution of the individual i to the log-likelihood

Definition at line 223 of file STK_IMixtureComposer.cpp.

224{
225#ifdef _OPENMP
226 CPointX lnComp_(tik_.cols());
227#endif
228 // get maximal element of ln(x_i,\theta_k) + ln(p_k)
229 for (int k=tik_.beginCols(); k< tik_.endCols(); k++)
230 { lnComp_[k] = std::log(pk_[k])+lnComponentProbability(i,k);}
231 if (lnComp_.isInfinite().any())
232 {
233#ifdef STK_MIXTURE_VERY_VERBOSE
234 stk_cout << _T("IMixtureComposer::eStep(") << i << _T(")\n");
235 stk_cout << _T("pk_ =") << pk_;
236 stk_cout << _T("tk_ =") << tk_;
237 stk_cout << _T("lnComp_ =") << lnComp_;
238 for (int k=tik_.beginCols(); k< tik_.endCols(); k++)
239 {
240 Real val = lnComponentProbability(i,k);
241 stk_cout << _T("lnComponentProbability(") << i << _T(",")<< k << _T(")=") << val << _T("\n");
242 }
244 stk_cout << _T("tik_ =\n") << tik_.transpose();
245#endif
246 throw(Clust::eStepFail_);
247 }
248 int kmax;
249 Real max = lnComp_.maxElt(kmax);
250 // set zi_
251 zi_[i] = kmax;
252 // return max + sum_k p_k exp{lnCom_k - lnComp_kmax}
253 Real sum = (tik_.row(i) = (lnComp_ - max).exp()).sum();
254 tik_.row(i) /= sum;
255 return max + std::log( sum );
256}
Type const maxElt(int &row, int &col) const
UnaryOperator< IsInfiniteOp< Type >, Derived > isInfinite() const
virtual void writeParameters(ostream &os) const
write the parameters of the model in the stream os.
virtual Real lnComponentProbability(int i, int k) const =0
CPointX pk_
The proportions of each mixtures.
hidden::SliceVisitorSelector< Derived, hidden::MaxVisitor, Arrays::by_col_ >::type_result max(Derived const &A)
If A is a row-vector or a column-vector then the function will return the usual maximal value of the ...
CArrayPoint< Real, UnknownSize, Arrays::by_col_ > CPointX

References _T, STK::Clust::eStepFail_, STK::ExprBase< Derived >::isInfinite(), lnComp_, STK::IMixtureStatModel::lnComponentProbability(), STK::max(), STK::ExprBase< Derived >::maxElt(), STK::IMixtureStatModel::pk_, STK::ICArray< Derived >::row(), stk_cout, STK::sum(), STK::IMixtureStatModel::tik_, STK::IMixtureStatModel::tk_, STK::ArrayBase< Derived >::transpose(), STK::IMixtureStatModel::writeParameters(), and STK::IMixtureStatModel::zi_.

◆ finalizeStep()

void STK::IMixtureComposer::finalizeStep ( )
virtual

Finalize the estimation of the model.

Compute obtained lnLikelihood and set state to finalized.

Reimplemented from STK::IMixtureStatModel.

Reimplemented in STK::MixtureComposer.

Definition at line 281 of file STK_IMixtureComposer.cpp.

282{
285}
void setState(Clust::modelState state)
set the state of the model : should be used by any strategy
@ modelFinalized_
the model is finalized

References STK::IMixtureStatModel::computeLnLikelihood(), STK::Clust::modelFinalized_, STK::IStatModelBase::setLnLikelihood(), and setState().

Referenced by STK::IMixtureAlgoPredict::predictBayesClassifier(), STK::EMPredict::run(), STK::SemiSEMPredict::run(), and STK::StrategyFacade::run().

◆ initializeMixtureParameters()

void STK::IMixtureComposer::initializeMixtureParameters ( )
protectedvirtual

Create the mixture model parameters pk_ and tik_.

Default implementation is to set pk_ and tik_ arrays to 1/K value.

Reimplemented in STK::MixtureSemiLearner.

Definition at line 290 of file STK_IMixtureComposer.cpp.

291{
292#ifdef STK_MIXTURE_VERBOSE
293 stk_cout << _T("Entering IMixtureComposer::initializeMixtureParameters\n");
294#endif
295 pk_ = 1./nbCluster_;
296 tik_ = 1./nbCluster_;
297 tk_ = sumByCol(tik_);
298}
int nbCluster_
number of cluster.
hidden::SliceVisitorSelector< Derived, hidden::SumVisitor, Arrays::by_col_ >::type_result sumByCol(Derived const &A)

References _T, STK::IMixtureStatModel::nbCluster_, STK::IMixtureStatModel::pk_, stk_cout, STK::sumByCol(), STK::IMixtureStatModel::tik_, and STK::IMixtureStatModel::tk_.

Referenced by initializeStep().

◆ initializeStep()

void STK::IMixtureComposer::initializeStep ( )
virtual

Initialize the model before its first use.

Initialize the values of the mixture parameters pk_ and @© tik_ using virtual method initializeMixtureParameters() and compute tk_ and zk_ using the virtual methods mapStep() and pStep().

Reimplemented from STK::IMixtureStatModel.

Definition at line 74 of file STK_IMixtureComposer.cpp.

75{
76#ifdef STK_MIXTURE_VERBOSE
77 stk_cout << _T("Entering IMixtureComposer::initializeStep\n");
78#endif
79 // (re)initialize the mixture parameters tik and pk. (virtual method)
81 // compute nk
83 // (re) initialize mixtures
85 // compute zi (virtual method)
86 mapStep();
87 // compute proportions pk_ (virtual method)
88 pStep();
89 // (re)initialize the likelihood
91 // compute the number of free parameters
93 // update state
95#ifdef STK_MIXTURE_VERBOSE
96 stk_cout << _T("IMixtureComposer::initializeStep done\n");
97#endif
98#ifdef STK_MIXTURE_VERY_VERBOSE
99 stk_cout << _T("Model parameters\n");
101#endif
102}
virtual void initializeMixtureParameters()
Create the mixture model parameters pk_ and tik_.
virtual void initializeStep()
Initialize the model before at its first use.
void setNbFreeParameter(int const &nbFreeParameter)
set the number of free parameters of the model
@ modelInitialized_
the model is initialized and its parameters are initialized to default values

References _T, STK::IMixtureStatModel::computeLnLikelihood(), STK::IMixtureStatModel::computeNbFreeParameters(), initializeMixtureParameters(), STK::IMixtureStatModel::initializeStep(), mapStep(), STK::Clust::modelInitialized_, pStep(), STK::IStatModelBase::setLnLikelihood(), STK::IStatModelBase::setNbFreeParameter(), setState(), stk_cout, STK::Stat::sumByCol(), STK::IMixtureStatModel::tik_, STK::IMixtureStatModel::tk_, and STK::IMixtureStatModel::writeParameters().

Referenced by STK::IMixtureAlgoPredict::burnStep(), STK::MixtureComposer::createComposer(), STK::IMixtureAlgoPredict::predictBayesClassifier(), randomClassInit(), randomFuzzyInit(), STK::MixtureComposer::randomInit(), STK::EMPredict::run(), STK::SemiSEMPredict::run(), STK::RandomInit::run(), STK::ClassInit::run(), and STK::FuzzyInit::run().

◆ mapStep() [1/2]

void STK::IMixtureComposer::mapStep ( )
virtual

Compute zi using the Map estimate.

Reimplemented in STK::MixtureSemiLearner.

Definition at line 205 of file STK_IMixtureComposer.cpp.

206{
207 for (int i = zi_.begin(); i< zi_.end(); ++i)
208 { mapStep(i); }
209}

References mapStep(), and STK::IMixtureStatModel::zi_.

Referenced by initializeStep(), mapStep(), STK::MixtureSemiLearner::mapStep(), STK::IMixtureAlgoPredict::predictBayesClassifier(), randomFuzzyInit(), STK::MixtureComposer::randomInit(), STK::EMPredict::run(), and STK::SemiSEMPredict::run().

◆ mapStep() [2/2]

void STK::IMixtureComposer::mapStep ( int  i)
protected

Compute zi using the Map estimate for i fixed.

Definition at line 259 of file STK_IMixtureComposer.cpp.

260{
261 int k;
262 tik_.row(i).maxElt(k);
263 zi_[i] = k;
264}

References STK::ICArray< Derived >::row(), STK::IMixtureStatModel::tik_, and STK::IMixtureStatModel::zi_.

◆ paramUpdateStep()

void STK::IMixtureComposer::paramUpdateStep ( )
pure virtual

Compute the proportions and the model parameters given the current tik mixture parameters.

Implemented in STK::MixtureComposer.

Definition at line 269 of file STK_IMixtureComposer.cpp.

270{ pStep();
271 /* implement specific parameters estimation in concrete class. */
272}

References pStep().

Referenced by randomClassInit(), randomFuzzyInit(), STK::EMAlgo::run(), STK::CEMAlgo::run(), STK::SEMAlgo::run(), and STK::SemiSEMAlgo::run().

◆ pStep()

void STK::IMixtureComposer::pStep ( )
virtual

Compute proportions using the ML estimates, default implementation.

Set as virtual in case we impose fixed proportions in derived Composer.

Reimplemented in STK::MixtureComposerFixedProp, and STK::MixtureSemiLearnerFixedProp.

Definition at line 275 of file STK_IMixtureComposer.cpp.

hidden::FunctorTraits< Derived, MeanOp >::Row meanByCol(Derived const &A)

References STK::Stat::meanByCol(), STK::IMixtureStatModel::pk_, and STK::IMixtureStatModel::tik_.

Referenced by initializeStep(), paramUpdateStep(), STK::EMAlgo::run(), STK::CEMAlgo::run(), STK::SEMAlgo::run(), and STK::SemiSEMAlgo::run().

◆ randomClassInit()

void STK::IMixtureComposer::randomClassInit ( )

Initialize randomly the labels zi of the model.

Initialize the model parameters using initializeStep() if it has not been already called. Simulate the zi, compute tik using eStep(), update the parameters using paramUpdateStep() and terminate using eStep().

Definition at line 105 of file STK_IMixtureComposer.cpp.

106{
107#ifdef STK_MIXTURE_VERY_VERBOSE
108 stk_cout << _T("Entering IMixtureComposer::randomClassInit()\n");
109#endif
110 // initialize mixture model if necessary
112 // generate random zi, compute tik and nk
114 // update parameters
116 // compute eStep()
117 eStep();
118 // model initialized
120#ifdef STK_MIXTURE_VERY_VERBOSE
121 stk_cout << _T("IMixtureComposer::randomClassInit() done\n");
122#endif
123}
virtual int randomZi()
generate random zi_
Clust::modelState state() const
@ modelParamInitialized_
The parameters of the model have been initialized.

References _T, eStep(), initializeStep(), STK::Clust::modelInitialized_, STK::Clust::modelParamInitialized_, paramUpdateStep(), STK::Clust::randomClassInitFail_, randomZi(), setState(), state(), and stk_cout.

Referenced by STK::ClassInit::run().

◆ randomFuzzyInit()

void STK::IMixtureComposer::randomFuzzyInit ( )

Initialize randomly the posterior probabilities tik of the model, then compute the zi values with mapStep, compute the initial parameter values using paramUpdateStep, and compute the tik using eStep().

Definition at line 126 of file STK_IMixtureComposer.cpp.

127{
128#ifdef STK_MIXTURE_VERBOSE
129 stk_cout << _T("Entering IMixtureComposer::randomFuzzyInit(). state= ") << state() << _T("\n");
130#endif
131 // initialize mixture model if necessary
133 // create random tik and compute nk
135 // compute zi
136 mapStep();
137 // update paramters values
139 // eStep
140 eStep();
141 // model intialized
143#ifdef STK_MIXTURE_VERBOSE
144 stk_cout << _T("IMixtureComposer::randomFuzzyInit() done\n");
145#endif
146}
virtual int randomTik()
generate random tik_

References _T, eStep(), initializeStep(), mapStep(), STK::Clust::modelInitialized_, STK::Clust::modelParamInitialized_, paramUpdateStep(), STK::Clust::randomFuzzyInitFail_, randomTik(), setState(), state(), and stk_cout.

Referenced by STK::FuzzyInit::run(), STK::FullStrategy::run(), STK::SimpleStrategy::run(), and STK::XemStrategy::run().

◆ randomInit()

virtual void STK::IMixtureComposer::randomInit ( )
pure virtual

initialize randomly the parameters of the components of the model

Implements STK::IMixtureStatModel.

Implemented in STK::MixtureComposer.

Referenced by STK::RandomInit::run().

◆ randomTik()

int STK::IMixtureComposer::randomTik ( )
protectedvirtual

generate random tik_

Reimplemented in STK::MixtureSemiLearner.

Definition at line 301 of file STK_IMixtureComposer.cpp.

302{
303 tk_ = 0.;
304 tik_.randUnif();
305 for (int i = tik_.beginRows(); i < tik_.endRows(); ++i)
306 {
307 // create a reference on the i-th row
308 CPointX tikRowi(tik_.row(i), true);
309 tikRowi = tikRowi * pk_;
310 tikRowi /= tikRowi.sum();
311 tk_ += tikRowi;
312 }
313 return tk_.minElt();
314}
Derived & randUnif()
set random values to this using a uniform law.
Type const sum() const

References STK::ExprBase< Derived >::minElt(), STK::IMixtureStatModel::pk_, STK::ArrayBase< Derived >::randUnif(), STK::ICArray< Derived >::row(), STK::IMixtureStatModel::tik_, and STK::IMixtureStatModel::tk_.

Referenced by randomFuzzyInit(), and STK::MixtureComposer::randomInit().

◆ randomZi()

int STK::IMixtureComposer::randomZi ( )
protectedvirtual

generate random zi_

Reimplemented in STK::MixtureSemiLearner.

Definition at line 317 of file STK_IMixtureComposer.cpp.

318{
319 Law::Categorical law(pk_);
320 zi_.rand(law);
321 return cStep();
322}
Derived & rand(Law::IUnivLaw< Type > const &law)
set random values to this using a distribution law given by the user.

References cStep(), STK::IMixtureStatModel::pk_, STK::ArrayBase< Derived >::rand(), and STK::IMixtureStatModel::zi_.

Referenced by randomClassInit().

◆ setState()

void STK::IMixtureComposer::setState ( Clust::modelState  state)
inline

set the state of the model : should be used by any strategy

Definition at line 105 of file STK_IMixtureComposer.h.

105{ state_ = state;}

References state(), and state_.

Referenced by finalizeStep(), initializeStep(), randomClassInit(), randomFuzzyInit(), STK::MixtureComposer::randomInit(), and STK::FullStrategy::run().

◆ sStep() [1/2]

int STK::IMixtureComposer::sStep ( )
virtual

Simulate zi accordingly to tik and replace tik by zik by calling cStep().

Returns
the minimal value of individuals in a class

Reimplemented in STK::MixtureSemiLearner.

Definition at line 167 of file STK_IMixtureComposer.cpp.

168{
169#ifdef STK_MIXTURE_VERY_VERBOSE
170 stk_cout << _T("Entering IMixtureComposer::sStep()\n");
171#endif
172 // simulate zi
173 for (int i = zi_.begin(); i< zi_.end(); ++i) { sStep(i);}
174#ifdef STK_MIXTURE_VERY_VERBOSE
175 stk_cout << _T("IMixtureComposer::sStep() done\n");
176#endif
177 return cStep();
178}

References _T, cStep(), sStep(), stk_cout, and STK::IMixtureStatModel::zi_.

Referenced by STK::IMixtureAlgoPredict::burnStep(), STK::SEMAlgo::run(), sStep(), and STK::MixtureSemiLearner::sStep().

◆ sStep() [2/2]

void STK::IMixtureComposer::sStep ( int  i)
protected

Simulate zi accordingly to tik.

Parameters
iindex of the the individual

Definition at line 214 of file STK_IMixtureComposer.cpp.

215{ zi_.elt(i) = Law::Categorical::rand(tik_.row(i));}
virtual int rand() const

References STK::Law::Categorical::rand(), STK::ICArray< Derived >::row(), STK::IMixtureStatModel::tik_, and STK::IMixtureStatModel::zi_.

◆ state()

Clust::modelState STK::IMixtureComposer::state ( ) const
inline
Returns
the state of the model

Definition at line 102 of file STK_IMixtureComposer.h.

102{ return state_;}

References state_.

Referenced by randomClassInit(), randomFuzzyInit(), STK::MixtureComposer::randomInit(), STK::SimpleStrategy::run(), STK::XemStrategy::run(), and setState().

Member Data Documentation

◆ lnComp_

CPointX STK::IMixtureComposer::lnComp_
private

Auxiliary array used in the eStep.

Definition at line 194 of file STK_IMixtureComposer.h.

Referenced by eStep().

◆ state_

Clust::modelState STK::IMixtureComposer::state_
private

state of the model

Definition at line 191 of file STK_IMixtureComposer.h.

Referenced by setState(), and state().


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