00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066 #include <string>
00067 #include <cstdlib>
00068 #include <iostream>
00069
00070 #include "TMatrix.h"
00071 #include "TObjString.h"
00072 #include "Riostream.h"
00073 #include "TMath.h"
00074
00075 #include "TMVA/ClassifierFactory.h"
00076 #include "TMVA/MethodCFMlpANN.h"
00077 #include "TMVA/MethodCFMlpANN_def.h"
00078 #include "TMVA/Tools.h"
00079
00080 REGISTER_METHOD(CFMlpANN)
00081
00082 ClassImp(TMVA::MethodCFMlpANN)
00083
00084
00085 namespace TMVA {
00086 Int_t MethodCFMlpANN_nsel = 0;
00087 }
00088
00089 TMVA::MethodCFMlpANN* TMVA::MethodCFMlpANN::fgThis = 0;
00090
00091
00092 TMVA::MethodCFMlpANN::MethodCFMlpANN( const TString& jobName,
00093 const TString& methodTitle,
00094 DataSetInfo& theData,
00095 const TString& theOption,
00096 TDirectory* theTargetDir ) :
00097 TMVA::MethodBase( jobName, Types::kCFMlpANN, methodTitle, theData, theOption, theTargetDir ),
00098 fData(0),
00099 fClass(0),
00100 fNlayers(0),
00101 fNcycles(0),
00102 fNodes(0),
00103 fYNN(0)
00104 {
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133 MethodCFMlpANN_Utils::SetLogger(&Log());
00134
00135 }
00136
00137
00138 TMVA::MethodCFMlpANN::MethodCFMlpANN( DataSetInfo& theData,
00139 const TString& theWeightFile,
00140 TDirectory* theTargetDir ):
00141 TMVA::MethodBase( Types::kCFMlpANN, theData, theWeightFile, theTargetDir ),
00142 fData(0),
00143 fClass(0),
00144 fNlayers(0),
00145 fNcycles(0),
00146 fNodes(0),
00147 fYNN(0)
00148 {
00149
00150 }
00151
00152
00153 Bool_t TMVA::MethodCFMlpANN::HasAnalysisType( Types::EAnalysisType type, UInt_t numberClasses, UInt_t )
00154 {
00155
00156 if (type == Types::kClassification && numberClasses == 2) return kTRUE;
00157 return kFALSE;
00158 }
00159
00160
00161 void TMVA::MethodCFMlpANN::DeclareOptions()
00162 {
00163
00164
00165
00166
00167 DeclareOptionRef( fNcycles =3000, "NCycles", "Number of training cycles" );
00168 DeclareOptionRef( fLayerSpec="N,N-1", "HiddenLayers", "Specification of hidden layer architecture" );
00169 }
00170
00171
00172 void TMVA::MethodCFMlpANN::ProcessOptions()
00173 {
00174
00175 fNodes = new Int_t[20];
00176 fNlayers = 2;
00177 Int_t currentHiddenLayer = 1;
00178 TString layerSpec(fLayerSpec);
00179 while(layerSpec.Length()>0) {
00180 TString sToAdd = "";
00181 if (layerSpec.First(',')<0) {
00182 sToAdd = layerSpec;
00183 layerSpec = "";
00184 }
00185 else {
00186 sToAdd = layerSpec(0,layerSpec.First(','));
00187 layerSpec = layerSpec(layerSpec.First(',')+1,layerSpec.Length());
00188 }
00189 Int_t nNodes = 0;
00190 if (sToAdd.BeginsWith("N") || sToAdd.BeginsWith("n")) { sToAdd.Remove(0,1); nNodes = GetNvar(); }
00191 nNodes += atoi(sToAdd);
00192 fNodes[currentHiddenLayer++] = nNodes;
00193 fNlayers++;
00194 }
00195 fNodes[0] = GetNvar();
00196 fNodes[fNlayers-1] = 2;
00197
00198 if (IgnoreEventsWithNegWeightsInTraining()) {
00199 Log() << kFATAL << "Mechanism to ignore events with negative weights in training not yet available for method: "
00200 << GetMethodTypeName()
00201 << " --> please remove \"IgnoreNegWeightsInTraining\" option from booking string."
00202 << Endl;
00203 }
00204
00205 Log() << kINFO << "Use configuration (nodes per layer): in=";
00206 for (Int_t i=0; i<fNlayers-1; i++) Log() << kINFO << fNodes[i] << ":";
00207 Log() << kINFO << fNodes[fNlayers-1] << "=out" << Endl;
00208
00209
00210 Log() << "Use " << fNcycles << " training cycles" << Endl;
00211
00212 Int_t nEvtTrain = Data()->GetNTrainingEvents();
00213
00214
00215 if (nEvtTrain>0) {
00216
00217
00218 fData = new TMatrix( nEvtTrain, GetNvar() );
00219 fClass = new vector<Int_t>( nEvtTrain );
00220
00221
00222
00223 UInt_t ivar;
00224 for (Int_t ievt=0; ievt<nEvtTrain; ievt++) {
00225 const Event * ev = GetEvent(ievt);
00226
00227
00228 (*fClass)[ievt] = DataInfo().IsSignal(ev) ? 1 : 2;
00229
00230
00231 for (ivar=0; ivar<GetNvar(); ivar++) {
00232 (*fData)( ievt, ivar ) = ev->GetValue(ivar);
00233 }
00234 }
00235
00236
00237
00238 }
00239
00240 }
00241
00242
00243 void TMVA::MethodCFMlpANN::Init( void )
00244 {
00245
00246
00247
00248 SetNormalised( kTRUE );
00249
00250
00251 fgThis = this;
00252
00253
00254 TMVA::MethodCFMlpANN_nsel = 0;
00255 }
00256
00257
00258 TMVA::MethodCFMlpANN::~MethodCFMlpANN( void )
00259 {
00260
00261 delete fData;
00262 delete fClass;
00263 delete[] fNodes;
00264
00265 if (fYNN!=0) {
00266 for (Int_t i=0; i<fNlayers; i++) delete[] fYNN[i];
00267 delete[] fYNN;
00268 fYNN=0;
00269 }
00270 }
00271
00272
00273 void TMVA::MethodCFMlpANN::Train( void )
00274 {
00275
00276
00277 Double_t dumDat(0);
00278 Int_t ntrain(Data()->GetNTrainingEvents());
00279 Int_t ntest(0);
00280 Int_t nvar(GetNvar());
00281 Int_t nlayers(fNlayers);
00282 Int_t *nodes = new Int_t[nlayers];
00283 Int_t ncycles(fNcycles);
00284
00285 for (Int_t i=0; i<nlayers; i++) nodes[i] = fNodes[i];
00286
00287 if (fYNN != 0) {
00288 for (Int_t i=0; i<fNlayers; i++) delete[] fYNN[i];
00289 delete[] fYNN;
00290 fYNN = 0;
00291 }
00292 fYNN = new Double_t*[nlayers];
00293 for (Int_t layer=0; layer<nlayers; layer++)
00294 fYNN[layer] = new Double_t[fNodes[layer]];
00295
00296
00297 #ifndef R__WIN32
00298 Train_nn( &dumDat, &dumDat, &ntrain, &ntest, &nvar, &nlayers, nodes, &ncycles );
00299 #else
00300 Log() << kWARNING << "<Train> sorry CFMlpANN does not run on Windows" << Endl;
00301 #endif
00302
00303 delete [] nodes;
00304 }
00305
00306
00307 Double_t TMVA::MethodCFMlpANN::GetMvaValue( Double_t* err, Double_t* errUpper )
00308 {
00309
00310 Bool_t isOK = kTRUE;
00311
00312 const Event* ev = GetEvent();
00313
00314
00315 vector<Double_t> inputVec( GetNvar() );
00316 for (UInt_t ivar=0; ivar<GetNvar(); ivar++) inputVec[ivar] = ev->GetValue(ivar);
00317
00318 Double_t myMVA = EvalANN( inputVec, isOK );
00319 if (!isOK) Log() << kFATAL << "EvalANN returns (!isOK) for event " << Endl;
00320
00321
00322 NoErrorCalc(err, errUpper);
00323
00324 return myMVA;
00325 }
00326
00327
00328 Double_t TMVA::MethodCFMlpANN::EvalANN( vector<Double_t>& inVar, Bool_t& isOK )
00329 {
00330
00331
00332
00333 Double_t* xeev = new Double_t[GetNvar()];
00334 for (UInt_t ivar=0; ivar<GetNvar(); ivar++) xeev[ivar] = inVar[ivar];
00335
00336
00337 isOK = kTRUE;
00338 for (UInt_t jvar=0; jvar<GetNvar(); jvar++) {
00339
00340 if (fVarn_1.xmax[jvar] < xeev[jvar]) xeev[jvar] = fVarn_1.xmax[jvar];
00341 if (fVarn_1.xmin[jvar] > xeev[jvar]) xeev[jvar] = fVarn_1.xmin[jvar];
00342 if (fVarn_1.xmax[jvar] == fVarn_1.xmin[jvar]) {
00343 isOK = kFALSE;
00344 xeev[jvar] = 0;
00345 }
00346 else {
00347 xeev[jvar] = xeev[jvar] - ((fVarn_1.xmax[jvar] + fVarn_1.xmin[jvar])/2);
00348 xeev[jvar] = xeev[jvar] / ((fVarn_1.xmax[jvar] - fVarn_1.xmin[jvar])/2);
00349 }
00350 }
00351
00352 NN_ava( xeev );
00353
00354 Double_t retval = 0.5*(1.0 + fYNN[fParam_1.layerm-1][0]);
00355
00356 delete [] xeev;
00357
00358 return retval;
00359 }
00360
00361
00362 void TMVA::MethodCFMlpANN::NN_ava( Double_t* xeev )
00363 {
00364
00365 for (Int_t ivar=0; ivar<fNeur_1.neuron[0]; ivar++) fYNN[0][ivar] = xeev[ivar];
00366
00367 for (Int_t layer=1; layer<fParam_1.layerm; layer++) {
00368 for (Int_t j=1; j<=fNeur_1.neuron[layer]; j++) {
00369
00370 Double_t x = Ww_ref(fNeur_1.ww, layer+1,j);
00371
00372 for (Int_t k=1; k<=fNeur_1.neuron[layer-1]; k++) {
00373 x += fYNN[layer-1][k-1]*W_ref(fNeur_1.w, layer+1, j, k);
00374 }
00375 fYNN[layer][j-1] = NN_fonc( layer, x );
00376 }
00377 }
00378 }
00379
00380
00381 Double_t TMVA::MethodCFMlpANN::NN_fonc( Int_t i, Double_t u ) const
00382 {
00383
00384 Double_t f(0);
00385
00386 if (u/fDel_1.temp[i] > 170) f = +1;
00387 else if (u/fDel_1.temp[i] < -170) f = -1;
00388 else {
00389 Double_t yy = TMath::Exp(-u/fDel_1.temp[i]);
00390 f = (1 - yy)/(1 + yy);
00391 }
00392
00393 return f;
00394 }
00395
00396
00397 void TMVA::MethodCFMlpANN::ReadWeightsFromStream( istream & istr )
00398 {
00399
00400 TString var;
00401
00402
00403 UInt_t nva(0), lclass(0);
00404 istr >> nva >> lclass;
00405
00406 if (GetNvar() != nva)
00407 Log() << kFATAL << "<ReadWeightsFromFile> mismatch in number of variables" << Endl;
00408
00409
00410 if (lclass != 2)
00411 Log() << kFATAL << "<ReadWeightsFromFile> mismatch in number of classes" << Endl;
00412
00413
00414 if (istr.eof( ))
00415 Log() << kFATAL << "<ReadWeightsFromStream> reached EOF prematurely " << Endl;
00416
00417
00418 for (UInt_t ivar=0; ivar<GetNvar(); ivar++)
00419 istr >> fVarn_1.xmax[ivar] >> fVarn_1.xmin[ivar];
00420
00421
00422 istr >> fParam_1.layerm;
00423
00424 if (fYNN != 0) {
00425 for (Int_t i=0; i<fNlayers; i++) delete[] fYNN[i];
00426 delete[] fYNN;
00427 fYNN = 0;
00428 }
00429 fYNN = new Double_t*[fParam_1.layerm];
00430 for (Int_t layer=0; layer<fParam_1.layerm; layer++) {
00431
00432
00433 istr >> fNeur_1.neuron[layer];
00434 fYNN[layer] = new Double_t[fNeur_1.neuron[layer]];
00435 }
00436
00437
00438 const Int_t nchar( 100 );
00439 char* dumchar = new char[nchar];
00440
00441
00442 for (Int_t layer=1; layer<=fParam_1.layerm-1; layer++) {
00443
00444 Int_t nq = fNeur_1.neuron[layer]/10;
00445 Int_t nr = fNeur_1.neuron[layer] - nq*10;
00446
00447 Int_t kk(0);
00448 if (nr==0) kk = nq;
00449 else kk = nq+1;
00450
00451 for (Int_t k=1; k<=kk; k++) {
00452 Int_t jmin = 10*k - 9;
00453 Int_t jmax = 10*k;
00454 if (fNeur_1.neuron[layer]<jmax) jmax = fNeur_1.neuron[layer];
00455 for (Int_t j=jmin; j<=jmax; j++) {
00456 istr >> Ww_ref(fNeur_1.ww, layer+1, j);
00457 }
00458 for (Int_t i=1; i<=fNeur_1.neuron[layer-1]; i++) {
00459 for (Int_t j=jmin; j<=jmax; j++) {
00460 istr >> W_ref(fNeur_1.w, layer+1, j, i);
00461 }
00462 }
00463
00464 istr.getline( dumchar, nchar );
00465 }
00466 }
00467
00468 for (Int_t layer=0; layer<fParam_1.layerm; layer++) {
00469
00470
00471 istr.getline( dumchar, nchar );
00472 istr.getline( dumchar, nchar );
00473
00474 istr >> fDel_1.temp[layer];
00475 }
00476
00477
00478 if ((Int_t)GetNvar() != fNeur_1.neuron[0]) {
00479 Log() << kFATAL << "<ReadWeightsFromFile> mismatch in zeroth layer:"
00480 << GetNvar() << " " << fNeur_1.neuron[0] << Endl;
00481 }
00482
00483 fNlayers = fParam_1.layerm;
00484 delete[] dumchar;
00485 }
00486
00487
00488 Int_t TMVA::MethodCFMlpANN::DataInterface( Double_t* , Double_t* ,
00489 Int_t* , Int_t* ,
00490 Int_t* , Int_t* nvar,
00491 Double_t* xpg, Int_t* iclass, Int_t* ikend )
00492 {
00493
00494
00495
00496 *ikend = 0;
00497
00498
00499 TMVA::MethodCFMlpANN* opt = TMVA::MethodCFMlpANN::This();
00500
00501
00502 if (0 == xpg) {
00503 Log() << kFATAL << "ERROR in MethodCFMlpANN_DataInterface zero pointer xpg" << Endl;
00504 }
00505 if (*nvar != (Int_t)opt->GetNvar()) {
00506 Log() << kFATAL << "ERROR in MethodCFMlpANN_DataInterface mismatch in num of variables: "
00507 << *nvar << " " << opt->GetNvar() << Endl;
00508 }
00509
00510
00511 *iclass = (int)opt->GetClass( TMVA::MethodCFMlpANN_nsel );
00512 for (UInt_t ivar=0; ivar<opt->GetNvar(); ivar++)
00513 xpg[ivar] = (double)opt->GetData( TMVA::MethodCFMlpANN_nsel, ivar );
00514
00515 ++TMVA::MethodCFMlpANN_nsel;
00516
00517 return 0;
00518 }
00519
00520
00521 void TMVA::MethodCFMlpANN::AddWeightsXMLTo( void* parent ) const
00522 {
00523
00524
00525 void *wght = gTools().AddChild(parent, "Weights");
00526 gTools().AddAttr(wght,"NVars",fParam_1.nvar);
00527 gTools().AddAttr(wght,"NClasses",fParam_1.lclass);
00528 gTools().AddAttr(wght,"NLayers",fParam_1.layerm);
00529 void* minmaxnode = gTools().AddChild(wght, "VarMinMax");
00530 stringstream s;
00531 s.precision( 16 );
00532 for (Int_t ivar=0; ivar<fParam_1.nvar; ivar++)
00533 s << std::scientific << fVarn_1.xmin[ivar] << " " << fVarn_1.xmax[ivar] << " ";
00534 gTools().AddRawLine( minmaxnode, s.str().c_str() );
00535 void* neurons = gTools().AddChild(wght, "NNeurons");
00536 stringstream n;
00537 n.precision( 16 );
00538 for (Int_t layer=0; layer<fParam_1.layerm; layer++)
00539 n << std::scientific << fNeur_1.neuron[layer] << " ";
00540 gTools().AddRawLine( neurons, n.str().c_str() );
00541 for (Int_t layer=1; layer<fParam_1.layerm; layer++) {
00542 void* layernode = gTools().AddChild(wght, "Layer"+gTools().StringFromInt(layer));
00543 gTools().AddAttr(layernode,"NNeurons",fNeur_1.neuron[layer]);
00544 void* neuronnode=NULL;
00545 for (Int_t neuron=0; neuron<fNeur_1.neuron[layer]; neuron++) {
00546 neuronnode = gTools().AddChild(layernode,"Neuron"+gTools().StringFromInt(neuron));
00547 stringstream weights;
00548 weights.precision( 16 );
00549 weights << std::scientific << Ww_ref(fNeur_1.ww, layer+1, neuron+1);
00550 for (Int_t i=0; i<fNeur_1.neuron[layer-1]; i++) {
00551 weights << " " << std::scientific << W_ref(fNeur_1.w, layer+1, neuron+1, i+1);
00552 }
00553 gTools().AddRawLine( neuronnode, weights.str().c_str() );
00554 }
00555 }
00556 void* tempnode = gTools().AddChild(wght, "LayerTemp");
00557 stringstream temp;
00558 temp.precision( 16 );
00559 for (Int_t layer=0; layer<fParam_1.layerm; layer++) {
00560 temp << std::scientific << fDel_1.temp[layer] << " ";
00561 }
00562 gTools().AddRawLine(tempnode, temp.str().c_str() );
00563 }
00564
00565 void TMVA::MethodCFMlpANN::ReadWeightsFromXML( void* wghtnode )
00566 {
00567
00568 gTools().ReadAttr( wghtnode, "NLayers",fParam_1.layerm );
00569 void* minmaxnode = gTools().GetChild(wghtnode);
00570 const char* minmaxcontent = gTools().GetContent(minmaxnode);
00571 std::stringstream content(minmaxcontent);
00572 for (UInt_t ivar=0; ivar<GetNvar(); ivar++)
00573 content >> fVarn_1.xmin[ivar] >> fVarn_1.xmax[ivar];
00574 if (fYNN != 0) {
00575 for (Int_t i=0; i<fNlayers; i++) delete[] fYNN[i];
00576 delete[] fYNN;
00577 fYNN = 0;
00578 }
00579 fYNN = new Double_t*[fParam_1.layerm];
00580 void *layernode=gTools().GetNextChild(minmaxnode);
00581 const char* neuronscontent = gTools().GetContent(layernode);
00582 stringstream ncontent(neuronscontent);
00583 for (Int_t layer=0; layer<fParam_1.layerm; layer++) {
00584
00585
00586 ncontent >> fNeur_1.neuron[layer];
00587 fYNN[layer] = new Double_t[fNeur_1.neuron[layer]];
00588 }
00589 for (Int_t layer=1; layer<fParam_1.layerm; layer++) {
00590 layernode=gTools().GetNextChild(layernode);
00591 void* neuronnode=NULL;
00592 neuronnode = gTools().GetChild(layernode);
00593 for (Int_t neuron=0; neuron<fNeur_1.neuron[layer]; neuron++) {
00594 const char* neuronweights = gTools().GetContent(neuronnode);
00595 stringstream weights(neuronweights);
00596 weights >> Ww_ref(fNeur_1.ww, layer+1, neuron+1);
00597 for (Int_t i=0; i<fNeur_1.neuron[layer-1]; i++) {
00598 weights >> W_ref(fNeur_1.w, layer+1, neuron+1, i+1);
00599 }
00600 neuronnode=gTools().GetNextChild(neuronnode);
00601 }
00602 }
00603 void* tempnode=gTools().GetNextChild(layernode);
00604 const char* temp = gTools().GetContent(tempnode);
00605 stringstream t(temp);
00606 for (Int_t layer=0; layer<fParam_1.layerm; layer++) {
00607 t >> fDel_1.temp[layer];
00608 }
00609 fNlayers = fParam_1.layerm;
00610 }
00611
00612
00613 void TMVA::MethodCFMlpANN::PrintWeights( std::ostream & o ) const
00614 {
00615
00616
00617
00618 o << "Number of vars " << fParam_1.nvar << endl;
00619 o << "Output nodes " << fParam_1.lclass << endl;
00620
00621
00622 for (Int_t ivar=0; ivar<fParam_1.nvar; ivar++)
00623 o << "Var " << ivar << " [" << fVarn_1.xmin[ivar] << " - " << fVarn_1.xmax[ivar] << "]" << endl;
00624
00625
00626 o << "Number of layers " << fParam_1.layerm << endl;
00627
00628 o << "Nodes per layer ";
00629 for (Int_t layer=0; layer<fParam_1.layerm; layer++)
00630
00631 o << fNeur_1.neuron[layer] << " ";
00632 o << endl;
00633
00634
00635 for (Int_t layer=1; layer<=fParam_1.layerm-1; layer++) {
00636
00637 Int_t nq = fNeur_1.neuron[layer]/10;
00638 Int_t nr = fNeur_1.neuron[layer] - nq*10;
00639
00640 Int_t kk(0);
00641 if (nr==0) kk = nq;
00642 else kk = nq+1;
00643
00644 for (Int_t k=1; k<=kk; k++) {
00645 Int_t jmin = 10*k - 9;
00646 Int_t jmax = 10*k;
00647 Int_t i, j;
00648 if (fNeur_1.neuron[layer]<jmax) jmax = fNeur_1.neuron[layer];
00649 for (j=jmin; j<=jmax; j++) {
00650
00651
00652 o << Ww_ref(fNeur_1.ww, layer+1, j) << " ";
00653
00654 }
00655 o << endl;
00656
00657 for (i=1; i<=fNeur_1.neuron[layer-1]; i++) {
00658 for (j=jmin; j<=jmax; j++) {
00659
00660 o << W_ref(fNeur_1.w, layer+1, j, i) << " ";
00661 }
00662 o << endl;
00663 }
00664
00665
00666 o << endl;
00667 }
00668 }
00669 for (Int_t layer=0; layer<fParam_1.layerm; layer++) {
00670 o << "Del.temp in layer " << layer << " : " << fDel_1.temp[layer] << endl;
00671 }
00672 }
00673
00674 TMVA::MethodCFMlpANN* TMVA::MethodCFMlpANN::This( void )
00675 {
00676
00677 return fgThis;
00678 }
00679 void TMVA::MethodCFMlpANN::MakeClassSpecific( std::ostream& fout, const TString& className ) const
00680 {
00681
00682 fout << " // not implemented for class: \"" << className << "\"" << endl;
00683 fout << "};" << endl;
00684 }
00685
00686
00687 void TMVA::MethodCFMlpANN::MakeClassSpecificHeader( std::ostream& , const TString& ) const
00688 {
00689
00690 }
00691
00692
00693 void TMVA::MethodCFMlpANN::GetHelpMessage() const
00694 {
00695
00696
00697
00698
00699 Log() << Endl;
00700 Log() << gTools().Color("bold") << "--- Short description:" << gTools().Color("reset") << Endl;
00701 Log() << Endl;
00702 Log() << "<None>" << Endl;
00703 Log() << Endl;
00704 Log() << gTools().Color("bold") << "--- Performance optimisation:" << gTools().Color("reset") << Endl;
00705 Log() << Endl;
00706 Log() << "<None>" << Endl;
00707 Log() << Endl;
00708 Log() << gTools().Color("bold") << "--- Performance tuning via configuration options:" << gTools().Color("reset") << Endl;
00709 Log() << Endl;
00710 Log() << "<None>" << Endl;
00711 }