1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 #include "AliMUONQADataMakerSim.h"
19 #include "AliMUONHit.h"
20 #include "AliMUONDigit.h"
21 #include "AliMUONVHitStore.h"
22 #include "AliMUONVDigitStore.h"
24 // --- AliRoot header files ---
26 #include "AliQAChecker.h"
28 // --- ROOT system ---
29 #include <TClonesArray.h>
36 //-----------------------------------------------------------------------------
37 /// \class AliMUONQADataMakerSim
39 /// MUON base class for quality assurance data (histo) maker
44 ClassImp(AliMUONQADataMakerSim)
47 //____________________________________________________________________________
48 AliMUONQADataMakerSim::AliMUONQADataMakerSim() :
49 AliQADataMakerSim(AliQA::GetDetName(AliQA::kMUON), "MUON Quality Assurance Data Maker"),
56 //____________________________________________________________________________
57 AliMUONQADataMakerSim::AliMUONQADataMakerSim(const AliMUONQADataMakerSim& qadm) :
65 fHitStore = static_cast<AliMUONVHitStore*>(qadm.fHitStore->Clone());
67 if ( qadm.fDigitStore )
69 fDigitStore = static_cast<AliMUONVDigitStore*>(qadm.fDigitStore->Clone());
71 SetName((const char*)qadm.GetName()) ;
72 SetTitle((const char*)qadm.GetTitle());
75 //__________________________________________________________________
76 AliMUONQADataMakerSim& AliMUONQADataMakerSim::operator = (const AliMUONQADataMakerSim& qadm )
79 this->~AliMUONQADataMakerSim();
80 new(this) AliMUONQADataMakerSim(qadm);
84 //__________________________________________________________________
85 AliMUONQADataMakerSim::~AliMUONQADataMakerSim()
92 //__________________________________________________________________
93 void AliMUONQADataMakerSim::InitHits()
95 /// Initialized hit spectra
96 TH1F* h0 = new TH1F("hHitDetElem", "DetElemId distribution in Hits", 1400, 100., 1500.);
99 TH1F* h1 = new TH1F("hHitPtot", "P distribution in Hits ", 300, 0., 300.);
104 //__________________________________________________________________
105 void AliMUONQADataMakerSim::InitSDigits()
107 /// Initialized SDigits spectra
108 TH1I* h0 = new TH1I("hSDigitsDetElem", "Detection element distribution in SDigits", 1400, 100, 1500);
109 Add2SDigitsList(h0, 0);
111 TH1F* h1 = new TH1F("hSDigitsCharge", "Charge distribution in SDigits", 4096, 0, 4095);
112 Add2SDigitsList(h1, 1);
116 //__________________________________________________________________
117 void AliMUONQADataMakerSim::InitDigits()
119 /// Initialized Digits spectra
120 TH1I* h0 = new TH1I("hDigitsDetElem", "Detection element distribution in Digits", 1400, 100, 1500);
121 Add2DigitsList(h0, 0);
123 TH1I* h1 = new TH1I("hDigitsADC", "ADC distribution in Digits", 4096, 0, 4095);
124 Add2DigitsList(h1, 1);
128 //__________________________________________________________________
129 void AliMUONQADataMakerSim::MakeHits(TTree* hitsTree)
131 /// makes data from Hits
133 fHitStore = AliMUONVHitStore::Create(*hitsTree);
134 fHitStore->Connect(*hitsTree, false);
135 hitsTree->GetEvent(0);
137 TIter next(fHitStore->CreateIterator());
139 AliMUONHit* hit = 0x0;
141 while ( ( hit = static_cast<AliMUONHit*>(next()) ) )
143 GetHitsData(0)->Fill(hit->DetElemId());
144 GetHitsData(1)->Fill(hit->Momentum());
150 //__________________________________________________________________
151 void AliMUONQADataMakerSim::MakeSDigits(TTree* sdigitsTree)
153 /// makes data from SDigits
155 fDigitStore = AliMUONVDigitStore::Create(*sdigitsTree);
156 fDigitStore->Connect(*sdigitsTree, false);
157 sdigitsTree->GetEvent(0);
159 TIter next(fDigitStore->CreateIterator());
161 AliMUONVDigit* dig = 0x0;
163 while ( ( dig = static_cast<AliMUONVDigit*>(next()) ) )
165 GetSDigitsData(0)->Fill(dig->DetElemId());
166 GetSDigitsData(1)->Fill(dig->Charge());
170 //__________________________________________________________________
171 void AliMUONQADataMakerSim::MakeDigits(TTree* digitsTree)
173 /// makes data from Digits
175 fDigitStore = AliMUONVDigitStore::Create(*digitsTree);
176 fDigitStore->Connect(*digitsTree, false);
177 digitsTree->GetEvent(0);
179 TIter next(fDigitStore->CreateIterator());
181 AliMUONVDigit* dig = 0x0;
183 while ( ( dig = static_cast<AliMUONVDigit*>(next()) ) )
185 GetDigitsData(0)->Fill(dig->DetElemId());
186 GetDigitsData(1)->Fill(dig->ADC());
190 //____________________________________________________________________________
191 void AliMUONQADataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray* list)
193 ///Detector specific actions at end of cycle
194 // do the QA checking
195 AliQAChecker::Instance()->Run(AliQA::kMUON, task, list) ;
199 //____________________________________________________________________________
200 void AliMUONQADataMakerSim::StartOfDetectorCycle()
202 /// Detector specific actions at start of cycle