]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONDigitizer.cxx
Correcting a bug (found by Sasha) that was leading to a pad being picked up twice...
[u/mrichter/AliRoot.git] / MUON / AliMUONDigitizer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2000, 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$ */
17
18 #include <assert.h>
19
20 #include "AliRun.h"
21 #include "AliRunDigitizer.h"
22 #include "AliRunLoader.h"
23
24 #include "AliMUONDigitizer.h"
25 #include "AliMUONConstants.h"
26 #include "AliMUONSegmentation.h"
27 #include "AliMUONHitMapA1.h"
28 #include "AliMUON.h"
29 #include "AliMUONLoader.h"
30 #include "AliMUONConstants.h"
31 #include "AliMUONTransientDigit.h"
32 #include "AliMUONTriggerDecision.h"
33 #include "AliLog.h"
34 #include "AliMUONGeometryTransformer.h"
35 #include "AliMUONGeometryModule.h"
36 #include "AliMUONGeometryStore.h"
37
38 /////////////////////////////////////////////////////////////////////////////////////
39 //
40 // AliMUONDigitizer should be base abstract class of all digitisers in the MUON 
41 // module. It implements the common functionality of looping over input streams
42 // filling the fTDList and writing the fTDList to the output stream. 
43 // Inheriting digitizers need to override certain methods to choose and initialize
44 // the correct input and output trees, apply the correct detector response if any
45 // and implement how the transient digits are generated from the input stream.
46 //
47 /////////////////////////////////////////////////////////////////////////////////////
48
49 ClassImp(AliMUONDigitizer)
50
51 //___________________________________________
52 AliMUONDigitizer::AliMUONDigitizer() : 
53         AliDigitizer(),
54         fHitMap(0),
55         fTDList(0),
56         fTDCounter(0),
57         fMask(0),
58         fSignal(0)
59 {
60 // Default constructor.
61 // Initializes all pointers to NULL.
62
63         fRunLoader = NULL;
64         fGime = NULL;
65         fMUON = NULL;
66         fMUONData = NULL;
67         fTrigDec = NULL;
68 }
69
70 //___________________________________________
71 AliMUONDigitizer::AliMUONDigitizer(AliRunDigitizer* manager) : 
72         AliDigitizer(manager),
73         fHitMap(0),
74         fTDList(0),
75         fTDCounter(0),
76         fMask(0),
77         fSignal(0)
78 {
79 // Constructor which should be used rather than the default constructor.
80 // Initializes all pointers to NULL.
81
82         fRunLoader = NULL;
83         fGime = NULL;
84         fMUON = NULL;
85         fMUONData = NULL;
86         fTrigDec = NULL;
87 }
88
89 //___________________________________________
90 AliMUONDigitizer::AliMUONDigitizer(const AliMUONDigitizer& rhs)
91   : AliDigitizer(rhs)
92 {
93 // Protected copy constructor
94
95   AliFatal("Not implemented.");
96 }
97
98 //___________________________________________
99 AliMUONDigitizer::~AliMUONDigitizer() 
100 {
101 // Destructor
102   if (fMUONData)
103     delete fMUONData;
104
105   if (fTrigDec)
106     delete fTrigDec;
107 }
108
109 //-------------------------------------------------------------------
110 AliMUONDigitizer&  
111 AliMUONDigitizer::operator=(const AliMUONDigitizer& rhs)
112 {
113 // Protected assignement operator
114
115   if (this == &rhs) return *this;
116
117   AliFatal("Not implemented.");
118     
119   return *this;  
120 }    
121          
122 //------------------------------------------------------------------------
123 Bool_t AliMUONDigitizer::Init()
124 {
125 // Does nothing. 
126         return kTRUE;
127 }
128
129 //------------------------------------------------------------------------
130 void AliMUONDigitizer::Exec(Option_t* /*option*/)
131 {
132 // The main work loop starts here. 
133 // The digitization process is broken up into two steps: 
134 // 1) Loop over input streams and create transient digits from the input.
135 //    Done in GenerateTransientDigits()
136 // 2) Loop over the generated transient digits and write them to the output
137 //    stream. Done in CreateDigits()
138
139         AliDebug(1, "Running digitiser.");
140
141         if (fManager->GetNinputs() == 0)
142         {
143                 AliWarning("No inputs set, nothing to do.");
144                 return;
145         }
146
147         if (!FetchLoaders(fManager->GetInputFolderName(0), fRunLoader, fGime) ) return;
148         if (! FetchGlobalPointers(fRunLoader) ) return;
149         if (! FetchTriggerPointer(fGime) ) return;
150
151         InitArrays();
152         
153         AliDebug(2, Form("Event Number is %d.", fManager->GetOutputEventNr()));
154
155         // Loop over files to merge and to digitize
156         fSignal = kTRUE;
157         for (Int_t inputFile = 0; inputFile < fManager->GetNinputs(); inputFile++)
158         {
159                 fMask = fManager->GetMask(inputFile);
160                 AliDebug(2, Form("Digitising folder %d, with fMask = %d: %s", inputFile, fMask,
161                              (const char*)fManager->GetInputFolderName(inputFile)));
162
163                 if (inputFile != 0)
164                         // If this is the first file then we already have the loaders loaded.
165                         if (! FetchLoaders(fManager->GetInputFolderName(inputFile), fRunLoader, fGime) )
166                                 continue;
167                 else
168                         // If this is not the first file then it is assumed to be background.
169                         fSignal = kFALSE;
170
171                 if (! InitInputData(fGime) ) continue;
172                 GenerateTransientDigits();
173                 CleanupInputData(fGime);
174         }
175
176         Bool_t ok = FetchLoaders(fManager->GetOutputFolderName(), fRunLoader, fGime);
177         if (ok) ok = InitOutputData(fGime);
178         if (ok) CreateDigits();
179         if (ok) CreateTrigger();
180         if (ok) CleanupOutputData(fGime);
181
182         CleanupArrays();
183         CleanupTriggerArrays();
184 }
185
186 //--------------------------------------------------------------------------
187 void AliMUONDigitizer::AddOrUpdateTransientDigit(AliMUONTransientDigit* mTD)
188 {
189 // Checks to see if the transient digit exists in the corresponding fHitMap.
190 // If it does then the digit is updated otherwise it is added. 
191
192         if (ExistTransientDigit(mTD)) 
193         {
194                 UpdateTransientDigit(mTD);
195                 delete mTD;  // The new digit can be deleted because it was not added.
196         }
197         else 
198                 AddTransientDigit(mTD);
199 }
200
201 //------------------------------------------------------------------------
202 void AliMUONDigitizer::UpdateTransientDigit(AliMUONTransientDigit* mTD)
203 {
204 // Update the transient digit that is already in the fTDList by adding the new
205 // transient digits charges and track lists to the existing one.
206
207         AliDebug(4,Form( "Updating transient digit 0x%X", (void*)mTD));
208         // Choosing the maping of the cathode plane of the chamber:
209         Int_t detElemId =  mTD->DetElemId();
210
211         Int_t iNchCpl= fNDetElemId[detElemId] + (mTD->Cathode()-1) * AliMUONConstants::NDetElem();
212
213         AliMUONTransientDigit *pdigit = 
214                 static_cast<AliMUONTransientDigit*>( fHitMap[iNchCpl]->GetHit(mTD->PadX(),mTD->PadY()) );
215
216         // update charge
217         pdigit->AddSignal( mTD->Signal() );
218         pdigit->AddPhysicsSignal( mTD->Physics() );  
219         
220         // update list of tracks
221         Int_t ntracks = mTD->GetNTracks();
222         if (ntracks > kMAXTRACKS)  // Truncate the number of tracks to kMAXTRACKS if we have to.
223         {
224                 AliDebug(1,Form( 
225                         "TransientDigit returned the number of tracks to be %d, which is bigger than kMAXTRACKS.",
226                         ntracks));
227                 AliDebug(1,Form( "Reseting the number of tracks to be %d.", kMAXTRACKS));
228                 ntracks = kMAXTRACKS;
229         }
230         
231         for (Int_t i = 0; i < ntracks; i++)
232         {
233                 pdigit->UpdateTrackList( mTD->GetTrack(i), mTD->GetCharge(i) );
234         }
235 }
236
237 //------------------------------------------------------------------------
238 void AliMUONDigitizer::AddTransientDigit(AliMUONTransientDigit* mTD)
239 {
240 // Adds the transient digit to the fTDList and sets the appropriate entry
241 // in the fHitMap arrays. 
242
243         AliDebug(4,Form( "Adding transient digit 0x%X", (void*)mTD));
244         // Choosing the maping of the cathode plane of the chamber:
245
246         Int_t detElemId =  mTD->DetElemId();
247         Int_t iNchCpl= fNDetElemId[detElemId] + (mTD->Cathode()-1) * AliMUONConstants::NDetElem();
248
249         fTDList->AddAtAndExpand(mTD, fTDCounter);
250         if (iNchCpl>-1 && iNchCpl<2*AliMUONConstants::NDetElem()) {
251           fHitMap[iNchCpl]->SetHit( mTD->PadX(), mTD->PadY(), fTDCounter);
252           fTDCounter++;
253         }
254 }
255
256 //------------------------------------------------------------------------
257 Bool_t AliMUONDigitizer::ExistTransientDigit(AliMUONTransientDigit* mTD)
258 {
259 // Checks if the transient digit already exists on the corresponding fHitMap.
260 // i.e. is there a transient digit on the same chamber, cathode and pad position 
261 // as mTD. If yes then kTRUE is returned else kFASLE is returned.
262
263         // Choosing the maping of the cathode plane of the chamber:
264         Int_t detElemId =  mTD->DetElemId();
265
266         Int_t iNchCpl= fNDetElemId[detElemId] + (mTD->Cathode()-1) *AliMUONConstants::NDetElem() ;
267
268         //      Int_t iNchCpl= mTD->Chamber() + (mTD->Cathode()-1) * AliMUONConstants::NCh();
269         if (iNchCpl>-1 && iNchCpl<2*AliMUONConstants::NDetElem())
270           return( fHitMap[iNchCpl]->TestHit(mTD->PadX(), mTD->PadY()) );
271         else return kFALSE;
272 }
273
274 //-----------------------------------------------------------------------
275 void AliMUONDigitizer::CreateDigits()
276 {
277 // Loops over the fTDList for each cathode, gets the correct signal for the 
278 // digit and adds the new digit to the output stream.
279
280         fTDList->Sort(); // sort by idDE
281         AliDebug(2, "Creating digits...");
282         //      for (Int_t icat = 0; icat < 2; icat++) {
283
284           Int_t digitindex[14] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
285
286           //
287           // Filling Digit List
288           Int_t nentries = fTDList->GetEntriesFast();
289           for (Int_t nent = 0; nent < nentries; nent++) {
290
291             AliMUONTransientDigit* td = (AliMUONTransientDigit*)fTDList->At(nent);
292             if (td == NULL) continue; 
293                         
294             // Must be the same cathode, otherwise we will fill a mixture
295             // of digits from both cathodes.
296             //if (icat != td->Cathode() - 1) continue;
297                         
298             AliDebug(3,Form( "Creating digit from transient digit 0x%X", (void*)td));
299
300             Int_t q = GetSignalFrom(td);
301             assert( 0 <= td->Chamber() && td->Chamber() <= 13 );
302             if (q > 0) AddDigit(td, q, digitindex[td->Chamber()]++);
303           }
304           FillOutputData();
305           //    }
306           fTDCounter = 0;
307 }
308
309 //------------------------------------------------------------------------
310 void AliMUONDigitizer::AddDigit(
311                 AliMUONTransientDigit* td, Int_t responseCharge,
312                 Int_t digitindex
313         )
314 {
315 // Prepares the digits, track and charge arrays in preparation for a call to
316 // AddDigit(Int_t, Int_t[kMAXTRACKS], Int_t[kMAXTRACKS], Int_t[6])
317 // This method is called by CreateDigits() whenever a new digit needs to be added
318 // to the output stream trees.
319 // The responseCharge value is used as the Signal of the new digit. 
320 // The OnWriteTransientDigit method is also called just before the adding the 
321 // digit to allow inheriting digitizers to be able to do some specific processing 
322 // at this point. 
323
324         Int_t tracks[kMAXTRACKS];
325         Int_t charges[kMAXTRACKS];
326         Int_t digits[7];
327       
328         digits[0] = td->PadX();
329         digits[1] = td->PadY();
330         digits[2] = td->Cathode() - 1;
331         digits[3] = responseCharge;
332         digits[4] = td->Physics();
333         digits[5] = td->Hit();
334         digits[6] =  td->DetElemId();
335
336         Int_t nptracks = td->GetNTracks();
337         if (nptracks > kMAXTRACKS) {
338
339           AliDebug(1, Form(
340                            "TransientDigit returned the number of tracks to be %d, which is bigger than kMAXTRACKS.",
341                            nptracks));
342           AliDebug(1, Form("Reseting the number of tracks to be %d.", kMAXTRACKS));
343           nptracks = kMAXTRACKS;
344         }
345         
346         for (Int_t i = 0; i < nptracks; i++) {
347
348           tracks[i]   = td->GetTrack(i);
349           charges[i]  = td->GetCharge(i);
350         }
351
352         // Sort list of tracks according to charge
353         SortTracks(tracks,charges,nptracks);
354
355         if (nptracks < kMAXTRACKS ) {
356
357           for (Int_t i = nptracks; i < kMAXTRACKS; i++) {
358             tracks[i]  = -1;
359             charges[i] = 0;
360           }
361         }
362
363         AliDebug(4,Form( "Adding digit with charge %d.", responseCharge));
364
365         OnWriteTransientDigit(td);
366         AddDigit(td->Chamber(), tracks, charges, digits);
367         AddDigitTrigger(td->Chamber(), tracks, charges, digits, digitindex);
368 }
369
370 //------------------------------------------------------------------------
371 void AliMUONDigitizer::OnCreateTransientDigit(AliMUONTransientDigit* /*digit*/, TObject* /*source_object*/)
372 {
373         // Does nothing.
374         //
375         // This is derived by Digitisers that want to trace which digits were made from
376         // which hits.
377 }
378
379 //------------------------------------------------------------------------
380 void AliMUONDigitizer::OnWriteTransientDigit(AliMUONTransientDigit* /*digit*/)
381 {
382         // Does nothing.
383         //
384         // This is derived by Digitisers that want to trace which digits were made from
385         // which hits.
386 }
387
388 //------------------------------------------------------------------------
389 Bool_t AliMUONDigitizer::FetchLoaders(const char* foldername, AliRunLoader*& runloader, AliMUONLoader*& muonloader)
390 {
391 // Fetches the run loader from the current folder, specified by 'foldername'. 
392 // The muon loader is then loaded from the fetched run loader. 
393 // kTRUE is returned if no error occurred otherwise kFALSE is returned. 
394
395         AliDebug(3, Form("Fetching run loader and muon loader from folder: %s", foldername));
396
397         runloader = AliRunLoader::GetRunLoader(foldername);
398         if (runloader == NULL)
399         {
400                 AliError(Form("RunLoader not found in folder: %s", foldername));
401                 return kFALSE; 
402         }                                                                                       
403         muonloader = (AliMUONLoader*) runloader->GetLoader("MUONLoader");
404         if (muonloader == NULL) 
405         {
406                 AliError(Form("MUONLoader not found in folder: %s", foldername));
407                 return kFALSE; 
408         }
409         return kTRUE;
410
411 }
412
413 //------------------------------------------------------------------------
414 Bool_t AliMUONDigitizer::FetchGlobalPointers(AliRunLoader* runloader)
415 {
416 // Fetches the AliRun object into the global gAlice pointer from the specified
417 // run loader. The AliRun object is loaded into memory using the run loader if 
418 // not yet loaded. The MUON module object is then loaded from gAlice and 
419 // AliMUONData fetched from the MUON module. 
420 // kTRUE is returned if no error occurred otherwise kFALSE is returned. 
421
422         AliDebug(3, Form("Fetching gAlice, MUON module and AliMUONData from runloader 0x%X.",
423                         (void*)runloader
424                     ));
425
426         if (runloader->GetAliRun() == NULL) runloader->LoadgAlice();
427         gAlice = runloader->GetAliRun();
428         if (gAlice == NULL)
429         {
430                 AliError(Form("Could not find the AliRun object in runloader 0x%X.", (void*)runloader));
431                 return kFALSE;
432         }
433         fMUON = (AliMUON*) gAlice->GetDetector("MUON");
434         if (fMUON == NULL)
435         {
436                 AliError(Form("Could not find the MUON module in runloader 0x%X.", (void*)runloader));
437                 return kFALSE;
438         }
439
440         AliMUONLoader *muonloader = (AliMUONLoader*) runloader->GetLoader("MUONLoader");
441         if (muonloader == NULL) 
442         {
443                 AliError( "MUONLoader not found ");
444                 return kFALSE; 
445         }
446
447
448         if (fMUONData == NULL) fMUONData = new AliMUONData(muonloader,"MUON","MUON");
449         if (fMUONData == NULL)
450         {
451                 AliError(Form("Could not find AliMUONData object in runloader 0x%X.", (void*)runloader));
452                 return kFALSE;
453         }
454
455         return kTRUE;
456 }
457 //-----------------------------------------------------------------------
458 Bool_t  AliMUONDigitizer::FetchTriggerPointer(AliMUONLoader* loader)
459 {
460   if (fMUONData == NULL) {
461     AliError("MUONData not found");
462     return kFALSE;
463   }
464
465   if (fTrigDec == NULL) 
466       fTrigDec = new AliMUONTriggerDecision(loader,0,fMUONData);
467   
468   return kTRUE;
469 }
470
471 //------------------------------------------------------------------------
472 void AliMUONDigitizer::InitArrays()
473 {
474 // Creates a new fTDList object. 
475 // Also creates an array of 2 * chamber_number AliMUONHitMapA1 objects
476 // in the fHitMaps array. Each one is set to a chamber and cathode
477 // specific segmentation model. 
478 //
479 // Note: the fTDList and fHitMap arrays must be NULL before calling this method.
480
481     AliDebug(2, "Initialising internal arrays.");
482     AliDebug(4, "Creating transient digits list.");
483     fTDList = new TObjArray;
484         
485     // Array of pointer of the AliMUONHitMapA1:
486     //  two HitMaps per chamber, or one HitMap per cahtode plane
487     fHitMap = new AliMUONHitMapA1* [2*AliMUONConstants::NDetElem()];
488     for (Int_t i=0; i<2*AliMUONConstants::NDetElem(); i++) fHitMap[i] = 0x0;
489
490     Int_t k = 0;
491     Int_t idDE;
492
493     for (Int_t i = 0; i < AliMUONConstants::NCh(); i++) {
494
495
496       AliDebug(4,Form( "Creating hit map for chamber %d, cathode 1.", i+1));
497       AliMUONSegmentation* segmentation = fMUON->GetSegmentation();
498       AliMUONGeometrySegmentation* c1Segmentation 
499         = segmentation->GetModuleSegmentation(i, 0); // Cathode plane 1
500       AliDebug(4,Form( "Creating hit map for chamber %d, cathode 2.", i+1));
501       AliMUONGeometrySegmentation* c2Segmentation 
502         = segmentation->GetModuleSegmentation(i, 1); // Cathode plane 2
503
504       const AliMUONGeometryTransformer* kGeometryTransformer 
505         = fMUON->GetGeometryTransformer();
506  
507       AliMUONGeometryStore* detElements 
508         = kGeometryTransformer->GetModuleTransformer(i)->GetDetElementStore();
509     
510
511     // Loop over detection elements
512       for (Int_t j=0; j<detElements->GetNofEntries(); j++) {
513        
514         idDE = detElements->GetEntry(j)->GetUniqueID();
515         fNDetElemId[idDE] = k;
516
517         fHitMap[k] = new AliMUONHitMapA1(idDE,c1Segmentation, fTDList); 
518      
519         fHitMap[k+AliMUONConstants::NDetElem()] = new AliMUONHitMapA1(idDE,c2Segmentation, fTDList);
520         k++;
521       }
522     }
523 }
524 //------------------------------------------------------------------------
525 void AliMUONDigitizer::CleanupArrays()
526 {
527 // The arrays fTDList and fHitMap are deleted and the pointers set to NULL.
528
529         AliDebug(2, "Deleting internal arrays.");
530         for(Int_t i = 0; i < 2*AliMUONConstants::NDetElem(); i++) {
531                 delete fHitMap[i];
532         }
533         delete [] fHitMap;
534         fHitMap = NULL;
535         
536         AliDebug(4, "Deleting transient digits list.");
537         fTDList->Delete();
538         delete fTDList;
539         fTDList = NULL;
540
541 }
542
543 //------------------------------------------------------------------------
544 void AliMUONDigitizer::SortTracks(Int_t *tracks, Int_t *charges, Int_t ntr) const
545 {
546 //
547 // Sort the list of tracks contributing to a given digit
548 // Only the 3 most significant tracks are actually sorted
549 //
550
551        if (ntr <= 1) return;
552
553        //
554        //  Loop over signals, only 3 times
555        //
556
557        Int_t qmax;
558        Int_t jmax;
559        Int_t idx[3] = {-2,-2,-2};
560        Int_t jch[3] = {-2,-2,-2};
561        Int_t jtr[3] = {-2,-2,-2};
562        Int_t i, j, imax;
563
564        if (ntr < 3) imax = ntr;
565        else imax=3;
566         
567        for(i = 0; i < imax; i++)
568          {
569            qmax=0;
570            jmax=0;
571
572            for(j = 0; j < ntr; j++)
573              {
574                if (     (i == 1 && j == idx[i-1]) || 
575                         (i == 2 && (j == idx[i-1] || j == idx[i-2]))
576                         ) 
577                  continue;
578
579                if(charges[j] > qmax) 
580                  {
581                    qmax = charges[j];
582                    jmax = j;
583                  }       
584              } 
585
586            if(qmax > 0)
587              {
588                idx[i] = jmax;
589                jch[i] = charges[jmax]; 
590                jtr[i] = tracks[jmax]; 
591              }
592
593          } 
594
595        for(i = 0; i < 3; i++)
596          {
597            if (jtr[i] == -2) 
598              {
599                charges[i] = 0;
600                tracks[i] = 0;
601              } 
602            else 
603              {
604                charges[i] = jch[i];
605                tracks[i] = jtr[i];
606              }
607          }
608 }