]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/TestRecPoints.C
Fixing the decoding of regional header bits.
[u/mrichter/AliRoot.git] / MUON / TestRecPoints.C
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
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 **************************************************************************/
15
16 /* $Id: TestRecPoints.C 23207 2007-12-20 09:59:20Z ivana $ */
17
18 /// \ingroup macros
19 /// \file TestRecPoints.C
20 /// \brief This macro is to be used to check the trigger and tracking chamber behaviour
21 /// through the RecPoints analysis
22 ///
23 /// \author D. Stocco (Torino)
24
25 #if !defined(__CINT__) || defined(__MAKECINT__)
26 // ROOT
27 #include "TString.h"
28 #include "TH1.h"
29 #include "TList.h"
30 #include "TFile.h"
31 #include "TDirectory.h"
32
33 #include "TTree.h"
34 #include "TClonesArray.h"
35
36 // STEER
37 #include "AliRun.h"
38 #include "AliRunLoader.h"
39 #include "AliLoader.h"
40 #include "AliCDBManager.h"
41
42 // tracker
43 #include "AliGeomManager.h"
44
45 // MUON
46 #include "AliMpCDB.h"
47 #include "AliMpDDL.h"
48 #include "AliMpDDLStore.h"
49 #include "AliMpConstants.h"
50
51 // trigger
52 #include "AliMUONVTriggerStore.h"
53 #include "AliMUONDigitStoreV1.h"
54 #include "AliMUONDigitMaker.h"
55 #include "AliMUONLocalTrigger.h"
56
57 // tracker
58 #include "AliMUONClusterStoreV2.h"
59 #include "AliMUONClusterFinderMLEM.h"
60 #include "AliMUONPreClusterFinder.h"
61 #include "AliMUONSimpleClusterServer.h"
62 #include "AliMUONGeometryTransformer.h"
63 #include "AliMUONVDigitStore.h"
64 #include "AliMUONDigitStoreV2R.h"
65 #include "AliMUONVCluster.h"
66
67 #endif
68
69 const Int_t kNtrigChambers = AliMpConstants::NofTriggerChambers();
70 const Int_t kNcathodes = AliMpConstants::NofCathodes();
71
72 const Int_t kNtrackChambers = AliMpConstants::NofTrackingChambers();
73
74 enum {kOnlyTracker, kOnlyTrigger, kTrackTrig};
75
76 Int_t GetPlane(Int_t ch, Int_t cath){return kNcathodes * ch + cath;}
77 void ClusterSize(TList&, AliMUONVDigit*, Int_t&, Int_t);
78
79 // Main Method
80 void TestRecPoints(TString baseDir=".", TString outDir=".", Float_t adcCut = 10., Int_t whatToTest=kTrackTrig, Int_t runNumber=0, TString cdbStorage="local://$ALICE_ROOT")
81 {
82   const Int_t kNplanes = kNtrigChambers * kNcathodes;
83   const Int_t kNslats = 18;
84   
85   TList histoListTrig;
86   
87   TH1F* trigStripMult[kNplanes];
88   TH1F* trigClusterSize[kNplanes];
89   TH1F* trigClusterMult[kNplanes];
90   TString histoName, histoTitle;
91   
92   if(whatToTest!=kOnlyTracker){
93     for(Int_t ch=0; ch<kNtrigChambers; ch++){
94       for(Int_t cath=0; cath<kNcathodes; cath++){
95         Int_t iplane = GetPlane(ch, cath);
96         histoName = "stripMult";
97         histoName.Append(Form("Ch%iCath%i",ch,cath));
98         histoTitle = Form("Strip multiplicity Ch %i Cathode %i", ch, cath);
99         trigStripMult[iplane] = new TH1F(histoName.Data(), histoTitle.Data(), kNslats, 0.-0.5, (Float_t)kNslats - 0.5);
100         trigStripMult[iplane]->SetXTitle("slat");
101         histoListTrig.Add(trigStripMult[iplane]);
102         
103         histoName = "clusterSize";
104         histoName.Append(Form("Ch%iCath%i",ch,cath));
105         histoTitle = Form("Cluster size Ch %i Cathode %i", ch, cath);
106         trigClusterSize[iplane] = new TH1F(histoName.Data(), histoTitle.Data(), 10, 0.-0.5, (Float_t)10 - 0.5);
107         trigClusterSize[iplane]->SetXTitle("cluster size");
108         histoListTrig.Add(trigClusterSize[iplane]);
109         
110         histoName = "clusterMult";
111         histoName.Append(Form("Ch%iCath%i",ch,cath));
112         histoTitle = Form("Cluster multiplicity Ch %i Cathode %i", ch, cath);
113         trigClusterMult[iplane] = new TH1F(histoName.Data(), histoTitle.Data(), kNslats, 0.-0.5, (Float_t)kNslats - 0.5);
114         trigClusterMult[iplane]->SetXTitle("slat");
115         histoListTrig.Add(trigClusterMult[iplane]);
116       }
117     }
118   }
119   
120   const Int_t kMaxDetElem = 26;
121   TList histoListTrack;
122   TH1F* trackClusterSize[kNtrackChambers];
123   TH1F* trackClusterMult[kNtrackChambers];
124   TH1F* trackClusterYield[kNtrackChambers];
125   
126   if(whatToTest!=kOnlyTrigger){
127     for(Int_t ch=0; ch<kNtrackChambers; ch++){
128       histoName = "clusterSize";
129       histoName.Append(Form("Ch%i",ch));
130       histoTitle = Form("Cluster size Ch %i", ch);
131       trackClusterSize[ch] = new TH1F(histoName.Data(), histoTitle.Data(), kMaxDetElem, 0.-0.5, (Float_t)kMaxDetElem - 0.5);
132       trackClusterSize[ch]->SetXTitle("cluster size");
133       histoListTrack.Add(trackClusterSize[ch]);
134       
135       histoName = "clusterMult";
136       histoName.Append(Form("Ch%i",ch));
137       histoTitle = Form("Cluster multiplicity vs. Charge Ch %i", ch);
138       trackClusterMult[ch] = new TH1F(histoName.Data(), histoTitle.Data(), 200, 0, 15000.);
139       trackClusterMult[ch]->SetXTitle("Charge (ADC)");
140       histoListTrack.Add(trackClusterMult[ch]);
141       
142       histoName = "clusterYield";
143       histoName.Append(Form("Ch%i",ch));
144       histoTitle = Form("Cluster yield vs. DetElem Ch %i", ch);
145       trackClusterYield[ch] = new TH1F(histoName.Data(), histoTitle.Data(), kMaxDetElem, 0.-0.5, (Float_t)kMaxDetElem - 0.5);
146       trackClusterYield[ch]->SetXTitle("Detector element #");
147       histoListTrack.Add(trackClusterYield[ch]);
148     }
149   }
150   
151   
152   // Creating Run Loader and opening RecPoints
153   TString filename = baseDir.Data();
154   filename.Append("/galice.root");
155   AliRunLoader * RunLoader = AliRunLoader::Open(filename.Data(),"MUONLoader","UPDATE");
156   if (RunLoader ==0x0) {
157     printf(">>> Error : Error Opening %s file \n",filename.Data());
158     return;
159   }
160   
161   // Loading MUON subsystem
162   AliLoader* MUONLoader = RunLoader->GetDetectorLoader("MUON");
163   if(whatToTest!=kOnlyTracker)  MUONLoader->LoadRecPoints("READ");
164   if(whatToTest!=kOnlyTrigger)  MUONLoader->LoadDigits("READ");
165   
166   TTree *treeR = 0x0, *treeD = 0x0;
167   AliMUONVTriggerStore* trigStore = 0x0;
168   AliMUONLocalTrigger* locTrg = 0x0;
169   
170   AliMUONVDigitStore* digitStoreTrack = 0x0;
171   AliMUONDigitStoreV2R digitStoreTrackCut;
172   AliMUONVCluster* cluster = 0x0;
173   
174   // Load segmentation
175   AliCDBManager::Instance()->SetDefaultStorage(cdbStorage.Data());
176   AliCDBManager::Instance()->SetRun(runNumber);
177   
178   AliMpCDB::LoadDDLStore();
179   
180   AliMUONGeometryTransformer* transformer = 0x0;
181   
182   
183   AliMUONClusterStoreV2 clusterStore;
184   AliMUONVClusterFinder* clusterFinder = 0x0;
185   
186   AliMUONSimpleClusterServer* clusterServer = 0x0;
187   
188   Int_t firstChamber(0);
189   Int_t lastChamber(9);
190   
191   if(whatToTest!=kOnlyTrigger){
192     // Import TGeo geometry
193     TString geoFilename = baseDir.Data();
194     geoFilename.Append("/geometry.root");
195     if ( ! AliGeomManager::GetGeometry() ) {
196       AliGeomManager::LoadGeometry(geoFilename);
197       if (! AliGeomManager::GetGeometry() ) {
198         printf(">>> Error : getting geometry from file %s failed\n", geoFilename.Data());
199         return;
200       }
201     }
202     transformer = new AliMUONGeometryTransformer();
203     // Load geometry data
204     transformer->LoadGeometryData();
205     clusterFinder = new AliMUONClusterFinderMLEM(kFALSE,new AliMUONPreClusterFinder);
206     clusterServer = new AliMUONSimpleClusterServer(clusterFinder,*transformer);
207   }
208   
209   AliMUONDigitStoreV1 digitStore;
210   AliMUONVDigit* mDigit;
211   
212   Int_t clusterSize;
213   
214   AliMUONDigitMaker digitMaker;
215   TList digitsList[kNplanes];
216   
217   Int_t nevents = RunLoader->GetNumberOfEvents();
218   Int_t analysisFrac = nevents/10+1;
219   
220   printf("\nNumber of events = %i\n\n",nevents);
221   
222   for(Int_t ievent=0; ievent<nevents; ievent++){
223     if(ievent%analysisFrac==0) printf("Analysing event = %i\n",ievent);
224     RunLoader->GetEvent(ievent);
225     if(whatToTest!=kOnlyTracker){
226       digitStore.Clear();
227       treeR = MUONLoader->TreeR();
228       trigStore = AliMUONVTriggerStore::Create(*treeR);
229       
230       if ( trigStore == 0x0 ) continue;
231       trigStore->Clear();
232       trigStore->Connect(*treeR);
233       treeR->GetEvent(0);
234       
235       TIter nextLocal(trigStore->CreateLocalIterator());
236       while ( (locTrg = static_cast<AliMUONLocalTrigger*>( nextLocal() )) != NULL )
237       {
238         TArrayS xyPattern[2];
239         locTrg->GetXPattern(xyPattern[0]);
240         locTrg->GetYPattern(xyPattern[1]);
241         
242         Int_t nBoard = locTrg->LoCircuit();
243         digitMaker.TriggerDigits(nBoard, xyPattern, digitStore);
244       }
245       
246       TIter next(digitStore.CreateIterator());
247       
248       while ( ( mDigit = static_cast<AliMUONVDigit*>(next()) ) )
249       {
250         Int_t detElemId = mDigit->DetElemId();
251         Int_t ch = detElemId/100 - 11;
252         Int_t cathode = mDigit->Cathode();
253         Int_t slat = detElemId%100;
254         Int_t iplane = GetPlane(ch, cathode);
255         trigStripMult[iplane]->Fill(slat);
256         digitsList[iplane].Add(mDigit);
257       } // loop on digits
258       
259       for(Int_t iplane=0; iplane<kNplanes; iplane++){
260         while(digitsList[iplane].GetEntries()){
261           clusterSize=1;
262           mDigit = (AliMUONVDigit*)digitsList[iplane].At(0);
263           digitsList[iplane].Remove(mDigit);
264           ClusterSize(digitsList[iplane],mDigit,clusterSize,1);
265           //if(clusterSize>1) printf("Cluster size = %i\n\n",clusterSize);
266           trigClusterSize[iplane]->Fill(clusterSize);
267           
268           Int_t detElemId = mDigit->DetElemId();
269           Int_t slat = detElemId%100;
270           trigClusterMult[iplane]->Fill(slat);
271         } // loop o sorted digits
272       } // loop on planes
273     } // trigger part
274     
275     if(whatToTest!=kOnlyTrigger){
276       clusterStore.Clear();
277       treeD = MUONLoader->TreeD();
278       digitStoreTrack = AliMUONVDigitStore::Create(*treeD);
279       
280       if ( digitStoreTrack == 0x0 ) continue;
281       digitStoreTrack->Clear();
282       digitStoreTrackCut.Clear();
283       digitStoreTrack->Connect(*treeD);
284       treeD->GetEvent(0);
285       
286       // Cut low charge channels (pedestal subtraction)
287       TIter nextDigitTrack(digitStoreTrack->CreateIterator());
288       while ( ( mDigit = static_cast<AliMUONVDigit*>(nextDigitTrack()) ) )
289       {
290         //printf("Digit class = %s",mDigit->ClassName());
291         Float_t charge = mDigit->Charge();
292         if(charge<adcCut) continue;
293         digitStoreTrackCut.Add(*mDigit, AliMUONVDigitStore::kDeny);
294       } // loop on digits
295       
296       TIter nextDigitTrackCut(digitStoreTrackCut.CreateIterator());
297       clusterServer->UseDigits(nextDigitTrackCut);
298       
299       for (Int_t ch = firstChamber; ch <= lastChamber; ++ch ){
300         clusterServer->Clusterize(ch, clusterStore, AliMpArea());
301       }
302       
303       TIter next(clusterStore.CreateIterator());
304       while ( ( cluster = static_cast<AliMUONVCluster*>(next()) ) )
305       {
306         Float_t charge = cluster->GetCharge();
307         if(charge==0.) continue;
308         Int_t ch = cluster->GetChamberId();
309         Int_t npads = cluster->GetNDigits();
310         Int_t detElemId = cluster->GetDetElemId()%100;
311         //printf("charge = %f   pads = %i   detElemId = %i\n",charge,npads,detElemId);
312         
313         trackClusterSize[ch]->Fill(npads);
314         trackClusterMult[ch]->Fill(charge);
315         trackClusterYield[ch]->Fill(detElemId);      
316       } // loop on clusters
317     } // tracker part
318   } // loop on events
319   
320   MUONLoader->UnloadRecPoints();
321   MUONLoader->UnloadDigits();
322   RunLoader->UnloadAll();
323   delete RunLoader;
324   TString outFileName = outDir.Data();
325   outFileName.Append("/outTestRecPoints.root");
326   TFile* outFile = new TFile(outFileName.Data(), "RECREATE");
327   TDirectory* dir = 0x0;
328   if(whatToTest!=kOnlyTracker){
329     outFile->cd();
330     dir = outFile->mkdir("Trigger");
331     dir->cd();
332     histoListTrig.Write();
333   }
334   if(whatToTest!=kOnlyTrigger){
335     outFile->cd();
336     dir = outFile->mkdir("Tracker");
337     dir->cd();
338     histoListTrack.Write();
339   }
340   outFile->Close();
341   printf("\nSee results in %s\n",outFileName.Data());
342 }
343
344 void ClusterSize(TList& list, AliMUONVDigit* refDigit, Int_t& clusterSize, Int_t depthLevel)
345 {
346   AliMUONVDigit* mDigit = 0x0;
347   for(Int_t idigit=0; idigit<list.GetEntries(); idigit++){
348     mDigit = (AliMUONVDigit*) list.At(idigit);
349     if(mDigit->DetElemId() != refDigit->DetElemId()) continue;
350     Int_t diffX = TMath::Abs(mDigit->PadX() - refDigit->PadX());
351     Int_t diffY = TMath::Abs(mDigit->PadY() - refDigit->PadY());
352     if(diffX>1) continue;
353     if(diffY>1) continue;
354     if(diffX*diffY != 0) continue;
355     clusterSize++;
356     list.Remove(mDigit);
357     //printf("DetElemId = %i   Level = %i   Ref. (%2i,%2i)    Found (%2i,%2i)\n",mDigit->DetElemId(),depthLevel,refDigit->PadX(),refDigit->PadY(),mDigit->PadX(),mDigit->PadY());
358     ClusterSize(list, mDigit, clusterSize, depthLevel+1);
359     Int_t val = idigit + depthLevel - clusterSize;
360     //printf("Level = %i   Current digit = %i\t",depthLevel,idigit);
361     idigit = TMath::Max(-1,val);
362     //printf("New digit = %i\n",idigit);
363   }
364 }