source: trunk/MagicSoft/Mars/mhbase/MHArray.cc@ 9374

Last change on this file since 9374 was 9153, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 20.2 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 07/2001 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26//
27// MHArray
28//
29// is a sequential collection of mars histograms. If the given index to
30// call the Fill function of the histogram excceeds the size of the
31// array by 1 a new entry is created.
32//
33// With Set/Inc/DecIndex you may specify the actual index of the histogram
34// wich should be filles by Fill.
35//
36// Use GetH to get the current histogram, the []-operator get the histogram
37// by its index.
38//
39// To access the histograms by a key instead of an index use SetIndexByKey
40// instead of Set/Inc/DecIndex. It will take the integerpart of the
41// floating point value (2 in case of 2.9). For each new key a new
42// index in the Mapping Table is created. So that you can access your
43// histograms by the key (eg in case of the Angle Theta=23.2deg use
44// SetIndexByKey(23.2)
45//
46// If the index is equal to the number of histograms in the array a call
47// to the Fill-member-function will create a new histogram.
48//
49// In the constructor istempl leads to two different behaviours of the
50// MHArray:
51//
52// - istempl=kTRUE tells MHArray to use the first histogram retrieved from
53// the Parameterlist by hname to be used as a template. New histograms
54// are not created using the root dictionary, but the New-member function
55// (see MParConatiner)
56// - In the case istempl=kFALSE new histograms are created using the root
57// dictionary while hname is the class name. For the creation their
58// default constructor is used.
59//
60//////////////////////////////////////////////////////////////////////////////
61#include "MHArray.h"
62
63#include <TH1.h>
64#include <TH2.h>
65#include <TH3.h>
66#include <TClass.h>
67#include <TStyle.h>
68#include <TGaxis.h>
69#include <TCanvas.h>
70#include <TLegend.h>
71#include <TPaveStats.h>
72
73#include "MLog.h"
74#include "MLogManip.h"
75
76#include "MParList.h"
77#include "MParContainer.h"
78
79#include "MBinning.h"
80
81ClassImp(MHArray);
82
83using namespace std;
84
85//////////////////////////////////////////////////////////////////////////////
86//
87// MMap
88//
89// This class maps a key-value to a given value. In its simple versions it
90// maps a key to an index.
91//
92//////////////////////////////////////////////////////////////////////////////
93#include <TArrayI.h>
94
95class MMap
96{
97private:
98 TArrayI fKeys;
99 TArrayI fValues;
100
101 Int_t K(Int_t i) const { return ((TArrayI)fKeys)[i]; }
102 Int_t V(Int_t i) const { return ((TArrayI)fValues)[i]; }
103
104public:
105 // --------------------------------------------------------------------------
106 //
107 // Return the size of the table
108 //
109 Int_t GetSize() const
110 {
111 return fKeys.GetSize();
112 }
113
114 // --------------------------------------------------------------------------
115 //
116 // Get the value which corresponds to the given key-value
117 //
118 Int_t GetValue(Int_t key) const
119 {
120 const Int_t n = fKeys.GetSize();
121 for (int i=0; i<n; i++)
122 {
123 if (K(i)==key)
124 return V(i);
125 }
126 return -1;
127 }
128
129 // --------------------------------------------------------------------------
130 //
131 // Get the key which corresponds to the given index
132 //
133 Int_t GetKey(Int_t value) const
134 {
135 const Int_t n = fKeys.GetSize();
136 for (int i=0; i<n; i++)
137 {
138 if (V(i)==value)
139 return K(i);
140 }
141 return -1;
142 }
143
144 // --------------------------------------------------------------------------
145 //
146 // Adds a new pair key-value. While the key is the key to the value.
147 // if the key already exists the pair is ignored.
148 //
149 void Add(Int_t key, Int_t value)
150 {
151 if (GetValue(key)>=0)
152 return;
153
154 const Int_t n = fKeys.GetSize();
155
156 fKeys.Set(n+1);
157 fValues.Set(n+1);
158
159 fKeys[n] = key;
160 fValues[n] = value;
161 }
162
163 // --------------------------------------------------------------------------
164 //
165 // Adds a new pair key-value. While the key is the key to the value.
166 // In this case the value is an automatically sequential created index.
167 // if the key already exists the pair is ignored.
168 //
169 Int_t Add(Int_t key)
170 {
171 const Int_t k = GetValue(key);
172 if (k>=0)
173 return k;
174
175 const Int_t n = fKeys.GetSize();
176
177 fKeys.Set(n+1);
178 fValues.Set(n+1);
179
180 fKeys[n] = key;
181 fValues[n] = n;
182
183 return n;
184 }
185};
186
187void MHArray::Init(const char *name)
188{
189 fName = name ? name : "MHArray";
190
191 fMapIdx = new MMap;
192
193 fArray = new TList;
194 fArray->SetOwner();
195}
196
197// --------------------------------------------------------------------------
198//
199// Can replace a constructor. Use the default constructor and afterwards
200// the Set function of your need.
201//
202void MHArray::Set(const TString hname, Bool_t istempl)
203{
204 if (fTemplate || fClass || fTemplateName!="<dummy>")
205 {
206 *fLog << warn << "WARNING - MHArray already setup... Set ignored." << endl;
207 return;
208 }
209
210 if (istempl)
211 {
212 fTemplateName = hname;
213 return;
214 }
215
216 //
217 // try to get class from root environment
218 //
219 fClass = gROOT->GetClass(hname);
220 if (!fClass)
221 {
222 //
223 // if class is not existing in the root environment
224 //
225 *fLog << err << dbginf << "Class '" << hname << "' not existing in dictionary." << endl;
226 }
227
228 //
229 // check for ineritance from MH
230 //
231 if (!fClass->InheritsFrom(MH::Class()))
232 {
233 //
234 // if class doesn't inherit from MH --> error
235 //
236 *fLog << err << dbginf << "Class '" << hname << "' doesn't inherit from MH." << endl;
237 fClass = NULL;
238 }
239}
240
241// --------------------------------------------------------------------------
242//
243// Can replace a constructor. Use the default constructor and afterwards
244// the Set function of your need.
245//
246void MHArray::Set(const MH *hist)
247{
248 fIdx=0;
249 fClass=NULL;
250 fTemplate=hist;
251 fTemplateName="<dummy>";
252}
253
254
255// --------------------------------------------------------------------------
256//
257// hname is the name of the histogram class which is in the array.
258//
259// istempl=kTRUE tells MHArray to use the first histogram retrieved from the
260// ParameterList by hname to be used as a template. New histograms are not
261// created using the root dictionary, but the New-member function (see
262// MParConatiner)
263// In the case istempl=kFALSE new histograms are created using the root
264// dictionary while hname is the class name. For the creation their
265// default constructor is used.
266//
267MHArray::MHArray(const TString hname, Bool_t istempl, const char *name, const char *title)
268 : fIdx(0), fClass(NULL), fTemplate(NULL)
269{
270 //
271 // set the name and title of this object
272 //
273 Init(name);
274 fTitle = title ? TString(title) : (TString("Base class for Mars histogram arrays:") + hname);
275
276 Set(hname, istempl);
277}
278
279// --------------------------------------------------------------------------
280//
281// Default constructor. Use MHArray::Set to setup the MHArray afterwards
282//
283MHArray::MHArray(const char *name, const char *title)
284 : fIdx(0), fClass(NULL), fTemplate(NULL), fTemplateName("<dummy>")
285{
286 //
287 // set the name and title of this object
288 //
289 Init(name);
290 fTitle = title ? title : "A Mars histogram array";
291}
292
293// --------------------------------------------------------------------------
294//
295// hname is the name of the histogram class which is in the array.
296//
297// istempl=kTRUE tells MHArray to use the first histogram retrieved from the
298// ParameterList by hname to be used as a template. New histograms are not
299// created using the root dictionary, but the New-member function (see
300// MParConatiner)
301// In the case istempl=kFALSE new histograms are created using the root
302// dictionary while hname is the class name. For the creation their
303// default constructor is used.
304//
305MHArray::MHArray(const MH *hist, const char *name, const char *title)
306 : fIdx(0), fClass(NULL), fTemplate(hist), fTemplateName("<dummy>")
307{
308 //
309 // set the name and title of this object
310 //
311 Init(name);
312 fTitle = title ? TString(title) : (TString("Base class for Mars histogram arrays:") + hist->GetName());
313}
314
315// --------------------------------------------------------------------------
316//
317// Destructor: Deleteing the array and all histograms which are part of the
318// array.
319//
320MHArray::~MHArray()
321{
322 fArray->Delete();
323 delete fArray;
324 delete fMapIdx;
325}
326
327// --------------------------------------------------------------------------
328//
329// Use this to access the histograms by a key. If you use values like
330// (in this order) 2.5, 7.2, 2.5, 9.3, 9.3, 3.3, 2.2, 1.1
331// it will be mapped to the following indices internally:
332// 0 1 0 2 2 3 4 5
333//
334// WARNING: Make sure that you don't create new histograms by setting
335// a new index (SetIndex or IncIndex) which is equal the size
336// of the array and create new histogram by CreateH. In this
337// case you will confuse the mapping completely.
338//
339void MHArray::SetIndexByKey(Double_t key)
340{
341 fIdx = fMapIdx->Add((Int_t)key);
342}
343
344// --------------------------------------------------------------------------
345//
346// Use this function to access a histogram by its index in the array.
347// Becarefull the index isn't checked!
348//
349MH &MHArray::operator[](Int_t i)
350{
351 return *(MH*)fArray->At(i);
352}
353
354// --------------------------------------------------------------------------
355//
356// Use this function to access a histogram by its index in the array.
357// Becarefull the index isn't checked!
358//
359MH *MHArray::At(Int_t i)
360{
361 return (MH*)fArray->At(i);
362}
363
364// --------------------------------------------------------------------------
365//
366// Use this function to access the histogram corresponding to the
367// currently set index (by Set/Inc/DecIndex or SetIndexByKey)
368// Becarefull the index set isn't checked!
369//
370MH *MHArray::GetH()
371{
372 return (MH*)fArray->At(fIdx);
373}
374
375// --------------------------------------------------------------------------
376//
377// Tries to create a new histogram, adds it as last entry to the array
378// and tries to call SetupFill for it. In case of success the last entry
379// in the array is the new histogram and kTRUE is returned. Otherwise
380// kFALSE is returned.
381//
382Bool_t MHArray::CreateH()
383{
384 TString cname = fClass ? fClass->GetName() : fTemplate->IsA()->GetName();
385
386 MH *hist = NULL;
387 if (fTemplate)
388 {
389 //
390 // create the parameter container as a clone of the existing
391 // template histogram.
392 //
393 hist = (MH*)fTemplate->New();
394 }
395 else
396 {
397 //
398 // create the parameter container of the the given class type
399 //
400 hist = (MH*)fClass->New();
401 }
402 if (!hist)
403 {
404 *fLog << err << dbginf << "Cannot create new instance of class '";
405 *fLog << cname << "' (Maybe no def. constructor)" << endl;
406 return kFALSE;
407 }
408
409 //
410 // Set the name of the container
411 //
412 if (!fTemplate)
413 {
414 TString name = TString(hist->GetName())+";";
415 name += fIdx;
416
417 hist->SetName(name);
418 }
419
420 //
421 // Try to setup filling for the histogram
422 //
423 if (!hist->SetupFill(fParList))
424 {
425 *fLog << err << dbginf << "SetupFill for new histogram of type '";
426 *fLog << cname << "' with Index #" << fIdx << " failed." << endl;
427 delete hist;
428 return kFALSE;
429 }
430
431 fArray->AddLast(hist);
432
433 return kTRUE;
434}
435
436// --------------------------------------------------------------------------
437//
438// Returns kFALSE if the class couldn't be found in the root dictionary or
439// if it doesn't inherit from MH.
440// The parameter list is remembert to be used for SetupFill in case a new
441// histogram is created.
442// The index is reset to 0
443//
444Bool_t MHArray::SetupFill(const MParList *pList)
445{
446 fParList = pList;
447 fIdx = 0;
448
449 if (fTemplate)
450 return kTRUE;
451
452 if (!fTemplateName.IsNull())
453 {
454 fTemplate = (MH*)pList->FindObject(fTemplateName, "MH");
455 return fTemplate ? kTRUE : kFALSE;
456 }
457
458 return fClass ? kTRUE : kFALSE;
459}
460
461// --------------------------------------------------------------------------
462//
463// Call Fill for the present histogram index. If the index is out of
464// bounds the event is skipped. If the index is the number of current
465// histograms in the array a new histogram is created and if creation was
466// successfull filled.
467//
468Int_t MHArray::Fill(const MParContainer *par, const Stat_t w)
469{
470 const Int_t n = fArray->GetSize();
471
472 if (fIdx<0 || fIdx>n)
473 {
474 *fLog << warn << "Histogram Index #" << fIdx << " out of bounds (>";
475 *fLog << n << ")... skipped." << endl;
476 return kERROR;
477 }
478
479 if (fIdx==n)
480 if (!CreateH())
481 return kERROR;
482
483 return GetH()->Fill(par, w);
484}
485
486Bool_t MHArray::AddHistogram()
487{
488 fIdx=fArray->GetSize();
489
490 return CreateH();
491}
492
493// --------------------------------------------------------------------------
494//
495// Calls Finalize for all histograms in the list. If at least one Finalize
496// fails kFALSE is returned.
497//
498Bool_t MHArray::Finalize()
499{
500 Bool_t rc = kTRUE;
501
502 TIter Next(fArray);
503 MH *hist = NULL;
504
505 while ((hist=(MH*)Next()))
506 if (!hist->Finalize())
507 rc = kFALSE;
508
509 return rc;
510}
511
512// --------------------------------------------------------------------------
513//
514// Print the number of entries in the array
515//
516void MHArray::Print(Option_t *option) const
517{
518 *fLog << all << GetDescriptor() << " contains " << fArray->GetSize();
519 *fLog << " histograms." << endl;
520
521 if (fMapIdx->GetSize()<=0)
522 return;
523
524 *fLog << " idx\t key" << endl;
525 for (int i=0; i<fMapIdx->GetSize(); i++)
526 *fLog << " " << i << "\t<--> " << fMapIdx->GetKey(i) << endl;
527 *fLog << endl;
528}
529
530// --------------------------------------------------------------------------
531//
532// Adds the given object to the given legend (if != NULL). The Legend
533// entry name is created from the key...
534//
535void MHArray::AddLegendEntry(TLegend *leg, TObject *obj, Int_t idx) const
536{
537 if (!leg)
538 return;
539
540 TString name = " ";
541 name += fMapIdx->GetKey(idx);
542 leg->AddEntry(obj, name, "lpf"); // l=line, p=polymarker, f=fill
543}
544
545
546// --------------------------------------------------------------------------
547//
548// The option is the name of the histogram, used to get a histogram from
549// the MH entries by calling their GetHist function.
550//
551void MHArray::Draw(Option_t *opt)
552{
553 if (!gPad)
554 MH::MakeDefCanvas(this);
555
556 const Int_t sstyle = gStyle->GetOptStat();
557 gStyle->SetOptStat(0);
558
559 //
560 // if the keymapping is used create a legend to identify the histograms
561 //
562 TLegend *leg = NULL;
563 if (fMapIdx->GetSize()>0)
564 {
565 leg = new TLegend(0.85, 0.80, 0.99, 0.99);
566 leg->SetBit(kCanDelete);
567 }
568
569 TIter Next(fArray);
570 MH *hist = (MH*)Next();
571
572 Int_t idx=0;
573 Double_t max=0;
574 Double_t min=0;
575
576 TH1 *h1=NULL;
577
578 //
579 // If the array has at least one entry:
580 // - find the starting boundaries
581 // - draw it and set its line color
582 //
583 if (hist)
584 {
585 if ((h1 = hist->GetHistByName(opt)))
586 {
587 h1->Draw();
588 h1->SetLineColor(idx+2);
589 max = h1->GetMaximum();
590 min = h1->GetMinimum();
591
592 AddLegendEntry(leg, h1, idx);
593 }
594 }
595
596 //
597 // For all following histograms:
598 // - update the boundaries
599 // - draw it and set its line color
600 //
601 while ((hist=(MH*)Next()))
602 {
603 TH1 *h=NULL;
604
605 if (!(h = hist->GetHistByName(opt)))
606 continue;
607
608 h->Draw("same");
609 h->SetLineColor(idx+2);
610 if (max<h->GetMaximum())
611 max = h->GetMaximum();
612 if (min>h->GetMinimum())
613 min = h->GetMinimum();
614
615 AddLegendEntry(leg, h, idx++);
616 }
617
618 //
619 // Now update the drawing region so that everything is displayed
620 //
621 if (h1)
622 {
623 h1->SetMinimum(min>0 ? min*0.95 : min*1.05);
624 h1->SetMaximum(max>0 ? max*1.05 : max*0.95);
625 }
626
627 if (leg)
628 leg->Draw();
629
630 gPad->Modified();
631 gPad->Update();
632
633 gStyle->SetOptStat(sstyle);
634}
635
636// --------------------------------------------------------------------------
637//
638// The option is the name of the histogram, used to get a histogram from
639// the MH entries by calling their GetHistByName function.
640// If the option also contains 'nonew' no new canvas is created.
641// The option "Scale=1" scales the area of all histogram to 1
642// The option "Scale=max" scales the maximum of all histogram to 1
643//
644TObject *MHArray::DrawClone(Option_t *opt) const
645{
646 TString o(opt);
647
648 TCanvas *c = NULL;
649
650 Int_t scale1 = o.Index("scale=1", TString::kIgnoreCase);
651 Int_t scalemax = o.Index("scale=max", TString::kIgnoreCase);
652 Int_t nonew = o.Index("nonew", TString::kIgnoreCase);
653
654 if (o.BeginsWith("scale=1", TString::kIgnoreCase))
655 scale1 = 0;
656 if (o.BeginsWith("scale=max", TString::kIgnoreCase))
657 scalemax = 0;
658 if (o.BeginsWith("nonew", TString::kIgnoreCase))
659 nonew = 0;
660
661 if (nonew>=0)
662 {
663 c = MH::MakeDefCanvas(this);
664
665 //
666 // This is necessary to get the expected bahviour of DrawClone
667 //
668 gROOT->SetSelectedPad(NULL);
669
670 o.Remove(nonew, 5);
671 }
672 if (scale1>=0)
673 o.Remove(scale1, 7);
674 if (scalemax>=0)
675 o.Remove(scalemax, 9);
676
677 const Int_t sstyle = gStyle->GetOptStat();
678 gStyle->SetOptStat(0);
679
680 //
681 // if the keymapping is used create a legend to identify the histograms
682 //
683 TLegend *leg = NULL;
684 if (fMapIdx->GetSize()>0)
685 {
686 leg = new TLegend(0.85, 0.80, 0.99, 0.99);
687 leg->SetBit(kCanDelete);
688 }
689
690 TIter Next(fArray);
691 MH *hist = (MH*)Next();
692
693 Int_t idx=0;
694 Double_t max=0;
695 Double_t min=0;
696
697 TH1 *h1=NULL;
698
699 //
700 // If the array has at least one entry:
701 // - find the starting boundaries
702 // - draw it and set its line color
703 //
704 if (hist)
705 {
706 if ((h1 = hist->GetHistByName(o)))
707 {
708 h1 = (TH1*)h1->DrawCopy();
709
710 if (scale1>=0)
711 h1->Scale(1./h1->Integral());
712 if (scalemax>=0)
713 h1->Scale(1./h1->GetMaximum());
714
715 h1->SetMarkerColor(idx);
716 h1->SetLineColor(idx+2);
717 h1->SetFillStyle(4000);
718 max = h1->GetMaximum();
719 min = h1->GetMinimum();
720
721 AddLegendEntry(leg, h1, idx++);
722 }
723 }
724
725 //
726 // For all following histograms:
727 // - update the boundaries
728 // - draw it and set its line color
729 //
730 while ((hist=(MH*)Next()))
731 {
732 TH1 *h=NULL;
733
734 if (!(h = hist->GetHistByName(o)))
735 continue;
736
737 h = (TH1*)h->DrawCopy("same");
738
739 if (scale1>=0)
740 h->Scale(1./h->Integral());
741 if (scalemax>=0)
742 h->Scale(1./h->GetMaximum());
743
744 h->SetMarkerColor(idx);
745 h->SetLineColor(idx+2);
746 h->SetFillStyle(4000); // transperent (why is this necessary?)
747 if (max<h->GetMaximum())
748 max = h->GetMaximum();
749 if (min>h->GetMinimum())
750 min = h->GetMinimum();
751
752 AddLegendEntry(leg, h, idx++);
753 }
754
755 //
756 // Now update the drawing region so that everything is displayed
757 //
758 if (h1)
759 {
760 h1->SetMinimum(min>0 ? min*0.95 : min*1.05);
761 h1->SetMaximum(max>0 ? max*1.05 : max*0.95);
762 }
763
764 if (leg)
765 leg->Draw();
766
767 gPad->Modified();
768 gPad->Update();
769
770 gStyle->SetOptStat(sstyle);
771
772 return c;
773}
Note: See TracBrowser for help on using the repository browser.