/* ======================================================================== *\ ! ! * ! * This file is part of MARS, the MAGIC Analysis and Reconstruction ! * Software. It is distributed to you in the hope that it can be a useful ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. ! * It is distributed WITHOUT ANY WARRANTY. ! * ! * Permission to use, copy, modify and distribute this software and its ! * documentation for any purpose is hereby granted without fee, ! * provided that the above copyright notice appear in all copies and ! * that both that copyright notice and this permission notice appear ! * in supporting documentation. It is provided "as is" without express ! * or implied warranty. ! * ! ! ! Author(s): Thomas Bretz 07/2001 ! ! Copyright: MAGIC Software Development, 2000-2002 ! ! \* ======================================================================== */ ////////////////////////////////////////////////////////////////////////////// // // MHArray // // is a sequential collection of mars histograms. If the given index to // call the Fill function of the histogram excceeds the size of the // array by 1 a new entry is created. // // With Set/Inc/DecIndex you may specify the actual index of the histogram // wich should be filles by Fill. // // Use GetH to get the current histogram, the []-operator get the histogram // by its index. // // To map a key to the index use/see MFillH, which automatically maps a // key-value to the array index. // // In the constructor istempl leads to two different behaviours of the // MHArray: // // - istempl=kTRUE tells MHArray to use the first histogram retrieved from // the Parameterlist by hname to be used as a template. New histograms // are not created using the root dictionary, but the New-member function // (see MParConatiner) // - In the case istempl=kFALSE new histograms are created using the root // dictionary while hname is the class name. For the creation their // default constructor is used. // ////////////////////////////////////////////////////////////////////////////// #include "MHArray.h" #include #include #include #include #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MParList.h" #include "MParContainer.h" #include "MBinning.h" ClassImp(MHArray); // -------------------------------------------------------------------------- // // Default Constructor. hname is the name of the histogram class which // is in the array. // // istempl=kTRUE tells MHArray to use the first histogram retrieved from the // ParameterList by hname to be used as a template. New histograms are not // created using the root dictionary, but the New-member function (see // MParConatiner) // In the case istempl=kFALSE new histograms are created using the root // dictionary while hname is the class name. For the creation their // default constructor is used. // MHArray::MHArray(const TString hname, Bool_t istempl, const char *name, const char *title) : fIdx(0), fClass(NULL), fTemplate(NULL) { // // set the name and title of this object // fName = name ? name : "MHArray"; fTitle = title ? TString(title) : (TString("Base class for Mars histogram arrays:") + hname); fArray = new TList; fArray->SetOwner(); if (istempl) { fTemplateName = hname; return; } // // try to get class from root environment // fClass = gROOT->GetClass(hname); if (!fClass) { // // if class is not existing in the root environment // *fLog << err << dbginf << "Class '" << hname << "' not existing in dictionary." << endl; } // // check for ineritance from MH // if (!fClass->InheritsFrom(MH::Class())) { // // if class doesn't inherit from MH --> error // *fLog << err << dbginf << "Class '" << hname << "' doesn't inherit from MH." << endl; fClass = NULL; } } // -------------------------------------------------------------------------- // // Destructor: Deleteing the array and all histograms which are part of the // array. // MHArray::~MHArray() { fArray->Delete(); delete fArray; } MH &MHArray::operator[](Int_t i) { return *(MH*)fArray->At(i); } MH *MHArray::At(Int_t i) { return (MH*)fArray->At(i); } MH *MHArray::GetH() { return (MH*)fArray->At(fIdx); } // -------------------------------------------------------------------------- // // Tries to create a new histogram, adds it as last entry to the array // and tries to call SetupFill for it. In case of success the last entry // in the array is the new histogram and kTRUE is returned. Otherwise // kFALSE is returned. // Bool_t MHArray::CreateH() { TString cname = fClass ? fClass->GetName() : fTemplate->IsA()->GetName(); MH *hist = NULL; if (fTemplate) { hist = (MH*)fTemplate->New(); } else { // // create the parameter container of the the given class type // hist = (MH*)fClass->New(); } if (!hist) { *fLog << err << dbginf << "Cannot create new instance of class '"; *fLog << cname << "' (Maybe no def. constructor)" << endl; return kFALSE; } // // Set the name of the container // if (!fTemplate) { TString name = TString(hist->GetName())+";"; name += fIdx; hist->SetName(name); } // // Try to setup filling for the histogram // if (!hist->SetupFill(fParList)) { *fLog << err << dbginf << "SetupFill for new histogram of type '"; *fLog << cname << "' with Index #" << fIdx << " failed." << endl; delete hist; return kFALSE; } fArray->AddLast(hist); return kTRUE; } // -------------------------------------------------------------------------- // // Returns kFALSE if the class couldn't be found in the root dictionary or // if it doesn't inherit from MH. // The parameter list is remembert to be used for SetupFill in case a new // histogram is created. // The index is reset to 0 // Bool_t MHArray::SetupFill(const MParList *pList) { fParList = pList; fIdx = 0; if (!fTemplateName.IsNull()) { fTemplate = (MH*)pList->FindObject(fTemplateName, "MH"); return fTemplate ? kTRUE : kFALSE; } return fClass ? kTRUE : kFALSE; } // -------------------------------------------------------------------------- // // Call Fill for the present histogram index. If the index is out of // bounds the event is skipped. If the index is the number of current // histograms in the array a new histogram is created and if creation was // successfull filled. // Bool_t MHArray::Fill(const MParContainer *par) { const Int_t n = fArray->GetSize(); if (fIdx <0 || fIdx>n) { *fLog << warn << "Histogram Index #" << fIdx << " out of bounds (>"; *fLog << n << ")... skipped." << endl; return kCONTINUE; } if (fIdx==n) if (!CreateH()) return kFALSE; return GetH()->Fill(par); } // -------------------------------------------------------------------------- // // Calls Finalize for all histograms in the list. If at least one Finalize // fails kFALSE is returned. // Bool_t MHArray::Finalize() { Bool_t rc = kTRUE; TIter Next(fArray); MH *hist = NULL; while ((hist=(MH*)Next())) if (!hist->Finalize()) rc = kFALSE; return rc; } // -------------------------------------------------------------------------- // // Print the number of entries in the array // void MHArray::Print(Option_t *option) const { *fLog << all << GetDescriptor() << " contains " << fArray->GetSize(); *fLog << " histograms." << endl; } // -------------------------------------------------------------------------- // // The option is the name of the histogram, used to get a histogram from // the MH entries by calling their GetHist function. // void MHArray::Draw(Option_t *opt) { if (!gPad) MH::MakeDefCanvas(this); const Stat_t sstyle = gStyle->GetOptStat(); gStyle->SetOptStat(0); TIter Next(fArray); MH *hist = (MH*)Next(); Int_t col=2; Double_t max=0; Double_t min=0; TH1 *h1=NULL; if (hist) { if ((h1 = hist->GetHistByName(opt))) { h1->Draw(); h1->SetLineColor(col++); max = h1->GetMaximum(); min = h1->GetMinimum(); } } while ((hist=(MH*)Next())) { TH1 *h=NULL; if (!(h = hist->GetHistByName(opt))) continue; h->Draw("same"); h->SetLineColor(col++); if (maxGetMaximum()) max = h->GetMaximum(); if (min>h->GetMinimum()) min = h->GetMinimum(); } if (h1) { h1->SetMinimum(min>0 ? min*0.95 : min*1.05); h1->SetMaximum(max>0 ? max*1.05 : max*0.95); } gPad->Modified(); gPad->Update(); gStyle->SetOptStat(sstyle); } // -------------------------------------------------------------------------- // // The option is the name of the histogram, used to get a histogram from // the MH entries by calling their GetHistByName function. // If the option also contains 'nonew' no new canvas is created. // TObject *MHArray::DrawClone(Option_t *opt) const { TString o(opt); TCanvas *c = NULL; Int_t nonew = o.Index("nonew", TString::kIgnoreCase); if (nonew>=0) { c = MH::MakeDefCanvas(this); // // This is necessary to get the expected bahviour of DrawClone // gROOT->SetSelectedPad(NULL); o.Remove(nonew, 5); } const Stat_t sstyle = gStyle->GetOptStat(); gStyle->SetOptStat(0); TIter Next(fArray); MH *hist = (MH*)Next(); Int_t col=2; Double_t max=0; Double_t min=0; TH1 *h1=NULL; if (hist) { if ((h1 = hist->GetHistByName(o))) { h1 = (TH1*)h1->DrawCopy(); h1->SetMarkerColor(col); h1->SetLineColor(col++); h1->SetFillStyle(4000); max = h1->GetMaximum(); min = h1->GetMinimum(); } } while ((hist=(MH*)Next())) { TH1 *h=NULL; if (!(h = hist->GetHistByName(o))) continue; h = (TH1*)h->DrawCopy("same"); h->SetMarkerColor(col); h->SetLineColor(col++); h->SetFillStyle(4000); if (maxGetMaximum()) max = h->GetMaximum(); if (min>h->GetMinimum()) min = h->GetMinimum(); } if (h1) { h1->SetMinimum(min>0 ? min*0.95 : min*1.05); h1->SetMaximum(max>0 ? max*1.05 : max*0.95); } gPad->Modified(); gPad->Update(); gStyle->SetOptStat(sstyle); return c; }