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