source: tags/Mars-V0.9.5/mtools/MTFillMatrix.cc

Last change on this file was 7413, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 11.7 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 Bretz, 12/2003 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2005
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MTFillMatrix
28//
29// Use this tool to fill eg trainings and test-matrices, while the matrix
30// to be filled can be a real matrix (MHMatrix) or a file (MWriteRootFile)
31// or both.
32//
33// First create a reference histogram (MH3). For more details see
34// MFEventSelector2 which is used to select events according to the
35// reference histogram.
36//
37// If no reference histogram is available the number of events are
38// randomly choosen from the sample with a probability which fits
39// the total destination number of events.
40//
41// Here is an example of how to choose 1000 events somehow distributed in
42// size from a file.
43// -----------------------------------------------------------------------
44//
45// MH3 ref("MHillas.fSize"); // choose a predefined size distribution
46// // Now setup the distribution
47//
48// MHMatrix matrix1; // create matrix to fill
49// matrix.AddColumn("MHillas.fWidth"); // setup columns of your matrix
50//
51// MReadMarsFile read("myfile.root"); // Setup your 'reader'
52// read.DisableAutoScheme(); // make sure everything is read
53//
54// MTFillMatrix fill(&ref); // setup MTFillMatrix
55// fill.SetNumDestEvents1(1000); // setup number of events to select
56// fill.SetDestMatrix1(&matrix1); // setup destination matrix
57// if (!fill.Process()) // check if everything worked
58// return;
59// fill.WriteMatrix1("myoutput.root"); // write matrix to file
60//
61//
62// To get two matrices instead of one (splitted randomly) you can add
63// the following before calling Process():
64// ------------------------------------------------------------------------
65//
66// MHMatrix matrix2;
67// fill.SetNumDestEvents2(500); // setup number of events to select
68// fill.SetDestMatrix2(&matrix2); // setup destination matrix
69// [...]
70// fill.WriteMatrix2("myoutput2.root");
71//
72//
73// To write both matrices into a single file use:
74// ----------------------------------------------
75//
76// fill.WriteMatrices("myfile.root");
77//
78/////////////////////////////////////////////////////////////////////////////
79#include "MTFillMatrix.h"
80
81#include <TFile.h>
82
83// environment
84#include "MHMatrix.h"
85
86#include "MLog.h"
87#include "MLogManip.h"
88
89// eventloop
90#include "MParList.h"
91#include "MTaskList.h"
92#include "MEvtLoop.h"
93
94// tasks
95#include "MRead.h"
96#include "MFillH.h"
97#include "MContinue.h"
98
99// filters
100#include "MF.h"
101#include "MFilterList.h"
102#include "MFEventSelector.h"
103#include "MFEventSelector2.h"
104
105ClassImp(MTFillMatrix);
106
107using namespace std;
108
109// --------------------------------------------------------------------------
110//
111// Print Size and contained columns.
112// Check whether the number of generated events is compatible with
113// the number of requested events.
114//
115Bool_t MTFillMatrix::CheckResult(MHMatrix *m, Int_t num) const
116{
117 if (!m || num==0)
118 return kTRUE;
119
120 m->Print("SizeCols");
121
122 const Int_t generated = m->GetM().GetNrows();
123 if (TMath::Abs(generated-num) <= TMath::Sqrt(9.0*num))
124 return kTRUE;
125
126 *fLog << warn << "WARNING - No. of generated events (";
127 *fLog << generated << ") incompatible with requested events (";
128 *fLog << num << ") for " << m->GetDescriptor() << endl;
129
130 return kFALSE;
131}
132
133// --------------------------------------------------------------------------
134//
135// Write the given MHMatrix with its name as default name to a
136// file fname.
137//
138Bool_t MTFillMatrix::WriteMatrix(MHMatrix *m, const TString &fname, Int_t i) const
139{
140 if (!m)
141 {
142 *fLog << "ERROR - Unable to write matrix #" << i << " (=NULL)... ignored." << endl;
143 return kFALSE;
144 }
145 if (fname.IsNull())
146 {
147 *fLog << "ERROR - Unable to write matrix, file name empty." << endl;
148 return kFALSE;
149 }
150
151 TFile file(fname, "RECREATE", m->GetTitle());
152 m->Write();
153 return kTRUE;
154}
155
156// --------------------------------------------------------------------------
157//
158// Constructor. Takes an MH3 as argument. This MH3 is the reference
159// distribution to fill the matrix. More information can be found
160// at MFEventSelector2 which is used to select the events.
161//
162// If no MH3 *ref is given the events are randomly selected from the
163// total sample - this may result in samples which don't have exactly
164// the predefined size, but it is much faster.
165//
166MTFillMatrix::MTFillMatrix(const MH3 *ref)
167: fReference(0), fReader(0), fDestMatrix1(0), fDestMatrix2(0),
168 fNumDestEvents1(0), fNumDestEvents2(0), fWriteFile1(0), fWriteFile2(0)
169{
170 fName = "MFillMatrix";
171 fTitle = "Tool to fill MHMatrix from file";
172
173 if (ref)
174 fReference = (MH3*)ref->Clone();
175}
176
177MTFillMatrix::~MTFillMatrix()
178{
179 if (fReference)
180 delete fReference;
181}
182
183//------------------------------------------------------------------------
184//
185// Add a cut which is used to fill the matrix, eg "MMcEvt.fOartId<1.5"
186// (The rule is applied, nit inverted: The matrix is filled with
187// the events fullfilling the condition)
188//
189void MTFillMatrix::AddPreCut(const char *rule)
190{
191 MFilter *f = new MF(rule);
192 f->SetBit(kCanDelete);
193 AddPreCut(f);
194}
195
196//------------------------------------------------------------------------
197//
198// Add a cut which is used to fill the matrix. If kCanDelete is set
199// MJOptimize takes the ownership.
200//
201void MTFillMatrix::AddPreCut(MFilter *f)
202{
203 fPreCuts.Add(f);
204}
205
206//------------------------------------------------------------------------
207//
208// Add all entries deriving from MFilter from list to PreCuts.
209// The ownership is not affected.
210//
211void MTFillMatrix::AddPreCuts(const TList &list)
212{
213 TIter Next(&list);
214 TObject *obj=0;
215 while ((obj=Next()))
216 if (obj->InheritsFrom(MFilter::Class()))
217 fPreCuts.Add(obj);
218}
219
220// --------------------------------------------------------------------------
221//
222// Fill the matrix (FIXME: Flow diagram missing)
223//
224Bool_t MTFillMatrix::Process(const MParList &parlist)
225{
226 if (!fReader)
227 {
228 *fLog << err << "ERROR - No task to read data was given... abort." << endl;
229 return kFALSE;
230 }
231
232 *fLog << inf;
233 fLog->Separator(GetDescriptor());
234 *fLog << "Fill " << fDestMatrix1->GetDescriptor() << " with " << fNumDestEvents1 << endl;
235 if (fDestMatrix2)
236 *fLog << "Fill " << fDestMatrix2->GetDescriptor() << " with " << fNumDestEvents2 << endl;
237 *fLog << "Distribution choosen ";
238 if (fReference && fReference->GetHist().GetEntries()>0)
239 *fLog << "from " << fReference->GetDescriptor();
240 else
241 *fLog << "randomly";
242 *fLog << endl;
243
244 //
245 // Create parameter list and task list, add tasklist to parlist
246 //
247 parlist.Print();
248 MParList plist(parlist);
249 MTaskList tlist;
250 plist.AddToList(&tlist);
251
252 //
253 // A selector to select a given number of events from a sample
254 //
255 // FIXME: Merge MFEventSelector and MFEventSelector2
256 MFilter *selector=0;
257 if (fNumDestEvents1>0 || fNumDestEvents2>0)
258 {
259 if (fReference)
260 {
261 // Case of a reference/nominal distribution
262 // The events must be read before selection
263 MFEventSelector2 *sel = new MFEventSelector2(*fReference);
264 sel->SetNumMax(fNumDestEvents1+fNumDestEvents2);
265 sel->SetInverted();
266
267 selector = sel;
268 }
269 else
270 {
271 // Case of a random distribution
272 // The events can be selected before reading
273 MFEventSelector *sel = new MFEventSelector;
274 sel->SetNumSelectEvts(fNumDestEvents1+fNumDestEvents2);
275 fReader->SetSelector(sel);
276
277 selector = sel;
278 }
279 }
280
281 //
282 // Continue for PreCuts
283 //
284 MFilterList list;
285 list.SetName("PreCuts");
286 if (!list.AddToList(fPreCuts))
287 *fLog << err << "ERROR - Calling MFilterList::AddToList for fPreCuts failed!" << endl;
288
289 MContinue cont0(&list);
290 cont0.SetInverted();
291
292 //
293 // Continue for all events which are not (SetInverted())
294 // selected by the 'selector'
295 //
296 MContinue cont(selector);
297
298 //
299 // Create a filter doing a random split
300 //
301 const Double_t prob = (Double_t)fNumDestEvents1/(fNumDestEvents1+fNumDestEvents2);
302 MFEventSelector split;
303 split.SetSelectionRatio(prob);
304
305 //
306 // Create the logical inverted filter for 'split'
307 //
308 MFilterList invsplit;
309 invsplit.AddToList(&split);
310 invsplit.SetInverted();
311
312 //
313 // The two tasks filling the two matrices
314 //
315 MFillH fill1(fDestMatrix1);
316 MFillH fill2(fDestMatrix2);
317 if (selector)
318 {
319 fill1.SetFilter(&split);
320 fill2.SetFilter(&invsplit);
321 }
322
323 // entries in MTaskList
324 tlist.AddToList(fReader); // Read events
325 if (fPreCuts.GetEntries()>0)
326 tlist.AddToList(&cont0); // PreCuts
327 if (fReference && selector)
328 tlist.AddToList(&cont); // select a sample of events
329 tlist.AddToList(&invsplit); // process invsplit (which implicitly processes split)
330 if (fDestMatrix1)
331 tlist.AddToList(&fill1); // fill matrix 1
332 if (fDestMatrix2)
333 tlist.AddToList(&fill2); // fill matrix 2
334 if (fWriteFile1)
335 {
336 fWriteFile1->SetFilter(&split);
337 tlist.AddToList(fWriteFile1);
338 }
339 if (fWriteFile2)
340 {
341 fWriteFile2->SetFilter(&invsplit);
342 tlist.AddToList(fWriteFile2);
343 }
344
345 //
346 // Execute the eventloop
347 //
348 MEvtLoop evtloop(fName);
349 evtloop.SetParList(&plist);
350 evtloop.SetDisplay(fDisplay);
351 evtloop.SetLogStream(fLog);
352
353 const Bool_t rc = evtloop.Eventloop();
354
355 if (selector)
356 delete selector;
357
358 if (!rc)
359 {
360 *fLog << err << GetDescriptor() << ": Failed." << endl;
361 return kFALSE;
362 }
363
364 // Check the result of filling...
365 CheckResult(fDestMatrix1, fNumDestEvents1);
366 CheckResult(fDestMatrix2, fNumDestEvents2);
367
368 *fLog << inf << GetDescriptor() << ": Done." << endl;
369 return kTRUE;
370}
371
372// --------------------------------------------------------------------------
373//
374// Write both matrices to a file. Return kFALSE if writing one of them
375// failed or both have the same name.
376//
377Bool_t MTFillMatrix::WriteMatrices(const TString &fname) const
378{
379 if (fDestMatrix1 && fDestMatrix2 &&
380 (TString)fDestMatrix1->GetName()==(TString)fDestMatrix2->GetName())
381 {
382 *fLog << "ERROR - Cannot write matrices (both have the same name)... ignored." << endl;
383 return kFALSE;
384 }
385
386 return WriteMatrix1(fname) && WriteMatrix2(fname);
387}
388
389Int_t MTFillMatrix::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
390{
391 Bool_t rc = kFALSE;
392 if (IsEnvDefined(env, prefix, "NumDestEvents1", print))
393 {
394 rc = kTRUE;
395 SetNumDestEvents1(GetEnvValue(env, prefix, "NumDestEvents1", fNumDestEvents1));
396 }
397 if (IsEnvDefined(env, prefix, "NumDestEvents2", print))
398 {
399 rc = kTRUE;
400 SetNumDestEvents2(GetEnvValue(env, prefix, "NumDestEvents2", fNumDestEvents2));
401 }
402 return rc;
403}
Note: See TracBrowser for help on using the repository browser.