library: libTMVA
#include "GeneticPopulation.h"

TMVA::GeneticPopulation


class description - header file - source file - inheritance tree (.pdf)

class TMVA::GeneticPopulation

Inheritance Chart:
TMVA::GeneticPopulation

    public:
GeneticPopulation() GeneticPopulation(const TMVA::GeneticPopulation&) virtual ~GeneticPopulation() void AddFactor(Double_t from, Double_t to) void AddPopulation(TMVA::GeneticPopulation* genePool) static TClass* Class() void ClearResults() void CreatePopulation(Int_t size) Double_t GetCounterFitness() const Double_t GetFitness(Int_t index) Double_t GetFitness() multimap<Double_t,GeneticGenes>* GetGenePool() const TMVA::GeneticGenes* GetGenes() TMVA::GeneticGenes* GetGenes(Int_t index) multimap<Double_t,GeneticGenes>* GetNewGenePool() const Int_t GetPopulationSize() const vector<TMVA::GeneticRange*>& GetRanges() void GiveHint(vector<Double_t> hint, Double_t fitness = 0) virtual TClass* IsA() const void MakeChildren() void MakeMutants(Double_t probability = 30, Bool_t near = kFALSE, Double_t spread = 0.1, Bool_t mirror = kFALSE) TMVA::GeneticGenes MakeSex(TMVA::GeneticGenes male, TMVA::GeneticGenes female) void Mutate(Double_t probability = 20, Int_t startIndex = 0, Bool_t near = kFALSE, Double_t spread = 0.1, Bool_t mirror = kFALSE) TMVA::GeneticPopulation& operator=(const TMVA::GeneticPopulation&) void Print(Int_t untilIndex = -1) void Print(ostream& out, Int_t utilIndex = -1) void Reset() Bool_t SetFitness(TMVA::GeneticGenes* g, Double_t fitness, Bool_t add = kTRUE) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b) void TrimPopulation() TH1F* VariableDistribution(Int_t varNumber, Int_t bins, Int_t min, Int_t max) vector<Double_t> VariableDistribution(Int_t varNumber)

Data Members

    private:
Double_t fCounterFitness Int_t fPopulationSize multimap<Double_t,GeneticGenes>* fGenePool multimap<Double_t,GeneticGenes>* fNewGenePool vector<TMVA::GeneticRange*,allocator<TMVA::GeneticRange*> > fRanges multimap<double,TMVA::GeneticGenes,less<double>,allocator<pair<const double,TMVA::GeneticGenes> > >::iterator fCounter TRandom* fRandomGenerator

Class Description

_______________________________________________________________________

 Population definition for genetic algorithm

_______________________________________________________________________
GeneticPopulation()
 Constructor
 create a randomGenerator for this population and set a seed
 create the genePools

void CreatePopulation( Int_t size )
 create a Population of individuals with the population size given
 in the parameter
--> every coefficient gets a random value within the constraints
 provided by the user

void AddPopulation( TMVA::GeneticPopulation *strangers )
 allows to connect two populations (using the same number of coefficients
 with the same ranges)
 this allows to calculate several populations on the same phase-space
 or on different parts of the same phase-space and combine them afterwards
 this improves the global outcome.

void TrimPopulation( )
 after adding another population or givingHint, the true size of the population
 may be bigger than the size which was given at createPopulation
 trimPopulation should be called (if necessary) after having checked the
 individuals fitness with calculateFitness

void MakeChildren()
 does what the name says,... it creates children out of members of the
 current generation
 children have a combination of the coefficients of their parents

TMVA::GeneticGenes MakeSex( TMVA::GeneticGenes male, TMVA::GeneticGenes female )
 this function takes two individuals and produces offspring by mixing (recombining) their
 coefficients

