star-travex
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
MuMcPrVKFV2012.C
Go to the documentation of this file.
1 /*
2  root.exe lMuDst.C MuMcPrV.C+
3 */
4 #ifndef __CINT__
5 #include <map>
6 #include <utility>
7 #include <cstdlib>
8 #include <string>
9 #include <vector>
10 #include <boost/algorithm/string/replace.hpp>
11 
12 #include "TROOT.h"
13 #include "TSystem.h"
14 #include "TH1D.h"
15 #include "TH2D.h"
16 #include "TTree.h"
17 #include "TNtuple.h"
18 #include "TChain.h"
19 #include "TFile.h"
20 #include "TString.h"
21 #include "TObjString.h"
22 #include "TArrayF.h"
23 #include "TArrayD.h"
24 #include "TMath.h"
25 #include "TMVA/Factory.h"
26 #include "TMVA/Tools.h"
27 #include "StMuDSTMaker/COMMON/StMuDst.h"
28 #include "StMuDSTMaker/COMMON/StMuEvent.h"
29 #include "StMuDSTMaker/COMMON/StMuTrack.h"
30 #include "StMuDSTMaker/COMMON/StMuPrimaryVertex.h"
31 #include "StMuDSTMaker/COMMON/StMuMcVertex.h"
32 #include "StMuDSTMaker/COMMON/StMuMcTrack.h"
33 #include "StMuDSTMaker/COMMON/StMuPrimaryTrackCovariance.h"
34 #include "utils.h"
35 #include "StMuDSTMaker/COMMON/StMuDstMaker.h"
36 #endif
37 
38 #ifdef __TMVA__
39 #include "TMVAClassification_BDT.class.C"
40 #endif /* __TMVA__ */
41 
42 struct data_t {Float_t beam, postx, prompt, cross, tof, notof, EEMC, noEEMC, chi2;};
43 const Char_t *vnames = "beam:postx:prompt:cross:tof:notof:EEMC:noEEMC:chi2";
44 
45 
46 
47 void FillData(data_t &data, const StMuPrimaryVertex *Vtx)
48 {
49  memset(&data.beam, 0, sizeof(data_t));
50  //memset(&data.postx, 0, sizeof(data_t)); //PPV
51  data.beam = Vtx->isBeamConstrained() ? 1 : 0;
52  data.postx = Vtx->nPostXtracks();
53  data.prompt = Vtx->nPromptTracks();
54  data.cross = Vtx->nCrossCentralMembrane();
55  data.tof = (Vtx->nCTBMatch() + Vtx->nBTOFMatch());
56  data.notof = (Vtx->nCTBNotMatch() + Vtx->nBTOFNotMatch());
57  data.EEMC = Vtx->nEEMCMatch();
58  data.noEEMC = Vtx->nEEMCNotMatch();
59  data.chi2 = Vtx->chiSquared();
60 }
61 
62 
63 //________________________________________________________________________________
64 Bool_t Accept(const StMuTrack *gTrack = 0)
65 {
66  if (! gTrack) return kFALSE;
67  if (! gTrack->idTruth()) return kFALSE;
68  if (! gTrack->charge()) return kFALSE;
69  if ( gTrack->flag() < 100 || gTrack->flag() % 100 == 11) return kFALSE; // bad fit or short track pointing to EEMC
70  if ( gTrack->flag() > 1000) return kFALSE; // pile up track in TPC
71  if ( gTrack->nHitsFit() < 10) return kFALSE;
72 
73  return kTRUE;
74 }
75 
76 
77 //________________________________________________________________________________
78 Bool_t AcceptVX(const StMuPrimaryVertex *Vtx = 0)
79 {
80  if (! Vtx) return kFALSE;
81  if (! Vtx->idTruth()) return kFALSE;
82  if ( Vtx->position().perp() > 3.0) return kFALSE;
83 
84  return kTRUE;
85 }
86 
87 
88 //________________________________________________________________________________
89 void ForceAnimate(unsigned int times = 0, int msecDelay = 0)
90 {
91  unsigned int counter = times;
92 
93  while ( (!times || counter) && !gSystem->ProcessEvents()) { --counter; if (msecDelay) gSystem->Sleep(msecDelay);}
94 }
95 
96 
127 void MuMcPrVKFV2012(Long64_t nevent, const char *file, const std::string& outFile, bool fillNtuple)
128 {
129 #ifdef __TMVA__
130  boost::replace_last(outFile, ".root", "");
131  outFile += ".TMVArank.root";
132 
133  // create a set of variables and declare them to the reader
134  // - the variable names must corresponds in name and type to
135  // those given in the weight file(s) that you use
136  TString separator(":");
137  TString Vnames(vnames);
138  TObjArray *array = Vnames.Tokenize(separator);
139 
140  std::vector<std::string> inputVars;
141  TIter next(array);
142  TObjString *objs;
143 
144  while ((objs = (TObjString *) next())) {
145  std::cout << objs->GetString() << std::endl;
146  }
147 
148  inputVars.push_back("beam");
149  inputVars.push_back("postx");
150  inputVars.push_back("prompt");
151  inputVars.push_back("cross");
152  inputVars.push_back("tof");
153  inputVars.push_back("notof");
154  inputVars.push_back("EEMC");
155  inputVars.push_back("noEEMC");
156  inputVars.push_back("chi2");
157 
158  std::vector<double> *inputVec = new std::vector<double>( inputVars.size() );
159  IClassifierReader *classReader = new ReadBDT( inputVars );
160 
161 #endif /* __TMVA__ */
162 
163  TFile *fOut = TFile::Open(outFile.c_str(), "recreate");
164  data_t data;
165 
166  // Book histograms
167  const int nMcRecMult = 75;
168  TArrayD xMult(nMcRecMult + 1);
169  xMult[0] = -0.5;
170 
171  for (int i = 1; i <= nMcRecMult; i++) {
172  if (xMult[i - 1] < 50) xMult[i] = xMult[i - 1] + 1; // 1 - 50
173  else if (xMult[i - 1] < 100) xMult[i] = xMult[i - 1] + 2; // 51 - 75
174  else if (xMult[i - 1] < 200) xMult[i] = xMult[i - 1] + 10; // 76 - 85
175  else xMult[i] = xMult[i - 1] + 100; // 86 -100
176  }
177 
178  TH1D *McRecMulT = new TH1D("McRecMulT", "Reconstructable multiplicity for trigger Mc Vertex", nMcRecMult, xMult.GetArray());
179  struct Name_t {
180  const Char_t *Name;
181  const Char_t *Title;
182  };
183  const Name_t HCases[3] = {
184  {"Any", "Any vertex matched with MC == 1"},
185  {"Good", "The best rank vertex with MC == 1"},
186  {"Bad", "The best rank vertex with MC != 1"}
187  };
188  const Name_t Plots[4] = {
189  {"Mult" , "the reconstructed (uncorrected) track multiplicity versus Reconstructable multiplicity"},
190  {"MultQA" , "the reconstructed (corrected for QA) track multiplicity versus Reconstructable multiplicity"},
191  {"McRecMul", "Reconstructable multiplicity"},
192  {"YvsX" , "Bad versus Good value"}
193  };
194  /* h p */
195  TH1 *hists[3][4];
196 
197  for (int h = 0; h < 3; h++) {
198  for (int p = 0; p < 4; p++) {
199  TString Name(Plots[p].Name); Name += HCases[h].Name;
200  TString Title(Plots[p].Title); Title += " for "; Title += HCases[h].Title; Title += " vertex";
201 
202  if (p < 2) hists[h][p] = new TH2D(Name, Title, nMcRecMult, xMult.GetArray(), nMcRecMult, xMult.GetArray());
203  else if (p == 2) hists[h][p] = new TH1D(Name, Title, nMcRecMult, xMult.GetArray());
204  }
205  }
206 
207  TNtuple *VertexG = new TNtuple("VertexG", "good vertex & global params info", vnames);
208  TNtuple *VertexB = new TNtuple("VertexB", "bad vertex & global params info", vnames);
209  // ----------------------------------------------
210  StMuDstMaker *maker = new StMuDstMaker(0, 0, "", file, "st:MuDst.root", 1e9); // set up maker in read mode
211  // 0,0 this mean read mode
212  // dir read all files in this directory
213  // file bla.lis real all file in this list, if (file!="") dir is ignored
214  // filter apply filter to filenames, multiple filters are separated by ':'
215  // 10 maximum number of file to read
216  maker->SetStatus("*", 0);
217 
218  std::vector<std::string> activeBranchNames = {
219  "MuEvent",
220  "PrimaryVertices",
221  "StStMuMcVertex",
222  "StStMuMcTrack"
223  };
224 
225  // Set Active braches
226  for (const auto& branchName : activeBranchNames)
227  maker->SetStatus(branchName.c_str(), 1);
228 
229  TChain *tree = maker->chain();
230  Long64_t nentries = tree->GetEntries();
231  nevent = TMath::Min(nevent, nentries);
232  std::cout << nentries << " events in chain " << nevent << " will be read." << std::endl;
233  tree->SetCacheSize(-1); //by setting the read cache to -1 we set it to the AutoFlush value when writing
234  tree->SetCacheLearnEntries(1); //one entry is sufficient to learn
235  tree->SetCacheEntryRange(0, nevent);
236 
237  for (Long64_t ev = 0; ev < nevent; ev++) {
238  if (maker->Make()) break;
239 
240  StMuDst *muDst = maker->muDst(); // get a pointer to the StMuDst class, the class that points to all the data
241  StMuEvent *muEvent = muDst->event(); // get a pointer to the class holding event-wise information
242  int referenceMultiplicity = muEvent->refMult(); // get the reference multiplicity
243 
244  TClonesArray *PrimaryVertices = muDst->primaryVertices();
245  int nPrimaryVertices = PrimaryVertices->GetEntriesFast();
246 
247  TClonesArray *MuMcVertices = muDst->mcArray(0);
248  int nMuMcVertices = MuMcVertices->GetEntriesFast();
249 
250  TClonesArray *MuMcTracks = muDst->mcArray(1);
251  int nMuMcTracks = MuMcTracks->GetEntriesFast();
252 
253  if ( nevent >= 10 && ev % int(nevent*0.1) == 0 )
254  {
255  std::cout << "Event #" << ev << "\tRun\t" << muEvent->runId()
256  << "\tId: " << muEvent->eventId()
257  << " refMult= " << referenceMultiplicity
258  << "\tPrimaryVertices " << nPrimaryVertices
259  << "\t" << " " << nMuMcVertices
260  << "\t" << " " << nMuMcTracks
261  << std::endl;
262  }
263 
264  // const Double_t field = muEvent->magneticField()*kilogauss;
265  if (! nMuMcVertices || ! nMuMcTracks || nPrimaryVertices <= 0) {
266  std::cout << "Ev. " << ev << " has no MC information ==> skip it" << std::endl;
267  std::cout << "OR no reconstructed verticies found" << std::endl;
268  continue;
269  }
270 
271  // Count number of MC tracks at a vertex with TPC reconstructable tracks
272  std::multimap<int, int> Mc2McHitTracks;
273 
274  for (int m = 0; m < nMuMcTracks; m++) {
275  StMuMcTrack *McTrack = (StMuMcTrack *) MuMcTracks->UncheckedAt(m);
276 
277  if (McTrack->No_tpc_hit() < 15) continue;
278 
279  Mc2McHitTracks.insert(std::pair<int, int>(McTrack->IdVx(), McTrack->Id()));
280  }
281 
282  // This is the "reconstructable" track multiplicity
283  int nMcTracksWithHits = Mc2McHitTracks.count(1);
284 
285  // Let's skip events in which we do not expect to reconstruct any tracks
286  // (and thus vertex) from the primary vertex
287  if (nMcTracksWithHits <= 0) continue;
288 
289  // This is our denominator histogram for efficiencies
290  McRecMulT->Fill(nMcTracksWithHits);
291 
292  // ============= Build map between Rc and Mc vertices
293  std::map<StMuPrimaryVertex *, StMuMcVertex *> reco2McVertices;
294  TArrayF vertexRanks(nPrimaryVertices);
295  int mcMatchedVertexIndex = -1; // any vertex with MC==1 and highest reconstrated multiplicity.
296  int vertexMaxMultiplicity = -1;
297 
298  // First loop over all verticies in this event. There is at least one
299  // must be available
300  for (int recoVertexIndex = 0; recoVertexIndex < nPrimaryVertices; recoVertexIndex++)
301  {
302  vertexRanks[recoVertexIndex] = -1e10;
303 
304  StMuPrimaryVertex *recoVertex = (StMuPrimaryVertex *) PrimaryVertices->UncheckedAt(recoVertexIndex);
305 
306  if ( !AcceptVX(recoVertex) ) continue;
307 
308  // Check Mc
309  if (recoVertex->idTruth() < 0 || recoVertex->idTruth() > nMuMcVertices) {
310  std::cout << "ERROR: Illegal idTruth " << recoVertex->idTruth() << " The track is ignored" << std::endl;
311  continue;
312  }
313 
314  StMuMcVertex *mcVertex = (StMuMcVertex *) MuMcVertices->UncheckedAt(recoVertex->idTruth() - 1);
315 
316  if (mcVertex->Id() != recoVertex->idTruth()) {
317  std::cout << "ERROR: Mismatched idTruth " << recoVertex->idTruth() << " and mcVertex Id " << mcVertex->Id()
318  << " The vertex is ignored" << std::endl;
319  continue;
320  }
321 
322  reco2McVertices[recoVertex] = mcVertex;
323  vertexRanks[recoVertexIndex] = recoVertex->ranking();
324 
325  if (recoVertex->idTruth() == 1 && recoVertex->noTracks() > vertexMaxMultiplicity)
326  {
327  mcMatchedVertexIndex = recoVertexIndex;
328  vertexMaxMultiplicity = recoVertex->noTracks();
329  }
330 
331  FillData(data, recoVertex);
332 
333 #ifdef __TMVA__
334  Float_t *dataArray = &data.beam;
335 
336  for (size_t j = 0; j < inputVec->size(); j++)
337  (*inputVec)[j] = dataArray[j];
338 
339  vertexRanks[recoVertexIndex] = classReader->GetMvaValue( *inputVec );
340 #endif
341  }
342 
343  // If we reconstructed a vertex which matches the MC one we fill the
344  // numerator of the "Any" efficiency histogram
345  if (mcMatchedVertexIndex != -1) {
346 
347  StMuPrimaryVertex *recoVertexMatchedMc = (StMuPrimaryVertex*) PrimaryVertices->UncheckedAt(mcMatchedVertexIndex);
348 
349  double nTracks = recoVertexMatchedMc->noTracks();
350  double nTracksQA = nTracks * recoVertexMatchedMc->qaTruth() / 100.;
351 
352  hists[0][0]->Fill(nMcTracksWithHits, nTracks);
353  hists[0][1]->Fill(nMcTracksWithHits, nTracksQA);
354  hists[0][2]->Fill(nMcTracksWithHits);
355  }
356 
357  // Now deal with the highest rank vertex
358  int maxRankVertexIndex = TMath::LocMax(nPrimaryVertices, vertexRanks.GetArray());
359 
360  StMuPrimaryVertex *recoVertexMaxRank = (StMuPrimaryVertex*) PrimaryVertices->UncheckedAt(maxRankVertexIndex);
361  StMuMcVertex *mcVertex = reco2McVertices[recoVertexMaxRank];
362 
363  double nTracks = recoVertexMaxRank->noTracks();
364  double nTracksQA = nTracks * recoVertexMaxRank->qaTruth() / 100.;
365 
366  // Fill numerator for "good" and "bad" efficiencies
367  int h = ( mcVertex && mcVertex->Id() == 1) ? 1 : 2;
368 
369  hists[h][0]->Fill(nMcTracksWithHits, nTracks);
370  hists[h][1]->Fill(nMcTracksWithHits, nTracksQA);
371  hists[h][2]->Fill(nMcTracksWithHits);
372 
373 
374  // Proceed with filling ntuple only if requested by the user
375  if ( !fillNtuple ) continue;
376 
377 
378  // Second loop over all verticies in this event
379  for (int recoVertexIndex = 0; recoVertexIndex < nPrimaryVertices; recoVertexIndex++)
380  {
381  StMuPrimaryVertex *recoVertex = (StMuPrimaryVertex *) PrimaryVertices->UncheckedAt(recoVertexIndex);
382 
383  if ( !AcceptVX(recoVertex) ) continue;
384 
385  StMuMcVertex *mcVertex = reco2McVertices[recoVertex];
386 
387  if ( !mcVertex ) {
388  std::cout << "No Match from RC to MC" << std::endl;
389  continue;
390  }
391 
392  if (vtxeval::gDebugFlag) {
393  std::cout << Form("Vx[%3i]", recoVertexIndex) << *recoVertex << " " << *mcVertex;
394  int nMcTracksWithHitsatL = Mc2McHitTracks.count(recoVertex->idTruth());
395  std::cout << Form("Number of McTkHit %4i rank %8.3f", nMcTracksWithHitsatL, vertexRanks[recoVertexIndex]);
396  }
397 
398  int IdPar = mcVertex->IdParTrk();
399 
400  if (IdPar > 0 && IdPar <= nMuMcTracks) {
401  StMuMcTrack *mcTrack = (StMuMcTrack *) MuMcTracks->UncheckedAt(IdPar - 1);
402 
403  if (mcTrack && vtxeval::gDebugFlag) std::cout << " " << mcTrack->GeName();
404  }
405 
406  FillData(data, recoVertex);
407 
408  double nTracks = recoVertex->noTracks();
409 
410  if (mcVertex->Id() == 1 && nTracks == vertexMaxMultiplicity) {// good
411  VertexG->Fill(&data.beam);
412  }
413  else { // bad
414  VertexB->Fill(&data.beam);
415  }
416  }
417 
418  if ( !gROOT->IsBatch() ) {
419  if (vtxeval::ask_user()) return;
420  }
421  else { vtxeval::gDebugFlag = false; }
422  }
423 
424  fOut->Write();
425 }
426 
427 
428 //+++++++++++++++++++++ from Jonathan +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
429 /**********************************************************************************
430  * Project : TMVA - a ROOT-integrated toolkit for multivariate data analysis *
431  * Package : TMVA *
432  * Root Macro: TMVAClassification *
433  * *
434  * This macro provides examples for the training and testing of the *
435  * TMVA classifiers. *
436  * *
437  * As input data is used a toy-MC sample consisting of four Gaussian-distributed *
438  * and linearly correlated input variables. *
439  * *
440  * The methods to be used can be switched on and off by means of booleans, or *
441  * via the prompt command, for example: *
442  * *
443  * root -l ./TMVAClassification.C\(\"Fisher,Likelihood\"\) *
444  * *
445  * (note that the backslashes are mandatory) *
446  * If no method given, a default set of classifiers is used. *
447  * *
448  * The output file "TMVA.root" can be analysed with the use of dedicated *
449  * macros (simply say: root -l <macro.C>), which can be conveniently *
450  * invoked through a GUI that will appear at the end of the run of this macro. *
451  * Launch the GUI via the command: *
452  * *
453  * root -l ./TMVAGui.C *
454  * *
455  **********************************************************************************/
456 void TMVAClassification( TString myMethodList = "")
457 {
458  TTree *signal = (TTree *)gDirectory->Get("VertexG");
459 
460  if (! signal) { std::cout << "No signal TTree" << std::endl; return;}
461 
462  TTree *background = (TTree *)gDirectory->Get("VertexB");
463 
464  if (! background) { std::cout << "No background TTree" << std::endl; return;}
465 
466  //---------------------------------------------------------------
467  // This loads the library
468  TMVA::Tools::Instance();
469 
470  // Default MVA methods to be trained + tested
471  std::map<std::string, int> Use;
472 
473  // --- Cut optimisation
474  Use["Cuts"] = 1;
475  Use["CutsD"] = 1;
476  Use["CutsPCA"] = 0;
477  Use["CutsGA"] = 0;
478  Use["CutsSA"] = 0;
479  //
480  // --- 1-dimensional likelihood ("naive Bayes estimator")
481  Use["Likelihood"] = 1;
482  Use["LikelihoodD"] = 0; // the "D" extension indicates decorrelated input variables (see option strings)
483  Use["LikelihoodPCA"] = 1; // the "PCA" extension indicates PCA-transformed input variables (see option strings)
484  Use["LikelihoodKDE"] = 0;
485  Use["LikelihoodMIX"] = 0;
486  //
487  // --- Mutidimensional likelihood and Nearest-Neighbour methods
488  Use["PDERS"] = 1;
489  Use["PDERSD"] = 0;
490  Use["PDERSPCA"] = 0;
491  Use["PDEFoam"] = 1;
492  Use["PDEFoamBoost"] = 0; // uses generalised MVA method boosting
493  Use["KNN"] = 1; // k-nearest neighbour method
494  //
495  // --- Linear Discriminant Analysis
496  Use["LD"] = 1; // Linear Discriminant identical to Fisher
497  Use["Fisher"] = 0;
498  Use["FisherG"] = 0;
499  Use["BoostedFisher"] = 0; // uses generalised MVA method boosting
500  Use["HMatrix"] = 0;
501  //
502  // --- Function Discriminant analysis
503  Use["FDA_GA"] = 1; // minimisation of user-defined function using Genetics Algorithm
504  Use["FDA_SA"] = 0;
505  Use["FDA_MC"] = 0;
506  Use["FDA_MT"] = 0;
507  Use["FDA_GAMT"] = 0;
508  Use["FDA_MCMT"] = 0;
509  //
510  // --- Neural Networks (all are feed-forward Multilayer Perceptrons)
511  Use["MLP"] = 0; // Recommended ANN
512  Use["MLPBFGS"] = 0; // Recommended ANN with optional training method
513  Use["MLPBNN"] = 1; // Recommended ANN with BFGS training method and bayesian regulator
514  Use["CFMlpANN"] = 0; // Depreciated ANN from ALEPH
515  Use["TMlpANN"] = 0; // ROOT's own ANN
516  //
517  // --- Support Vector Machine
518  Use["SVM"] = 1;
519  //
520  // --- Boosted Decision Trees
521  Use["BDT"] = 1; // uses Adaptive Boost
522  Use["BDTG"] = 0; // uses Gradient Boost
523  Use["BDTB"] = 0; // uses Bagging
524  Use["BDTD"] = 0; // decorrelation + Adaptive Boost
525  Use["BDTF"] = 0; // allow usage of fisher discriminant for node splitting
526  Use["myBDTD"] = 1; // mine
527  //
528  // --- Friedman's RuleFit method, ie, an optimised series of cuts ("rules")
529  Use["RuleFit"] = 1;
530  // ---------------------------------------------------------------
531 
532  std::cout << std::endl;
533  std::cout << "==> Start TMVAClassification" << std::endl;
534 
535  // Select methods (don't look at this code - not of interest)
536  if (myMethodList != "") {
537  for (std::map<std::string, int>::iterator it = Use.begin(); it != Use.end(); it++) it->second = 0;
538 
539  std::vector<TString> mlist = TMVA::gTools().SplitString( myMethodList, ',' );
540 
541  for (size_t i = 0; i < mlist.size(); i++) {
542  std::string regMethod(mlist[i]);
543 
544  if (Use.find(regMethod) == Use.end()) {
545  std::cout << "Method \"" << regMethod << "\" not known in TMVA under this name. Choose among the following:" << std::endl;
546 
547  for (std::map<std::string, int>::iterator it = Use.begin(); it != Use.end(); it++) std::cout << it->first << " ";
548 
549  std::cout << std::endl;
550  return;
551  }
552 
553  Use[regMethod] = 1;
554  }
555  }
556 
557  // --------------------------------------------------------------------------------------------------
558 
559  // --- Here the preparation phase begins
560 
561  // Create a ROOT output file where TMVA will store ntuples, histograms, etc.
562  TString outfileName( "TMVA.root" );
563  TFile *outputFile = TFile::Open( outfileName, "RECREATE" );
564 
565  // Create the factory object. Later you can choose the methods
566  // whose performance you'd like to investigate. The factory is
567  // the only TMVA object you have to interact with
568  //
569  // The first argument is the base of the name of all the
570  // weightfiles in the directory weight/
571  //
572  // The second argument is the output file for the training results
573  // All TMVA output can be suppressed by removing the "!" (not) in
574  // front of the "Silent" argument in the option string
575  TMVA::Factory *factory = new TMVA::Factory( "TMVAClassification", outputFile,
576  "!V:!Silent:Color:DrawProgressBar:Transformations=I;D;P;G,D:AnalysisType=Classification" );
577 
578  // If you wish to modify default settings
579  // (please check "src/Config.h" to see all available global options)
580  // (TMVA::gConfig().GetVariablePlotting()).fTimesRMS = 8.0;
581  // (TMVA::gConfig().GetIONames()).fWeightFileDir = "myWeightDirectory";
582 
583  // load the signal and background event samples from ROOT trees
584 
585  std::cout << " starts ... " << std::endl;
586  // global event weights per tree (see below for setting event-wise weights)
587  // Float_t w;
588  double signalWeight = 1.0;
589  double backgroundWeight = 1.0;
590 
591  std::cout << " signalWeight = " << signalWeight << " backWeight = " << backgroundWeight << std::endl;
592  factory->AddSignalTree( signal, signalWeight );
593  factory->AddBackgroundTree( background, backgroundWeight );
594 
595  TString separator(":");
596  TString Vnames(vnames);
597  TObjArray *array = Vnames.Tokenize(separator);
598 
599  std::vector<std::string> inputVars;
600  TIter next(array);
601  TObjString *objs;
602 
603  while ((objs = (TObjString *) next())) {
604  // std::cout << objs->GetString() << std::endl;
605  TString name(objs->GetString());
606 
607  if (name == "BEMC") continue;
608 
609  if (name == "noBEMC") continue;
610 
611  factory->AddVariable(name, 'F');
612  }
613 
614  // This would set individual event weights (the variables defined in the
615  // expression need to exist in the original TTree)
616  // for signal : factory->SetSignalWeightExpression("weight1*weight2");
617  // for background: factory->SetBackgroundWeightExpression("weight1*weight2");
618  // commented JB : 04/26 ??
619  //factory->dSetBackgroundWeightExpression("weight");
620 
621  // Apply additional cuts on the signal and background samples (can be different)
622  TCut mycuts = "";
623  TCut mycutb = "";
624 
625  // Tell the factory how to use the training and testing events
626  //
627  // If no numbers of events are given, half of the events in the tree are used
628  // for training, and the other half for testing:
629  // factory->PrepareTrainingAndTestTree( mycut, "SplitMode=random:!V" );
630  // To also specify the number of testing events, use:
631  //factory->PrepareTrainingAndTestTree( mycuts,mycutb,"NSigTrain=9000:NBkgTrain=50000:NSigTest=9000:NBkgTest=50000:SplitMode=Random:!V" );
632  factory->PrepareTrainingAndTestTree( mycuts, mycutb, "nTrain_Signal=4900:nTrain_Background=49000:nTest_Signal=4900:nTest_Background=49000:SplitMode=Random:!V"); // for KFVertex
633  // factory->PrepareTrainingAndTestTree( mycuts, mycutb,"nTrain_Signal=20000:nTrain_Background=40000:nTest_Signal=20000:nTest_Background=40000:SplitMode=Random:!V"); // for PPV
634 
635  // ---- Book MVA methods
636  //
637  // Please lookup the various method configuration options in the corresponding cxx files, eg:
638  // src/MethoCuts.cxx, etc, or here: http://tmva.sourceforge.net/optionRef.html
639  // it is possible to preset ranges in the option string in which the cut optimisation should be done:
640  // "...:CutRangeMin[2]=-1:CutRangeMax[2]=1"...", where [2] is the third input variable
641 
642  // Cut optimisation
643  if (Use["Cuts"])
644  factory->BookMethod( TMVA::Types::kCuts, "Cuts",
645  "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart" );
646 
647  if (Use["CutsD"])
648  factory->BookMethod( TMVA::Types::kCuts, "CutsD",
649  "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart:VarTransform=Decorrelate" );
650 
651  if (Use["CutsPCA"])
652  factory->BookMethod( TMVA::Types::kCuts, "CutsPCA",
653  "!H:!V:FitMethod=MC:EffSel:SampleSize=200000:VarProp=FSmart:VarTransform=PCA" );
654 
655  if (Use["CutsGA"])
656  factory->BookMethod( TMVA::Types::kCuts, "CutsGA",
657  "H:!V:FitMethod=GA:CutRangeMin[0]=-10:CutRangeMax[0]=10:VarProp[1]=FMax:EffSel:Steps=30:Cycles=3:PopSize=400:SC_steps=10:SC_rate=5:SC_factor=0.95" );
658 
659  if (Use["CutsSA"])
660  factory->BookMethod( TMVA::Types::kCuts, "CutsSA",
661  "!H:!V:FitMethod=SA:EffSel:MaxCalls=150000:KernelTemp=IncAdaptive:InitialTemp=1e+6:MinTemp=1e-6:Eps=1e-10:UseDefaultScale" );
662 
663  // Likelihood ("naive Bayes estimator")
664  if (Use["Likelihood"])
665  factory->BookMethod( TMVA::Types::kLikelihood, "Likelihood",
666  "H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmoothBkg[1]=10:NSmooth=1:NAvEvtPerBin=50" );
667 
668  // Decorrelated likelihood
669  if (Use["LikelihoodD"])
670  factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodD",
671  "!H:!V:TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmooth=5:NAvEvtPerBin=50:VarTransform=Decorrelate" );
672 
673  // PCA-transformed likelihood
674  if (Use["LikelihoodPCA"])
675  factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodPCA",
676  "!H:!V:!TransformOutput:PDFInterpol=Spline2:NSmoothSig[0]=20:NSmoothBkg[0]=20:NSmooth=5:NAvEvtPerBin=50:VarTransform=PCA" );
677 
678  // Use a kernel density estimator to approximate the PDFs
679  if (Use["LikelihoodKDE"])
680  factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodKDE",
681  "!H:!V:!TransformOutput:PDFInterpol=KDE:KDEtype=Gauss:KDEiter=Adaptive:KDEFineFactor=0.3:KDEborder=None:NAvEvtPerBin=50" );
682 
683  // Use a variable-dependent mix of splines and kernel density estimator
684  if (Use["LikelihoodMIX"])
685  factory->BookMethod( TMVA::Types::kLikelihood, "LikelihoodMIX",
686  "!H:!V:!TransformOutput:PDFInterpolSig[0]=KDE:PDFInterpolBkg[0]=KDE:PDFInterpolSig[1]=KDE:PDFInterpolBkg[1]=KDE:PDFInterpolSig[2]=Spline2:PDFInterpolBkg[2]=Spline2:PDFInterpolSig[3]=Spline2:PDFInterpolBkg[3]=Spline2:KDEtype=Gauss:KDEiter=Nonadaptive:KDEborder=None:NAvEvtPerBin=50" );
687 
688  // Test the multi-dimensional probability density estimator
689  // here are the options strings for the MinMax and RMS methods, respectively:
690  // "!H:!V:VolumeRangeMode=MinMax:DeltaFrac=0.2:KernelEstimator=Gauss:GaussSigma=0.3" );
691  // "!H:!V:VolumeRangeMode=RMS:DeltaFrac=3:KernelEstimator=Gauss:GaussSigma=0.3" );
692  if (Use["PDERS"])
693  factory->BookMethod( TMVA::Types::kPDERS, "PDERS",
694  "!H:!V:NormTree=T:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600" );
695 
696  if (Use["PDERSD"])
697  factory->BookMethod( TMVA::Types::kPDERS, "PDERSD",
698  "!H:!V:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600:VarTransform=Decorrelate" );
699 
700  if (Use["PDERSPCA"])
701  factory->BookMethod( TMVA::Types::kPDERS, "PDERSPCA",
702  "!H:!V:VolumeRangeMode=Adaptive:KernelEstimator=Gauss:GaussSigma=0.3:NEventsMin=400:NEventsMax=600:VarTransform=PCA" );
703 
704  // Multi-dimensional likelihood estimator using self-adapting phase-space binning
705  if (Use["PDEFoam"])
706  factory->BookMethod( TMVA::Types::kPDEFoam, "PDEFoam",
707  "!H:!V:SigBgSeparate=F:TailCut=0.001:VolFrac=0.0666:nActiveCells=500:nSampl=2000:nBin=5:Nmin=100:Kernel=None:Compress=T" );
708 
709  if (Use["PDEFoamBoost"])
710  factory->BookMethod( TMVA::Types::kPDEFoam, "PDEFoamBoost",
711  "!H:!V:Boost_Num=30:Boost_Transform=linear:SigBgSeparate=F:MaxDepth=4:UseYesNoCell=T:DTLogic=MisClassificationError:FillFoamWithOrigWeights=F:TailCut=0:nActiveCells=500:nBin=20:Nmin=400:Kernel=None:Compress=T" );
712 
713  // K-Nearest Neighbour classifier (KNN)
714  if (Use["KNN"])
715  factory->BookMethod( TMVA::Types::kKNN, "KNN",
716  "H:nkNN=20:ScaleFrac=0.8:SigmaFact=1.0:Kernel=Gaus:UseKernel=F:UseWeight=T:!Trim" );
717 
718  // H-Matrix (chi2-squared) method
719  if (Use["HMatrix"])
720  factory->BookMethod( TMVA::Types::kHMatrix, "HMatrix", "!H:!V:VarTransform=None" );
721 
722  // Linear discriminant (same as Fisher discriminant)
723  if (Use["LD"])
724  factory->BookMethod( TMVA::Types::kLD, "LD", "H:!V:VarTransform=None:CreateMVAPdfs:PDFInterpolMVAPdf=Spline2:NbinsMVAPdf=50:NsmoothMVAPdf=10" );
725 
726  // Fisher discriminant (same as LD)
727  if (Use["Fisher"])
728  factory->BookMethod( TMVA::Types::kFisher, "Fisher", "H:!V:Fisher:VarTransform=None:CreateMVAPdfs:PDFInterpolMVAPdf=Spline2:NbinsMVAPdf=50:NsmoothMVAPdf=10" );
729 
730  // Fisher with Gauss-transformed input variables
731  if (Use["FisherG"])
732  factory->BookMethod( TMVA::Types::kFisher, "FisherG", "H:!V:VarTransform=Gauss" );
733 
734  // Composite classifier: ensemble (tree) of boosted Fisher classifiers
735  if (Use["BoostedFisher"])
736  factory->BookMethod( TMVA::Types::kFisher, "BoostedFisher",
737  "H:!V:Boost_Num=20:Boost_Transform=log:Boost_Type=AdaBoost:Boost_AdaBoostBeta=0.2:!Boost_DetailedMonitoring" );
738 
739  // Function discrimination analysis (FDA) -- test of various fitters - the recommended one is Minuit (or GA or SA)
740  if (Use["FDA_MC"])
741  factory->BookMethod( TMVA::Types::kFDA, "FDA_MC",
742  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=MC:SampleSize=100000:Sigma=0.1" );
743 
744  if (Use["FDA_GA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options])
745  factory->BookMethod( TMVA::Types::kFDA, "FDA_GA",
746  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=GA:PopSize=300:Cycles=3:Steps=20:Trim=True:SaveBestGen=1" );
747 
748  if (Use["FDA_SA"]) // can also use Simulated Annealing (SA) algorithm (see Cuts_SA options])
749  factory->BookMethod( TMVA::Types::kFDA, "FDA_SA",
750  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=SA:MaxCalls=15000:KernelTemp=IncAdaptive:InitialTemp=1e+6:MinTemp=1e-6:Eps=1e-10:UseDefaultScale" );
751 
752  if (Use["FDA_MT"])
753  factory->BookMethod( TMVA::Types::kFDA, "FDA_MT",
754  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=MINUIT:ErrorLevel=1:PrintLevel=-1:FitStrategy=2:UseImprove:UseMinos:SetBatch" );
755 
756  if (Use["FDA_GAMT"])
757  factory->BookMethod( TMVA::Types::kFDA, "FDA_GAMT",
758  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=GA:Converger=MINUIT:ErrorLevel=1:PrintLevel=-1:FitStrategy=0:!UseImprove:!UseMinos:SetBatch:Cycles=1:PopSize=5:Steps=5:Trim" );
759 
760  if (Use["FDA_MCMT"])
761  factory->BookMethod( TMVA::Types::kFDA, "FDA_MCMT",
762  "H:!V:Formula=(0)+(1)*x0+(2)*x1+(3)*x2+(4)*x3:ParRanges=(-1,1);(-10,10);(-10,10);(-10,10);(-10,10):FitMethod=MC:Converger=MINUIT:ErrorLevel=1:PrintLevel=-1:FitStrategy=0:!UseImprove:!UseMinos:SetBatch:SampleSize=20" );
763 
764  // TMVA ANN: MLP (recommended ANN) -- all ANNs in TMVA are Multilayer Perceptrons
765  if (Use["MLP"])
766  factory->BookMethod( TMVA::Types::kMLP, "MLP", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:!UseRegulator" );
767 
768  if (Use["MLPBFGS"])
769  factory->BookMethod( TMVA::Types::kMLP, "MLPBFGS", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:TrainingMethod=BFGS:!UseRegulator" );
770 
771  if (Use["MLPBNN"])
772  factory->BookMethod( TMVA::Types::kMLP, "MLPBNN", "H:!V:NeuronType=tanh:VarTransform=N:NCycles=600:HiddenLayers=N+5:TestRate=5:TrainingMethod=BFGS:UseRegulator" ); // BFGS training with bayesian regulators
773 
774  // CF(Clermont-Ferrand)ANN
775  if (Use["CFMlpANN"])
776  factory->BookMethod( TMVA::Types::kCFMlpANN, "CFMlpANN", "!H:!V:NCycles=2000:HiddenLayers=N+1,N" ); // n_cycles:#nodes:#nodes:...
777 
778  // Tmlp(Root)ANN
779  if (Use["TMlpANN"])
780  factory->BookMethod( TMVA::Types::kTMlpANN, "TMlpANN", "!H:!V:NCycles=200:HiddenLayers=N+1,N:LearningMethod=BFGS:ValidationFraction=0.3" ); // n_cycles:#nodes:#nodes:...
781 
782  // Support Vector Machine
783  if (Use["SVM"])
784  factory->BookMethod( TMVA::Types::kSVM, "SVM", "Gamma=0.25:Tol=0.001:VarTransform=Norm" );
785 
786  // Boosted Decision Trees
787  if (Use["BDTG"]) // Gradient Boost
788  factory->BookMethod( TMVA::Types::kBDT, "BDTG",
789  "!H:!V:NTrees=1000:BoostType=Grad:Shrinkage=0.10:UseBaggedGrad:GradBaggingFraction=0.5:nCuts=20:NNodesMax=5" );
790 
791  if (Use["BDT"]) // Adaptive Boost
792  factory->BookMethod( TMVA::Types::kBDT, "BDT",
793  "!H:!V:NTrees=850:nEventsMin=150:MaxDepth=3:BoostType=AdaBoost:AdaBoostBeta=0.5:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning" );
794 
795 
796  if (Use["BDTB"]) // Bagging
797  factory->BookMethod( TMVA::Types::kBDT, "BDTB",
798  "!H:!V:NTrees=400:BoostType=Bagging:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning" );
799 
800  if (Use["BDTD"]) // Decorrelation + Adaptive Boost
801  factory->BookMethod( TMVA::Types::kBDT, "BDTD",
802  "!H:!V:NTrees=400:nEventsMin=400:MaxDepth=3:BoostType=AdaBoost:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning:VarTransform=Decorrelate" );
803 
804  if (Use["myBDTD"]) // Decorrelation + Adaptive Boost
805  factory->BookMethod( TMVA::Types::kBDT, "BDTDTEST",
806  "!H:!V:NTrees=1000:nEventsMin=400:MaxDepth=6:BoostType=AdaBoost:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning:VarTransform=Decorrelate" );
807 
808  if (Use["BDTF"]) // Allow Using Fisher discriminant in node splitting for (strong) linearly correlated variables
809  factory->BookMethod( TMVA::Types::kBDT, "BDTMitFisher",
810  "!H:!V:NTrees=50:nEventsMin=150:UseFisherCuts:MaxDepth=3:BoostType=AdaBoost:AdaBoostBeta=0.5:SeparationType=GiniIndex:nCuts=20:PruneMethod=NoPruning" );
811 
812  // RuleFit -- TMVA implementation of Friedman's method
813  if (Use["RuleFit"])
814  factory->BookMethod( TMVA::Types::kRuleFit, "RuleFit",
815  "H:!V:RuleFitModule=RFTMVA:Model=ModRuleLinear:MinImp=0.001:RuleMinDist=0.001:NTrees=20:fEventsMin=0.01:fEventsMax=0.5:GDTau=-1.0:GDTauPrec=0.01:GDStep=0.01:GDNSteps=10000:GDErrScale=1.02" );
816 
817  // For an example of the category classifier usage, see: TMVAClassificationCategory
818 
819  // TMVA::IMethod* category = factory->BookMethod( TMVA::Types::kCategory,"Category","" );
820 
821  // --------------------------------------------------------------------------------------------------
822 
823  // ---- Now you can optimize the setting (configuration) of the MVAs using the set of training events
824 #if 0
825  factory->OptimizeAllMethods("SigEffAt001", "Scan");
826  factory->OptimizeAllMethods("ROCIntegral", "GA");
827 #endif
828  // --------------------------------------------------------------------------------------------------
829 
830  // ---- Now you can tell the factory to train, test, and evaluate the MVAs
831 
832  // Train MVAs using the set of training events
833  factory->TrainAllMethods();
834 
835  // ---- Evaluate all MVAs using the set of test events
836  factory->TestAllMethods();
837 
838  // ----- Evaluate and compare performance of all configured MVAs
839  factory->EvaluateAllMethods();
840 
841  // --------------------------------------------------------------
842 
843  // Save the output
844  outputFile->Close();
845 
846  std::cout << "==> Wrote root file: " << outputFile->GetName() << std::endl;
847  std::cout << "==> TMVAClassification is done!" << std::endl;
848 
849  delete factory;
850 }
Bool_t Accept(const StMuTrack *gTrack=0)
Bool_t AcceptVX(const StMuPrimaryVertex *Vtx=0)
Float_t cross
void ForceAnimate(unsigned int times=0, int msecDelay=0)
Float_t EEMC
void MuMcPrVKFV2012(Long64_t nevent, const char *file, const std::string &outFile, bool fillNtuple)
Float_t noEEMC
void FillData(data_t &data, const StMuPrimaryVertex *Vtx)
Float_t prompt
Float_t tof
Float_t chi2
void TMVAClassification(TString myMethodList="")
Float_t postx
Float_t notof
const Char_t * vnames
Float_t beam