]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONDigitizer.cxx
Update rawdata format for trigger (Christian)
[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 "AliMUONDigitizer.h"
19 #include "AliMUONConstants.h"
20 #include "AliMUONSegmentation.h"
21 #include "AliMUONHitMapA1.h"
22 #include "AliMUON.h"
23 #include "AliMUONLoader.h"
24 #include "AliMUONConstants.h"
25 #include "AliMUONTransientDigit.h"
26 #include "AliMUONTriggerDecision.h"
27 #include "AliLog.h"
28 #include "AliMUONGeometryTransformer.h"
29 #include "AliMUONGeometryModule.h"
30 #include "AliMUONGeometryStore.h"
31
32 #include "AliRun.h"
33 #include "AliRunDigitizer.h"
34 #include "AliRunLoader.h"
35
36 /////////////////////////////////////////////////////////////////////////////////////
37 //
38 // AliMUONDigitizer should be base abstract class of all digitisers in the MUON 
39 // module. It implements the common functionality of looping over input streams
40 // filling the fTDList and writing the fTDList to the output stream. 
41 // Inheriting digitizers need to override certain methods to choose and initialize
42 // the correct input and output trees, apply the correct detector response if any
43 // and implement how the transient digits are generated from the input stream.
44 //
45 /////////////////////////////////////////////////////////////////////////////////////
46
47 ClassImp(AliMUONDigitizer)
48
49 //___________________________________________
50 AliMUONDigitizer::AliMUONDigitizer() : 
51         AliDigitizer(),
52         fHitMap(0),
53         fTDList(0),
54         fTDCounter(0),
55         fMask(0),
56         fSignal(0)
57 {
58 /// Default constructor.
59 /// Initializes all pointers to NULL.
60
61         fRunLoader = NULL;
62         fGime = NULL;
63         fMUON = NULL;
64         fMUONData = NULL;
65         fTrigDec = NULL;
66 }
67
68 //___________________________________________
69 AliMUONDigitizer::AliMUONDigitizer(AliRunDigitizer* manager) : 
70         AliDigitizer(manager),
71         fHitMap(0),
72         fTDList(0),
73         fTDCounter(0),
74         fMask(0),
75         fSignal(0)
76 {
77 /// Constructor which should be used rather than the default constructor.
78 /// Initializes all pointers to NULL.
79
80         fRunLoader = NULL;
81         fGime = NULL;
82         fMUON = NULL;
83         fMUONData = NULL;
84         fTrigDec = NULL;
85 }
86
87 //___________________________________________
88 AliMUONDigitizer::AliMUONDigitizer(const AliMUONDigitizer& rhs)
89   : AliDigitizer(rhs)
90 {
91 /// Protected copy constructor
92
93   AliFatal("Not implemented.");
94 }
95
96 //___________________________________________
97 AliMUONDigitizer::~AliMUONDigitizer() 
98 {
99 /// Destructor
100
101   if (fMUONData)
102     delete fMUONData;
103
104   if (fTrigDec)
105     delete fTrigDec;
106 }
107
108 //-------------------------------------------------------------------
109 AliMUONDigitizer&  
110 AliMUONDigitizer::operator=(const AliMUONDigitizer& rhs)
111 {
112 /// Protected assignement operator
113
114   if (this == &rhs) return *this;
115
116   AliFatal("Not implemented.");
117     
118   return *this;  
119 }    
120          
121 //------------------------------------------------------------------------
122 Bool_t AliMUONDigitizer::Init()
123 {
124 /// Initialize - Does nothing. 
125
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             if (q > 0) {
302               Int_t chamber = td->Chamber();
303               if (0 <= chamber && chamber <= 13 )
304                 AddDigit(td, q, digitindex[chamber]++);
305               else
306                 AliError(Form("Invalid chamber %d\n",chamber));
307             }
308           }
309           FillOutputData();
310           //    }
311           fTDCounter = 0;
312 }
313
314 //------------------------------------------------------------------------
315 void AliMUONDigitizer::AddDigit(
316                 AliMUONTransientDigit* td, Int_t responseCharge,
317                 Int_t digitindex
318         )
319 {
320 /// Prepares the digits, track and charge arrays in preparation for a call to
321 /// AddDigit(Int_t, Int_t[kMAXTRACKS], Int_t[kMAXTRACKS], Int_t[6])
322 /// This method is called by CreateDigits() whenever a new digit needs to be added
323 /// to the output stream trees.
324 /// The responseCharge value is used as the Signal of the new digit. 
325 /// The OnWriteTransientDigit method is also called just before the adding the 
326 /// digit to allow inheriting digitizers to be able to do some specific processing 
327 /// at this point. 
328
329         Int_t tracks[kMAXTRACKS];
330         Int_t charges[kMAXTRACKS];
331         Int_t digits[7];
332       
333         digits[0] = td->PadX();
334         digits[1] = td->PadY();
335         digits[2] = td->Cathode() - 1;
336         digits[3] = responseCharge;
337         digits[4] = td->Physics();
338         digits[5] = td->Hit();
339         digits[6] =  td->DetElemId();
340
341         Int_t nptracks = td->GetNTracks();
342         if (nptracks > kMAXTRACKS) {
343
344           AliDebug(1, Form(
345                            "TransientDigit returned the number of tracks to be %d, which is bigger than kMAXTRACKS.",
346                            nptracks));
347           AliDebug(1, Form("Reseting the number of tracks to be %d.", kMAXTRACKS));
348           nptracks = kMAXTRACKS;
349         }
350         
351         for (Int_t i = 0; i < nptracks; i++) {
352
353           tracks[i]   = td->GetTrack(i);
354           charges[i]  = td->GetCharge(i);
355         }
356
357         // Sort list of tracks according to charge
358         SortTracks(tracks,charges,nptracks);
359
360         if (nptracks < kMAXTRACKS ) {
361
362           for (Int_t i = nptracks; i < kMAXTRACKS; i++) {
363             tracks[i]  = -1;
364             charges[i] = 0;
365           }
366         }
367
368         AliDebug(4,Form( "Adding digit with charge %d.", responseCharge));
369
370         OnWriteTransientDigit(td);
371         AddDigit(td->Chamber(), tracks, charges, digits);
372         AddDigitTrigger(td->Chamber(), tracks, charges, digits, digitindex);
373 }
374
375 //------------------------------------------------------------------------
376 void AliMUONDigitizer::OnCreateTransientDigit(AliMUONTransientDigit* /*digit*/, TObject* /*source_object*/)
377 {
378 /// Does nothing.
379 ///
380 /// This is derived by Digitisers that want to trace which digits were made from
381 /// which hits.
382 }
383
384 //------------------------------------------------------------------------
385 void AliMUONDigitizer::OnWriteTransientDigit(AliMUONTransientDigit* /*digit*/)
386 {
387 /// Does nothing.
388 ///
389 /// This is derived by Digitisers that want to trace which digits were made from
390 /// which hits.
391 }
392
393 //------------------------------------------------------------------------
394 Bool_t AliMUONDigitizer::FetchLoaders(const char* foldername, AliRunLoader*& runloader, AliMUONLoader*& muonloader)
395 {
396 /// Fetches the run loader from the current folder, specified by 'foldername'. 
397 /// The muon loader is then loaded from the fetched run loader. 
398 /// kTRUE is returned if no error occurred otherwise kFALSE is returned. 
399
400         AliDebug(3, Form("Fetching run loader and muon loader from folder: %s", foldername));
401
402         runloader = AliRunLoader::GetRunLoader(foldername);
403         if (runloader == NULL)
404         {
405                 AliError(Form("RunLoader not found in folder: %s", foldername));
406                 return kFALSE; 
407         }                                                                                       
408         muonloader = (AliMUONLoader*) runloader->GetLoader("MUONLoader");
409         if (muonloader == NULL) 
410         {
411                 AliError(Form("MUONLoader not found in folder: %s", foldername));
412                 return kFALSE; 
413         }
414         return kTRUE;
415
416 }
417
418 //------------------------------------------------------------------------
419 Bool_t AliMUONDigitizer::FetchGlobalPointers(AliRunLoader* runloader)
420 {
421 /// Fetches the AliRun object into the global gAlice pointer from the specified
422 /// run loader. The AliRun object is loaded into memory using the run loader if 
423 /// not yet loaded. The MUON module object is then loaded from gAlice and 
424 /// AliMUONData fetched from the MUON module. 
425 /// kTRUE is returned if no error occurred otherwise kFALSE is returned. 
426
427         AliDebug(3, Form("Fetching gAlice, MUON module and AliMUONData from runloader 0x%X.",
428                         (void*)runloader
429                     ));
430
431         if (runloader->GetAliRun() == NULL) runloader->LoadgAlice();
432         gAlice = runloader->GetAliRun();
433         if (gAlice == NULL)
434         {
435                 AliError(Form("Could not find the AliRun object in runloader 0x%X.", (void*)runloader));
436                 return kFALSE;
437         }
438         fMUON = (AliMUON*) gAlice->GetDetector("MUON");
439         if (fMUON == NULL)
440         {
441                 AliError(Form("Could not find the MUON module in runloader 0x%X.", (void*)runloader));
442                 return kFALSE;
443         }
444
445         AliMUONLoader *muonloader = (AliMUONLoader*) runloader->GetLoader("MUONLoader");
446         if (muonloader == NULL) 
447         {
448                 AliError( "MUONLoader not found ");
449                 return kFALSE; 
450         }
451
452
453         if (fMUONData == NULL) fMUONData = new AliMUONData(muonloader,"MUON","MUON");
454         if (fMUONData == NULL)
455         {
456                 AliError(Form("Could not find AliMUONData object in runloader 0x%X.", (void*)runloader));
457                 return kFALSE;
458         }
459
460         return kTRUE;
461 }
462 //-----------------------------------------------------------------------
463 Bool_t  AliMUONDigitizer::FetchTriggerPointer(AliMUONLoader* loader)
464 {
465 /// \todo add description
466
467   if (fMUONData == NULL) {
468     AliError("MUONData not found");
469     return kFALSE;
470   }
471
472   if (fTrigDec == NULL) 
473       fTrigDec = new AliMUONTriggerDecision(loader,0,fMUONData);
474   
475   return kTRUE;
476 }
477
478 //------------------------------------------------------------------------
479 void AliMUONDigitizer::InitArrays()
480 {
481 /// Creates a new fTDList object. 
482 /// Also creates an array of 2 * chamber_number AliMUONHitMapA1 objects
483 /// in the fHitMaps array. Each one is set to a chamber and cathode
484 /// specific segmentation model. 
485 ///
486 /// Note: the fTDList and fHitMap arrays must be NULL before calling this method.
487
488     AliDebug(2, "Initialising internal arrays.");
489     AliDebug(4, "Creating transient digits list.");
490     fTDList = new TObjArray;
491         
492     // Array of pointer of the AliMUONHitMapA1:
493     //  two HitMaps per chamber, or one HitMap per cahtode plane
494     fHitMap = new AliMUONHitMapA1* [2*AliMUONConstants::NDetElem()];
495     for (Int_t i=0; i<2*AliMUONConstants::NDetElem(); i++) fHitMap[i] = 0x0;
496
497     Int_t k = 0;
498     Int_t idDE;
499
500     for (Int_t i = 0; i < AliMUONConstants::NCh(); i++) {
501
502
503       AliDebug(4,Form( "Creating hit map for chamber %d, cathode 1.", i+1));
504       AliMUONSegmentation* segmentation = fMUON->GetSegmentation();
505       AliMUONGeometrySegmentation* c1Segmentation 
506         = segmentation->GetModuleSegmentation(i, 0); // Cathode plane 1
507       AliDebug(4,Form( "Creating hit map for chamber %d, cathode 2.", i+1));
508       AliMUONGeometrySegmentation* c2Segmentation 
509         = segmentation->GetModuleSegmentation(i, 1); // Cathode plane 2
510
511       const AliMUONGeometryTransformer* kGeometryTransformer 
512         = fMUON->GetGeometryTransformer();
513  
514       AliMUONGeometryStore* detElements 
515         = kGeometryTransformer->GetModuleTransformer(i)->GetDetElementStore();
516     
517
518     // Loop over detection elements
519       for (Int_t j=0; j<detElements->GetNofEntries(); j++) {
520        
521         idDE = detElements->GetEntry(j)->GetUniqueID();
522         fNDetElemId[idDE] = k;
523
524         Int_t npx1 = c1Segmentation->Npx(idDE)+1;
525         Int_t npy1 = c1Segmentation->Npy(idDE)+1;
526         fHitMap[k] = new AliMUONHitMapA1(npx1, npy1, fTDList); 
527      
528         Int_t npx2 = c2Segmentation->Npx(idDE)+1;
529         Int_t npy2 = c2Segmentation->Npy(idDE)+1;
530         fHitMap[k+AliMUONConstants::NDetElem()] = new AliMUONHitMapA1(npx2, npy2, fTDList);
531         k++;
532       }
533     }
534 }
535 //------------------------------------------------------------------------
536 void AliMUONDigitizer::CleanupArrays()
537 {
538 /// The arrays fTDList and fHitMap are deleted and the pointers set to NULL.
539
540         AliDebug(2, "Deleting internal arrays.");
541         for(Int_t i = 0; i < 2*AliMUONConstants::NDetElem(); i++) {
542                 delete fHitMap[i];
543         }
544         delete [] fHitMap;
545         fHitMap = NULL;
546         
547         AliDebug(4, "Deleting transient digits list.");
548         fTDList->Delete();
549         delete fTDList;
550         fTDList = NULL;
551
552 }
553
554 //------------------------------------------------------------------------
555 void AliMUONDigitizer::SortTracks(Int_t *tracks, Int_t *charges, Int_t ntr) const
556 {
557 /// Sort the list of tracks contributing to a given digit
558 /// Only the 3 most significant tracks are actually sorted
559
560        if (ntr <= 1) return;
561
562        //
563        //  Loop over signals, only 3 times
564        //
565
566        Int_t qmax;
567        Int_t jmax;
568        Int_t idx[3] = {-2,-2,-2};
569        Int_t jch[3] = {-2,-2,-2};
570        Int_t jtr[3] = {-2,-2,-2};
571        Int_t i, j, imax;
572
573        if (ntr < 3) imax = ntr;
574        else imax=3;
575         
576        for(i = 0; i < imax; i++)
577          {
578            qmax=0;
579            jmax=0;
580
581            for(j = 0; j < ntr; j++)
582              {
583                if (     (i == 1 && j == idx[i-1]) || 
584                         (i == 2 && (j == idx[i-1] || j == idx[i-2]))
585                         ) 
586                  continue;
587
588                if(charges[j] > qmax) 
589                  {
590                    qmax = charges[j];
591                    jmax = j;
592                  }       
593              } 
594
595            if(qmax > 0)
596              {
597                idx[i] = jmax;
598                jch[i] = charges[jmax]; 
599                jtr[i] = tracks[jmax]; 
600              }
601
602          } 
603
604        for(i = 0; i < 3; i++)
605          {
606            if (jtr[i] == -2) 
607              {
608                charges[i] = 0;
609                tracks[i] = 0;
610              } 
611            else 
612              {
613                charges[i] = jch[i];
614                tracks[i] = jtr[i];
615              }
616          }
617 }