void MakeMutants( Double_t probability, Bool_t near, Double_t spread, Bool_t mirror )
 produces offspring which is are mutated versions of their parents
 Parameters:
         double probability : gives the probability (in percent) of a mutation of a coefficient
         bool near : if true, the mutation will produce a new coefficient which is "near" the old one
                     (gaussian around the current value)
         double spread : if near==true, spread gives the sigma of the gaussian
         bool mirror : if the new value obtained would be outside of the given constraints
                    the value is mapped between the constraints again. This can be done either
                    by a kind of periodic boundary conditions or mirrored at the boundary.
                    (mirror = true seems more "natural")

void Mutate( Double_t probability , Int_t startIndex, Bool_t near, Double_t spread, Bool_t mirror )
 mutates the individuals in the genePool
 Parameters:
         double probability : gives the probability (in percent) of a mutation of a coefficient
         int startIndex : leaves unchanged (without mutation) the individuals which are better ranked
                     than indicated by "startIndex". This means: if "startIndex==3", the first (and best)
                     three individuals are not mutaded. This allows to preserve the best result of the
                     current Generation for the next generation.
         bool near : if true, the mutation will produce a new coefficient which is "near" the old one
                     (gaussian around the current value)
         double spread : if near==true, spread gives the sigma of the gaussian
         bool mirror : if the new value obtained would be outside of the given constraints
                    the value is mapped between the constraints again. This can be done either
                    by a kind of periodic boundary conditions or mirrored at the boundary.
                    (mirror = true seems more "natural")

void AddFactor( Double_t from, Double_t to )
 adds a new coefficient to the individuals.
 Parameters:
          double from : minimum value of the coefficient
          double to : maximum value of the coefficient

TMVA::GeneticGenes* GetGenes( Int_t index )
 gives back the "Genes" of the population with the given index.

Double_t GetFitness( Int_t index )
 gives back the calculated fitness of the individual with the given index
 (after the evaluation of the fitness ["calculateFitness"] index==0
 is the best individual.

void ClearResults()
 delete the results of the last calculation of the fitnesses of the
 population.
 (to prepare a new Generation)

TMVA::GeneticGenes* GetGenes()
 get the Genes of where an internal pointer is pointing to in the population

Double_t GetFitness()
 gives back the currently calculated fitness

void Reset()
 prepare for a new generation

Bool_t SetFitness( TMVA::GeneticGenes *g, Double_t fitness, Bool_t add )
 set the fitness of "g" to the value "fitness".
 add==true indicates, that this individual is created newly in this generation
 if add==false, this is a reavaluation of the fitness of the individual.

void GiveHint( vector< Double_t > hint, Double_t fitness )
 if there is some good configuration of coefficients one might give this Hint to
 the genetic algorithm.
 Parameters:
       vector< double > hint : is the collection of coefficients
       double fitness : is the fitness this collection has got

void Print( Int_t untilIndex)
 make a little printout of the individuals up to index "untilIndex"
 this means, .. write out the best "untilIndex" individuals.

void Print( ostream & out, Int_t untilIndex )
 make a little printout to the stream "out" of the individuals up to index "untilIndex"
 this means, .. write out the best "untilIndex" individuals.

TH1F* VariableDistribution( Int_t varNumber, Int_t bins, Int_t min, Int_t max )
 give back a histogram with the distribution of the coefficients
 parameters:
          int bins : number of bins of the histogram
          int min : histogram minimum
          int max : maximum value of the histogram

vector<Double_t> VariableDistribution( Int_t varNumber )
 gives back all the values of coefficient "varNumber" of the current generation

~GeneticPopulation()
 destructor
GeneticPopulation()
Double_t GetCounterFitness()
Int_t GetPopulationSize()
std::multimap<Double_t, GeneticGenes >* GetGenePool()
std::multimap<Double_t, GeneticGenes >* GetNewGenePool()
std::vector< TMVA::GeneticRange* >& GetRanges()

Author: Peter Speckmayer
Last update: root/tmva $Id: GeneticPopulation.cxx,v 1.4 2006/05/31 14:01:33 rdm Exp $
Copyright (c) 2005: *


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.