source: trunk/MagicSoft/Mars/mtemp/mpadova/macros/RanForestPD.C@ 6547

Last change on this file since 6547 was 5665, checked in by MAGIC, 20 years ago
*** empty log message ***
File size: 8.3 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Hengstebeck 3/2003 <mailto:hengsteb@alwa02.physik.uni-siegen.de>!
19!
20! Copyright: MAGIC Software Development, 2000-2003
21!
22!
23\* ======================================================================== */
24
25#include "MMcEvt.hxx"
26
27void RanForestPD(Int_t minsize=0)
28{
29 //
30 // Create a empty Parameter List and an empty Task List
31 // The tasklist is identified in the eventloop by its name
32 //
33 MParList plist;
34
35 MTaskList tlist;
36 plist.AddToList(&tlist);
37
38 MReadMarsFile read("Events", "star_gamma_train.root");
39 Float_t numgammas = read.GetEntries();
40
41 read.AddFile("star_proton_train.root");
42 read.AddFile("star_helium_train.root");
43
44 Float_t numhadrons = read.GetEntries() - numgammas;
45
46 // Fraction of gammas to be used in training:
47 // (less than total sample, to speed up execution)
48 //
49 Float_t gamma_frac = 2.*(numhadrons / numgammas);
50
51
52 read.DisableAutoScheme();
53
54 tlist.AddToList(&read);
55
56 MFParticleId fhadrons("MMcEvt", '!', MMcEvt::kGAMMA);
57 tlist.AddToList(&fhadrons);
58
59 MHMatrix matrixg("MatrixGammas");
60
61 // matrixg.AddColumn("floor(0.5+(cos(MMcEvt.fTelescopeTheta)/0.01))");
62 // matrixg.AddColumn("MMcEvt.fTelescopePhi");
63 // matrixg.AddColumn("MSigmabar.fSigmabar");
64
65 matrixg.AddColumn("MHillas.fSize");
66 matrixg.AddColumn("MHillasSrc.fDist");
67
68 // matrixg.AddColumn("MHillas.fWidth");
69 // matrixg.AddColumn("MHillas.fLength");
70
71 // Scaled Width and Length:
72 //
73 matrixg.AddColumn("MHillas.fWidth/(-13.1618+(10.0492*log10(MHillas.fSize)))");
74 matrixg.AddColumn("MHillas.fLength/(-57.3784+(32.6131*log10(MHillas.fSize)))");
75
76 // matrixg.AddColumn("(MHillas.fSize/(MHillas.fWidth*MHillas.fLength))");
77 matrixg.AddColumn("MNewImagePar.fConc");
78
79 // matrixg.AddColumn("MNewImagePar.fConc1");
80 // matrixg.AddColumn("MNewImagePar.fAsym");
81 // matrixg.AddColumn("MHillasExt.fM3Trans");
82 // matrixg.AddColumn("abs(MHillasSrc.fHeadTail)");
83 // matrixg.AddColumn("MNewImagePar.fLeakage1");
84
85 //
86 // Third moment, with sign referred to source position:
87 //
88 matrixg.AddColumn("MHillasExt.fM3Long*MHillasSrc.fCosDeltaAlpha/abs(MHillasSrc.fCosDeltaAlpha)");
89
90 //# matrixg.AddColumn("MHillasSrc.fAlpha");
91
92 //
93 // SIZE cut:
94 //
95 TString dummy("(MHillas.fSize<");
96 dummy += minsize;
97 dummy += ")";
98
99 // AM, useful FOR High Energy only:
100 // dummy += " || (MImagePar.fNumSatPixelsLG>0)";
101
102 MContinue SizeCut(dummy);
103 tlist.AddToList(&SizeCut);
104
105 // AM TEST!!!, no muons in train
106 // dummy = "(MMcEvt.fMuonCphFraction>0.1)";
107 // MContinue MuonsCut(dummy);
108 // tlist.AddToList(&MuonsCut);
109
110 plist.AddToList(&matrixg);
111
112 MHMatrix matrixh("MatrixHadrons");
113 matrixh.AddColumns(matrixg.GetColumns());
114 plist.AddToList(&matrixh);
115
116 MFillH fillmath("MatrixHadrons");
117 fillmath.SetFilter(&fhadrons);
118 tlist.AddToList(&fillmath);
119
120 MContinue skiphad(&fhadrons);
121 tlist.AddToList(&skiphad);
122
123 MFEventSelector reduce_training_gammas;
124 reduce_training_gammas.SetSelectionRatio(gamma_frac);
125 tlist.AddToList(&reduce_training_gammas);
126
127 MFillH fillmatg("MatrixGammas");
128 fillmatg.SetFilter(&reduce_training_gammas);
129 tlist.AddToList(&fillmatg);
130
131 //
132 // Create and setup the eventloop
133 //
134 MEvtLoop evtloop;
135 evtloop.SetParList(&plist);
136
137 MProgressBar bar;
138 evtloop.SetProgressBar(&bar);
139 //
140 // Execute your analysis
141 //
142 if (!evtloop.Eventloop())
143 return;
144
145 tlist.PrintStatistics();
146
147 // ---------------------------------------------------------------
148 // ---------------------------------------------------------------
149 // second event loop: the trees of the random forest are grown,
150 // the event loop is now actually a tree loop (loop of training
151 // process)
152 // ---------------------------------------------------------------
153 // ---------------------------------------------------------------
154
155 MTaskList tlist2;
156 plist.Replace(&tlist2);
157
158 MRanForestGrow rfgrow2;
159 rfgrow2.SetNumTrees(100);
160 rfgrow2.SetNumTry(3);
161 rfgrow2.SetNdSize(10);
162
163 MWriteRootFile rfwrite2("RF.root");
164 rfwrite2.AddContainer("MRanTree","Tree"); //save all trees
165 MFillH fillh2("MHRanForestGini");
166
167 tlist2.AddToList(&rfgrow2);
168 tlist2.AddToList(&rfwrite2);
169 tlist2.AddToList(&fillh2);
170
171 // gRandom is accessed from MRanForest (-> bootstrap aggregating)
172 // and MRanTree (-> random split selection) and should be initialized
173 // here if you want to set a certain random number generator
174 if(gRandom)
175 delete gRandom;
176 // gRandom = new TRandom3(0); // Takes seed from the computer clock
177 gRandom = new TRandom3(); // Uses always same seed
178 //
179 // Execute tree growing (now the eventloop is actually a treeloop)
180 //
181
182 evtloop.SetProgressBar(&bar);
183
184 if (!evtloop.Eventloop())
185 return;
186
187 tlist2.PrintStatistics();
188
189 plist.FindObject("MHRanForestGini")->DrawClone();
190
191 // ---------------------------------------------------------------
192 // ---------------------------------------------------------------
193 // third event loop: the control sample (star2.root) is processed
194 // through the previously grown random forest,
195 //
196 // the histograms MHHadronness (quality of g/h-separation) and
197 // MHRanForest are created and displayed.
198 // MHRanForest shows the convergence of the classification error
199 // as function of the number of grown (and combined) trees
200 // and tells the user how many trees are actually needed in later
201 // classification tasks.
202 // ---------------------------------------------------------------
203 // ---------------------------------------------------------------
204
205 MTaskList tlist3;
206
207 plist.Replace(&tlist3);
208
209 MReadMarsFile read3("Events", "star_gamma_test.root");
210
211 read3.AddFile("star_proton_test.root");
212 read3.AddFile("star_helium_test.root");
213
214 read3.DisableAutoScheme();
215 tlist3.AddToList(&read3);
216
217 tlist3.AddToList(&SizeCut);
218
219 MRanForestCalc calc;
220 tlist3.AddToList(&calc);
221
222 // parameter containers you want to save
223 //
224
225 TString outfname = "star_S";
226 outfname += minsize;
227 outfname += "_all.root";
228 MWriteRootFile write(outfname, "recreate");
229
230 write.AddContainer("MMcEvt", "Events");
231 write.AddContainer("MHillas", "Events");
232 write.AddContainer("MHillasExt", "Events");
233 write.AddContainer("MImagePar", "Events");
234 write.AddContainer("MNewImagePar", "Events");
235 write.AddContainer("MHillasSrc", "Events");
236 write.AddContainer("MHadronness", "Events");
237
238 write.AddContainer("MRawRunHeader", "RunHeaders");
239 write.AddContainer("MSrcPosCam", "RunHeaders");
240 write.AddContainer("MMcRunHeader", "RunHeaders");
241
242 /*
243 write.AddContainer("MMcTrig", "Events");
244 write.AddContainer("MRawEvtData", "Events");
245 write.AddContainer("MRawEvtHeader", "Events");
246 */
247
248
249 MFillH fillh3a("MHHadronness");
250 MFillH fillh3b("MHRanForest");
251
252 tlist3.AddToList(&fillh3a);
253 tlist3.AddToList(&fillh3b);
254
255
256 tlist3.AddToList(&write);
257
258 evtloop.SetProgressBar(&bar);
259
260 //
261 // Execute your analysis
262 //
263 if (!evtloop.Eventloop())
264 return;
265
266 tlist3.PrintStatistics();
267
268 plist.FindObject("MHRanForest")->DrawClone();
269
270 TCanvas *had = new TCanvas;
271 had->cd();
272 plist.FindObject("MHHadronness")->Draw();
273
274 TString gifname = "had_S";
275 gifname += minsize;
276 gifname += ".gif";
277
278 had->SaveAs(gifname);
279 had->Close();
280
281 plist.FindObject("MHHadronness")->DrawClone();
282 plist.FindObject("MHHadronness")->Print();//*/
283
284 return;
285}
Note: See TracBrowser for help on using the repository browser.