#include <iostream>
#include <iomanip>
#include "TMVA/GeneticPopulation.h"
#include "TMVA/GeneticGenes.h"
#include "Riostream.h"
#include "Rstrstream.h"
#include "TSystem.h"
#include "TRandom3.h"
#include "TH1.h"
ClassImp(TMVA::GeneticPopulation)
TMVA::GeneticPopulation::GeneticPopulation()
: fLogger( "GeneticPopulation" )
{
fRandomGenerator = new TRandom3( 100 );
fRandomGenerator->Uniform(0.,1.);
fGenePool = new multimap<Double_t, TMVA::GeneticGenes>();
fNewGenePool = new multimap<Double_t, TMVA::GeneticGenes>();
fCounterFitness = 0;
}
void TMVA::GeneticPopulation::SetRandomSeed( UInt_t seed )
{
fRandomGenerator->SetSeed( seed );
}
void TMVA::GeneticPopulation::CreatePopulation( Int_t size )
{
fPopulationSize = size;
fGenePool->clear();
fNewGenePool->clear();
vector< TMVA::GeneticRange* >::iterator rIt;
vector< Double_t > newEntry;
for (Int_t i=0; i<fPopulationSize; i++) {
newEntry.clear();
for (rIt = fRanges.begin(); rIt<fRanges.end(); rIt++) {
newEntry.push_back( (*rIt)->Random() );
}
entry e(0, TMVA::GeneticGenes( newEntry) );
fGenePool->insert( e );
}
fCounter = fGenePool->begin();
}
void TMVA::GeneticPopulation::AddPopulation( TMVA::GeneticPopulation *strangers )
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
for (it = strangers->fGenePool->begin(); it != strangers->fGenePool->end(); it++ ) {
GiveHint( it->second.GetFactors(), it->first );
}
}
void TMVA::GeneticPopulation::AddPopulation( TMVA::GeneticPopulation& strangers )
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
for (it = strangers.fGenePool->begin(); it != strangers.fGenePool->end(); it++ ) {
GiveHint( it->second.GetFactors(), it->first );
}
}
void TMVA::GeneticPopulation::TrimPopulation( )
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it = fGenePool->begin() ;
for (Int_t i=0; i<fPopulationSize; i++ ) it++;
fGenePool->erase( it, fGenePool->end()-- );
}
void TMVA::GeneticPopulation::NextGeneration()
{
fGenePool->swap( (*fNewGenePool) );
fNewGenePool->clear();
}
void TMVA::GeneticPopulation::MakeChildren()
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
multimap<Double_t, TMVA::GeneticGenes >::iterator it2;
Int_t pos = 0;
Int_t n = 0;
for (it = fGenePool->begin(); it != fGenePool->end(); it++) {
if (n< (Int_t)(fGenePool->size()/2)) {
fNewGenePool->insert( entry(0, it->second) );
pos = (Int_t)fRandomGenerator->Integer( fGenePool->size()/2 );
it2 = fGenePool->begin();
for (Int_t i=0; i<pos; i++) it2++;
fNewGenePool->insert( entry( 0, MakeSex( it->second, it2->second ) ) );
}
else continue;
n++;
}
}
TMVA::GeneticGenes TMVA::GeneticPopulation::MakeSex( TMVA::GeneticGenes male,
TMVA::GeneticGenes female )
{
vector< Double_t > child;
vector< Double_t >::iterator itM;
vector< Double_t >::iterator itF = female.GetFactors().begin();
for (itM = male.GetFactors().begin(); itM < male.GetFactors().end(); itM++) {
if (fRandomGenerator->Integer( 2 ) == 0) {
child.push_back( (*itM) );
}else{
child.push_back( (*itF) );
}
itF++;
}
return TMVA::GeneticGenes( child );
}
void TMVA::GeneticPopulation::MakeMutants( Double_t probability, Bool_t near,
Double_t spread, Bool_t mirror )
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
int n=0;
for (it = fGenePool->begin(); it != fGenePool->end(); it++) {
if (n< (fPopulationSize/2)) {
fNewGenePool->insert( entry(0, GeneticPopulation::Mutate( it->second, probability, near, spread, mirror ) ) );
fNewGenePool->insert( entry(1, GeneticPopulation::Mutate( it->second, probability, near, spread, mirror ) ) );
}
else continue;
n++;
}
}
void TMVA::GeneticPopulation::MakeCopies( int number )
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
Int_t n = 0;
for (it = fGenePool->begin(); it != fGenePool->end(); it++) {
if (n< number) {
fNewGenePool->insert( entry(0, it->second) );
}
else continue;
n++;
}
}
TMVA::GeneticGenes TMVA::GeneticPopulation::Mutate( TMVA::GeneticGenes individual,
Double_t probability,
Bool_t near,
Double_t spread,
Bool_t mirror ){
vector< Double_t > child;
vector< Double_t >::iterator vec;
vector< TMVA::GeneticRange* >::iterator vecRange;
double val;
vecRange = fRanges.begin();
for (vec = individual.GetFactors().begin(); vec < (individual.GetFactors()).end(); vec++) {
if (fRandomGenerator->Uniform( 100 ) <= probability) {
val = (*vecRange)->Random( near, (*vec), spread, mirror );
child.push_back( val );
}
vecRange++;
}
return TMVA::GeneticGenes( child );
}
void TMVA::GeneticPopulation::Mutate( Double_t probability , Int_t startIndex,
Bool_t near, Double_t spread, Bool_t mirror )
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
Int_t index = 0;
vector< Double_t >::iterator vec;
vector< TMVA::GeneticRange* >::iterator vecRange;
for (it = fGenePool->begin(); it != fGenePool->end(); it++) {
if (index >= startIndex) {
vecRange = fRanges.begin();
for (vec = (it->second.GetFactors()).begin(); vec < (it->second.GetFactors()).end(); vec++) {
if (fRandomGenerator->Uniform( 100 ) <= probability) {
(*vec) = (*vecRange)->Random( near, (*vec), spread, mirror );
}
vecRange++;
}
}
index++;
}
}
void TMVA::GeneticPopulation::AddFactor( TMVA::Interval *interval )
{
fRanges.push_back( new TMVA::GeneticRange( fRandomGenerator, interval ) );
}
TMVA::GeneticGenes* TMVA::GeneticPopulation::GetGenes( Int_t index )
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it = fGenePool->begin();
for (Int_t i=0; i<index; i++) it++;
return &(it->second);
}
Double_t TMVA::GeneticPopulation::GetFitness( Int_t index )
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it = fGenePool->begin();
for (Int_t i=0; i<index; i++) it++;
return it->first;
}
void TMVA::GeneticPopulation::ClearResults()
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
for (it = fGenePool->begin(); it!=fGenePool->end(); it++) {
it->second.ClearResults();
}
}
TMVA::GeneticGenes* TMVA::GeneticPopulation::GetGenes()
{
TMVA::GeneticGenes *g;
if (fCounter == fGenePool->end() ) {
g = new TMVA::GeneticGenes();
return g;
}
g = &(fCounter->second);
fCounterFitness = fCounter->first;
return g;
}
Double_t TMVA::GeneticPopulation::GetFitness()
{
if (fCounter == fGenePool->end() ) {
Reset();
return -1.;
}
return fCounter->first;
}
void TMVA::GeneticPopulation::Reset()
{
fCounter = fGenePool->begin();
fNewGenePool->clear();
}
Bool_t TMVA::GeneticPopulation::SetFitness( TMVA::GeneticGenes *g, Double_t fitness, Bool_t add )
{
if (add) g->GetResults().push_back( fitness );
fNewGenePool->insert( entry( fitness, *g) );
fCounter++;
if (fCounter == fGenePool->end()) {
fGenePool->swap( (*fNewGenePool) );
fCounter = fGenePool->begin();
Reset();
return kFALSE;
}
return kTRUE;
}
void TMVA::GeneticPopulation::GiveHint( vector< Double_t >& hint, Double_t fitness )
{
TMVA::GeneticGenes g;
g.GetFactors().assign( hint.begin(), hint.end() );
fGenePool->insert( entry( fitness, g ) );
}
void TMVA::GeneticPopulation::Print( Int_t untilIndex )
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
Int_t n;
for (it = fGenePool->begin(); it != fGenePool->end(); it++) {
if (untilIndex >= -1 ) {
if (untilIndex == -1 ) return;
untilIndex--;
}
n = 0;
fLogger << "fitness: " << it->first << " ";
for (vector< Double_t >::iterator vec = it->second.GetFactors().begin();
vec < it->second.GetFactors().end(); vec++ ) {
fLogger << "f_" << n++ << ": " << (*vec) << " ";
}
fLogger << Endl;
}
}
void TMVA::GeneticPopulation::Print( ostream & out, Int_t untilIndex )
{
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
Int_t n;
for (it = fGenePool->begin(); it != fGenePool->end(); it++) {
if (untilIndex > -1 ) {
untilIndex--;
if (untilIndex == -1 ) return;
}
n = 0;
out << "fitness: " << it->first << " ";
for (vector< Double_t >::iterator vec = it->second.GetFactors().begin();
vec < it->second.GetFactors().end(); vec++) {
out << "f_" << n++ << ": " << (*vec) << " ";
}
out << endl;
}
}
TH1F* TMVA::GeneticPopulation::VariableDistribution( Int_t varNumber, Int_t bins,
Int_t min, Int_t max )
{
std::stringstream histName;
histName.clear();
histName.str("v");
histName << varNumber;
TH1F *hist = new TH1F( histName.str().c_str(),histName.str().c_str(), bins,min,max );
hist->SetBit(TH1::kCanRebin);
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
for (it = fGenePool->begin(); it != fGenePool->end(); it++) {
hist->Fill( it->second.GetFactors().at(varNumber));
}
return hist;
}
vector<Double_t> TMVA::GeneticPopulation::VariableDistribution( Int_t varNumber )
{
vector< Double_t > varDist;
multimap<Double_t, TMVA::GeneticGenes >::iterator it;
for (it = fGenePool->begin(); it != fGenePool->end(); it++) {
varDist.push_back( it->second.GetFactors().at( varNumber ) );
}
return varDist;
}
TMVA::GeneticPopulation::~GeneticPopulation()
{
if (fRandomGenerator != NULL) delete fRandomGenerator;
if (fGenePool != NULL) delete fGenePool;
if (fNewGenePool != NULL) delete fNewGenePool;
std::vector<GeneticRange*>::iterator it = fRanges.begin();
for (;it!=fRanges.end(); it++) delete *it;
}
Last change: Wed Jun 25 08:48:12 2008
Last generated: 2008-06-25 08:48
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.