1 /**************************************************************************
2 * Copyright(c) 1998-2000, 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 Revision 1.9 2002/03/13 07:04:11 jchudoba
19 Connect only MUON branches when reading the event to speed up digitisation.
21 Revision 1.8 2002/02/22 12:14:21 morsch
22 Validate pad hit before digitization.
24 Revision 1.7 2001/11/22 11:15:41 jchudoba
25 Proper deletion of arrays (thanks to Rene Brun)
27 Revision 1.6 2001/11/02 12:55:45 jchudoba
28 cleanup of the code, add const to Get methods
30 Revision 1.5 2001/10/31 16:35:07 jchudoba
31 some functionality move to AliMUONTransientDigit class
33 Revision 1.4 2001/03/20 13:36:11 egangler
34 TFile memory leak and "too many files open" problem solved (?):
35 the backgroud file is open only once forever, and not once for each event.
37 Revision 1.3 2001/03/05 23:57:44 morsch
38 Writing of digit tree moved to macro.
40 Revision 1.2 2001/03/05 08:40:25 morsch
41 Method SortTracks(..) imported from AliMUON.
43 Revision 1.1 2001/02/02 14:11:53 morsch
44 AliMUONMerger prototype to be called by the merge manager.
48 #include <TDirectory.h>
50 #include <TObjArray.h>
54 #include "AliHitMap.h"
56 #include "AliMUONChamber.h"
57 #include "AliMUONConstants.h"
58 #include "AliMUONDigit.h"
59 #include "AliMUONHit.h"
60 #include "AliMUONHitMapA1.h"
61 #include "AliMUONMerger.h"
62 #include "AliMUONPadHit.h"
63 #include "AliMUONTransientDigit.h"
66 ClassImp(AliMUONMerger)
68 //___________________________________________
69 AliMUONMerger::AliMUONMerger()
71 // Default constructor
86 //------------------------------------------------------------------------
87 AliMUONMerger::~AliMUONMerger()
90 if (fTrH1) delete fTrH1;
91 if (fHitsBgr) delete fHitsBgr;
92 if (fPadHitsBgr) delete fPadHitsBgr;
93 if (fHitMap) delete [] fHitMap;
94 if (fList) delete fList;
95 if (fBgrFile) delete fBgrFile;
98 //------------------------------------------------------------------------
99 Bool_t AliMUONMerger::Exists(const AliMUONPadHit *padhit) const
101 // test if the given padhit was already fired
102 return (fHitMap[fNch]->TestHit(padhit->PadX(),padhit->PadY()));
105 //------------------------------------------------------------------------
106 void AliMUONMerger::Update(AliMUONPadHit *padhit)
108 // add new contribution to the fired padhit
109 AliMUONTransientDigit *pdigit =
110 static_cast<AliMUONTransientDigit*>(
111 fHitMap[fNch]->GetHit(padhit->PadX(),padhit->PadY()));
115 Int_t iqpad = padhit->QPad(); // charge per pad
116 pdigit->AddSignal(iqpad);
117 pdigit->AddPhysicsSignal(iqpad);
119 // update list of tracks
129 pdigit->UpdateTrackList(track,charge);
132 //------------------------------------------------------------------------
133 void AliMUONMerger::CreateNew(AliMUONPadHit *padhit)
135 // add new transient digit to the list, update hit map
136 fList->AddAtAndExpand(
137 new AliMUONTransientDigit(fNch,fDigits),fCounter);
138 fHitMap[fNch]->SetHit(padhit->PadX(),padhit->PadY(),fCounter);
139 AliMUONTransientDigit* pdigit =
140 static_cast<AliMUONTransientDigit*>(fList->Last());
145 charge = padhit->QPad();
150 pdigit->AddToTrackList(track,charge);
155 //------------------------------------------------------------------------
156 void AliMUONMerger::Init()
159 // open only once the background file !!
160 if (fMerge && !fBgrFile) fBgrFile = InitBgr();
165 //------------------------------------------------------------------------
166 TFile* AliMUONMerger::InitBgr()
168 // Initialise background event
169 TFile *file = new TFile(fFnBgr);
170 // add error checking later
171 printf("\n AliMUONMerger has opened %s file with background event \n", fFnBgr);
172 fHitsBgr = new TClonesArray("AliMUONHit",1000);
173 fPadHitsBgr = new TClonesArray("AliMUONPadHit",1000);
177 //------------------------------------------------------------------------
178 void AliMUONMerger::Digitise()
181 // keep galice.root for signal and name differently the file for
182 // background when add! otherwise the track info for signal will be lost !
184 AliMUONChamber* iChamber;
185 AliSegmentation* segmentation;
187 fList = new TObjArray;
189 AliMUON *pMUON = (AliMUON *) gAlice->GetModule("MUON");
190 fHitMap= new AliHitMap* [AliMUONConstants::NCh()];
194 // Get Hits Tree header from file
195 //if(fHitsBgr) fHitsBgr->Clear(); // Useless because line 327
196 //if(fPadHitsBgr) fPadHitsBgr->Clear(); // Useless because line 328
197 if(fTrH1) delete fTrH1;
201 sprintf(treeName,"TreeH%d",fEvNrBgr);
202 fTrH1 = (TTree*)gDirectory->Get(treeName);
204 printf("ERROR: cannot find Hits Tree for event:%d\n",fEvNrBgr);
207 // Set branch addresses
210 sprintf(branchname,"%s",pMUON->GetName());
211 if (fTrH1 && fHitsBgr) {
212 branch = fTrH1->GetBranch(branchname);
213 if (branch) branch->SetAddress(&fHitsBgr);
215 if (fTrH1 && fPadHitsBgr) {
216 branch = fTrH1->GetBranch("MUONCluster");
217 if (branch) branch->SetAddress(&fPadHitsBgr);
221 // loop over cathodes
224 for (int icat = 0; icat < 2; icat++) {
226 for (Int_t i = 0; i < AliMUONConstants::NCh(); i++) {
227 iChamber = &(pMUON->Chamber(i));
228 if (iChamber->Nsec() == 1 && icat == 1) {
231 segmentation = iChamber->SegmentationModel(icat+1);
233 fHitMap[i] = new AliMUONHitMapA1(segmentation, fList);
240 TTree *treeH = gAlice->TreeH();
241 Int_t ntracks = (Int_t) treeH->GetEntries();
242 treeH->SetBranchStatus("*",0); // switch off all branches
243 treeH->SetBranchStatus("MUON*",1); // switch on only MUON
245 for (fTrack = 0; fTrack < ntracks; fTrack++) {
247 treeH->GetEntry(fTrack,0);
250 for(AliMUONHit* mHit = (AliMUONHit*)pMUON->FirstHit(-1);
252 mHit = (AliMUONHit*)pMUON->NextHit())
254 fNch = mHit->Chamber()-1; // chamber number
255 if (fNch > AliMUONConstants::NCh()-1) continue;
256 iChamber = &(pMUON->Chamber(fNch));
259 // Loop over pad hits
260 for (AliMUONPadHit* mPad =
261 (AliMUONPadHit*)pMUON->FirstPad(mHit,pMUON->PadHits());
263 mPad = (AliMUONPadHit*)pMUON->NextPad(pMUON->PadHits()))
265 Int_t cathode = mPad->Cathode(); // cathode number
266 if (cathode != (icat+1)) continue;
267 Int_t iqpad = Int_t(mPad->QPad()); // charge per pad
268 // segmentation = iChamber->SegmentationModel(cathode);
269 fDigits[0] = mPad->PadX();
270 fDigits[1] = mPad->PadY();
271 if (!(fHitMap[fNch]->ValidateHit(fDigits[0], fDigits[1]))) continue;
275 if (mHit->Particle() == kMuonPlus ||
276 mHit->Particle() == kMuonMinus) {
277 fDigits[5] = mPad->HitNumber();
278 } else fDigits[5] = -1;
280 // build the list of fired pads and update the info
287 } //end loop over clusters
291 // open the file with background
295 ntracks = (Int_t)fTrH1->GetEntries();
299 for (fTrack = 0; fTrack < ntracks; fTrack++) {
301 if (fHitsBgr) fHitsBgr->Clear();
302 if (fPadHitsBgr) fPadHitsBgr->Clear();
304 fTrH1->GetEvent(fTrack);
308 for(Int_t i = 0; i < fHitsBgr->GetEntriesFast(); ++i)
310 mHit = (AliMUONHit*) (*fHitsBgr)[i];
311 fNch = mHit->Chamber()-1; // chamber number
312 iChamber = &(pMUON->Chamber(fNch));
314 // Loop over pad hits
315 for (AliMUONPadHit* mPad =
316 (AliMUONPadHit*)pMUON->FirstPad(mHit,fPadHitsBgr);
318 mPad = (AliMUONPadHit*)pMUON->NextPad(fPadHitsBgr))
320 Int_t cathode = mPad->Cathode(); // cathode number
321 Int_t ipx = mPad->PadX(); // pad number on X
322 Int_t ipy = mPad->PadY(); // pad number on Y
323 Int_t iqpad = Int_t(mPad->QPad()); // charge per pad
324 if (!(fHitMap[fNch]->ValidateHit(ipx, ipy))) continue;
326 if (cathode != (icat+1)) continue;
334 // build the list of fired pads and update the info
340 } //end loop over clusters
344 TTree *treeK = gAlice->TreeK();
347 if (treeK) file = treeK->GetCurrentFile();
351 Int_t tracks[kMAXTRACKS];
352 Int_t charges[kMAXTRACKS];
353 Int_t nentries = fList->GetEntriesFast();
355 for (Int_t nent = 0; nent < nentries; nent++) {
356 AliMUONTransientDigit *address = (AliMUONTransientDigit*)fList->At(nent);
357 if (address == 0) continue;
358 Int_t ich = address->Chamber();
359 Int_t q = address->Signal();
360 iChamber = &(pMUON->Chamber(ich));
362 // Digit Response (noise, threshold, saturation, ...)
363 AliMUONResponse * response = iChamber->ResponseModel();
364 q = response->DigitResponse(q);
368 fDigits[0] = address->PadX();
369 fDigits[1] = address->PadY();
370 fDigits[2] = address->Cathode();
372 fDigits[4] = address->Physics();
373 fDigits[5] = address->Hit();
375 Int_t nptracks = address->GetNTracks();
377 if (nptracks > kMAXTRACKS) {
378 printf("\n Attention - nptracks > kMAXTRACKS %d \n", nptracks);
379 nptracks = kMAXTRACKS;
382 printf("Attention - nptracks > 2 %d \n",nptracks);
383 printf("cat,ich,ix,iy,q %d %d %d %d %d \n",icat,ich,fDigits[0],fDigits[1],q);
385 for (Int_t tr = 0; tr < nptracks; tr++) {
386 tracks[tr] = address->GetTrack(tr);
387 charges[tr] = address->GetCharge(tr);
388 } //end loop over list of tracks for one pad
389 // Sort list of tracks according to charge
391 SortTracks(tracks,charges,nptracks);
393 if (nptracks < kMAXTRACKS ) {
394 for (Int_t i = nptracks; i < kMAXTRACKS; i++) {
401 pMUON->AddDigits(ich,tracks,charges,fDigits);
403 gAlice->TreeD()->Fill();
404 pMUON->ResetDigits();
408 for(Int_t ii = 0; ii < AliMUONConstants::NCh(); ++ii) {
414 } //end loop over cathodes
418 // no need to delete ... and it makes a crash also
419 // if (fHitsBgr) fHitsBgr->Delete();
420 // if (fPadHitsBgr) fPadHitsBgr->Delete();
421 // gObjectTable->Print();
426 void AliMUONMerger::SortTracks(Int_t *tracks,Int_t *charges,Int_t ntr)
429 // Sort the list of tracks contributing to a given digit
430 // Only the 3 most significant tracks are acctually sorted
434 // Loop over signals, only 3 times
439 Int_t idx[3] = {-2,-2,-2};
440 Int_t jch[3] = {-2,-2,-2};
441 Int_t jtr[3] = {-2,-2,-2};
452 if((i == 1 && j == idx[i-1])
453 ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
455 if(charges[j] > qmax) {
463 jch[i]=charges[jmax];