source: trunk/MagicSoft/Mars/mtools/MTFillMatrix.cc@ 8785

Last change on this file since 8785 was 8644, checked in by tbretz, 17 years ago
*** empty log message ***
File size: 13.5 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 "MFDataPhrase.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
156void MTFillMatrix::Init(const char *name, const char *title)
157{
158 fName = name ? name : "MFillMatrix";
159 fTitle = title ? title : "Tool to fill MHMatrix from file";
160}
161
162MTFillMatrix::MTFillMatrix(const char *name, const char *title)
163: fReference(0), fReader(0), fDestMatrix1(0), fDestMatrix2(0),
164 fNumDestEvents1(0), fNumDestEvents2(0), fWriteFile1(0), fWriteFile2(0)
165{
166 Init(name, title);
167}
168
169// --------------------------------------------------------------------------
170//
171// Constructor. Takes an MH3 as argument. This MH3 is the reference
172// distribution to fill the matrix. More information can be found
173// at MFEventSelector2 which is used to select the events.
174//
175// If no MH3 *ref is given the events are randomly selected from the
176// total sample - this may result in samples which don't have exactly
177// the predefined size, but it is much faster.
178//
179MTFillMatrix::MTFillMatrix(const MH3 *ref, const char *name, const char *title)
180: fReference(0), fReader(0), fDestMatrix1(0), fDestMatrix2(0),
181 fNumDestEvents1(0), fNumDestEvents2(0), fWriteFile1(0), fWriteFile2(0)
182{
183 Init(name, title);
184 if (ref)
185 fReference = (MH3*)ref->Clone();
186}
187
188MTFillMatrix::~MTFillMatrix()
189{
190 if (fReference)
191 delete fReference;
192}
193
194//------------------------------------------------------------------------
195//
196// Function serving AddPreCuts, AddPreTasks and AddPostTasks
197//
198void MTFillMatrix::Add(const TList &src, const TClass *cls, TList &dest)
199{
200 TIter Next(&src);
201 TObject *obj=0;
202 while ((obj=Next()))
203 if (obj->InheritsFrom(cls))
204 dest.Add(obj);
205}
206
207//------------------------------------------------------------------------
208//
209// Add a cut which is used to fill the matrix, eg "MMcEvt.fOartId<1.5"
210// (The rule is applied, nit inverted: The matrix is filled with
211// the events fullfilling the condition)
212//
213void MTFillMatrix::AddPreCut(const char *rule)
214{
215 MFilter *f = new MFDataPhrase(rule);
216 f->SetBit(kCanDelete);
217 AddPreCut(f);
218}
219
220//------------------------------------------------------------------------
221//
222// Add a cut which is used to fill the matrix. If kCanDelete is set
223// MJOptimize takes the ownership.
224//
225void MTFillMatrix::AddPreCut(MFilter *f)
226{
227 fPreCuts.Add(f);
228}
229
230//------------------------------------------------------------------------
231//
232// Add all entries deriving from MFilter from list to PreCuts.
233// The ownership is not affected.
234//
235void MTFillMatrix::AddPreCuts(const TList &list)
236{
237 Add(list, MFilter::Class(), fPreCuts);
238}
239
240//------------------------------------------------------------------------
241//
242// Add a task which is executed before the precuts. If kCanDelete is set
243// MJOptimize takes the ownership.
244//
245void MTFillMatrix::AddPreTask(MTask *t)
246{
247 fPreTasks.Add(t);
248}
249
250//------------------------------------------------------------------------
251//
252// Add all entries deriving from MTask from list to PreTasks.
253// The ownership is not affected.
254//
255void MTFillMatrix::AddPreTasks(const TList &list)
256{
257 Add(list, MTask::Class(), fPreTasks);
258}
259
260//------------------------------------------------------------------------
261//
262// Add a task which is executed after the precuts. If kCanDelete is set
263// MJOptimize takes the ownership.
264//
265void MTFillMatrix::AddPostTask(MTask *t)
266{
267 fPostTasks.Add(t);
268}
269
270//------------------------------------------------------------------------
271//
272// Add all entries deriving from MTask from list to PostTasks.
273// The ownership is not affected.
274//
275void MTFillMatrix::AddPostTasks(const TList &list)
276{
277 Add(list, MTask::Class(), fPostTasks);
278}
279
280// --------------------------------------------------------------------------
281//
282// Fill the matrix (FIXME: Flow diagram missing)
283//
284Bool_t MTFillMatrix::Process(const MParList &parlist)
285{
286 if (!fReader)
287 {
288 *fLog << err << "ERROR - No task to read data was given... abort." << endl;
289 return kFALSE;
290 }
291
292 *fLog << inf;
293 fLog->Separator(GetDescriptor());
294 *fLog << "Fill " << fDestMatrix1->GetDescriptor() << " with " << fNumDestEvents1 << endl;
295 if (fDestMatrix2)
296 *fLog << "Fill " << fDestMatrix2->GetDescriptor() << " with " << fNumDestEvents2 << endl;
297 *fLog << "Distribution choosen ";
298 if (fReference && fReference->GetHist().GetEntries()>0)
299 *fLog << "from " << fReference->GetDescriptor();
300 else
301 *fLog << "randomly";
302 *fLog << endl;
303
304 //
305 // Create parameter list and task list, add tasklist to parlist
306 //
307 parlist.Print();
308 MParList plist(parlist);
309 MTaskList tlist;
310 plist.AddToList(&tlist);
311
312 //
313 // A selector to select a given number of events from a sample
314 //
315 // FIXME: Merge MFEventSelector and MFEventSelector2
316 MFilter *selector=0;
317 if (fNumDestEvents1>0 || fNumDestEvents2>0)
318 {
319 if (fReference)
320 {
321 // Case of a reference/nominal distribution
322 // The events must be read before selection
323 MFEventSelector2 *sel = new MFEventSelector2(*fReference);
324 sel->SetNumMax(fNumDestEvents1+fNumDestEvents2);
325 sel->SetInverted();
326
327 selector = sel;
328 }
329 else
330 {
331 // Case of a random distribution
332 // The events can be selected before reading
333 MFEventSelector *sel = new MFEventSelector;
334 sel->SetNumSelectEvts(fNumDestEvents1+fNumDestEvents2);
335 fReader->SetSelector(sel);
336
337 selector = sel;
338 }
339 }
340
341 //
342 // Continue for PreCuts
343 //
344 MFilterList list;
345 list.SetName("PreCuts");
346 if (!list.AddToList(fPreCuts))
347 *fLog << err << "ERROR - Calling MFilterList::AddToList for fPreCuts failed!" << endl;
348
349 MContinue cont0(&list);
350 cont0.SetInverted();
351
352 //
353 // Continue for all events which are not (SetInverted())
354 // selected by the 'selector'
355 //
356 MContinue cont(selector);
357
358 //
359 // Create a filter doing a random split
360 //
361 const Double_t prob = (Double_t)fNumDestEvents1/(fNumDestEvents1+fNumDestEvents2);
362 MFEventSelector split;
363 split.SetSelectionRatio(prob);
364
365 //
366 // Create the logical inverted filter for 'split'
367 //
368 MFilterList invsplit;
369 invsplit.AddToList(&split);
370 invsplit.SetInverted();
371
372 //
373 // The two tasks filling the two matrices
374 //
375 MFillH fill1(fDestMatrix1);
376 MFillH fill2(fDestMatrix2);
377 fill1.SetFilter(&split);
378 fill2.SetFilter(&invsplit);
379
380 // entries in MTaskList
381 tlist.AddToList(fReader); // Read events
382 tlist.AddToList(fPreTasks, 0); // PreTasks
383 if (fPreCuts.GetEntries()>0)
384 tlist.AddToList(&cont0); // PreCuts
385 if (fReference && selector)
386 tlist.AddToList(&cont); // select a sample of events
387 tlist.AddToList(&invsplit); // process invsplit (which implicitly processes split)
388 tlist.AddToList(fPostTasks, 0); // PostTasks
389 if (fDestMatrix1)
390 tlist.AddToList(&fill1); // fill matrix 1
391 if (fDestMatrix2)
392 tlist.AddToList(&fill2); // fill matrix 2
393 if (fWriteFile1)
394 {
395 fWriteFile1->SetFilter(&split);
396 tlist.AddToList(fWriteFile1);
397 }
398 if (fWriteFile2)
399 {
400 fWriteFile2->SetFilter(&invsplit);
401 tlist.AddToList(fWriteFile2);
402 }
403
404 tlist.SetAccelerator(MTask::kAccDontReset|MTask::kAccDontTime);
405
406 //
407 // Execute the eventloop
408 //
409 MEvtLoop evtloop(fName);
410 evtloop.SetParList(&plist);
411 evtloop.SetDisplay(fDisplay);
412 evtloop.SetLogStream(fLog);
413
414 const Bool_t rc = evtloop.Eventloop();
415
416 if (selector)
417 delete selector;
418
419 if (!rc)
420 {
421 *fLog << err << GetDescriptor() << ": Failed." << endl;
422 return kFALSE;
423 }
424
425 // Check the result of filling...
426 CheckResult(fDestMatrix1, fNumDestEvents1);
427 CheckResult(fDestMatrix2, fNumDestEvents2);
428
429 *fLog << inf << GetDescriptor() << ": Done." << endl;
430 return kTRUE;
431}
432
433// --------------------------------------------------------------------------
434//
435// Write both matrices to a file. Return kFALSE if writing one of them
436// failed or both have the same name.
437//
438Bool_t MTFillMatrix::WriteMatrices(const TString &fname) const
439{
440 if (fDestMatrix1 && fDestMatrix2 &&
441 (TString)fDestMatrix1->GetName()==(TString)fDestMatrix2->GetName())
442 {
443 *fLog << "ERROR - Cannot write matrices (both have the same name)... ignored." << endl;
444 return kFALSE;
445 }
446
447 return WriteMatrix1(fname) && WriteMatrix2(fname);
448}
449
450Int_t MTFillMatrix::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
451{
452 Bool_t rc = kFALSE;
453 if (IsEnvDefined(env, prefix, "NumDestEvents1", print))
454 {
455 rc = kTRUE;
456 SetNumDestEvents1(GetEnvValue(env, prefix, "NumDestEvents1", fNumDestEvents1));
457 }
458 if (IsEnvDefined(env, prefix, "NumDestEvents2", print))
459 {
460 rc = kTRUE;
461 SetNumDestEvents2(GetEnvValue(env, prefix, "NumDestEvents2", fNumDestEvents2));
462 }
463 return rc;
464}
Note: See TracBrowser for help on using the repository browser.