Transition to NewIO
[u/mrichter/AliRoot.git] / MUON / AliMUONDigitizerv1.cxx
1
2
3 #include <Riostream.h>
4 #include <TDirectory.h>
5 #include <TFile.h>
6 #include <TObjArray.h>
7 #include <TPDGCode.h>
8 #include <TTree.h> 
9 #include <TMath.h>
10
11 #include "AliMUON.h"
12 #include "AliMUONChamber.h"
13 #include "AliMUONConstants.h"
14 #include "AliMUONDigit.h"
15 #include "AliMUONDigitizerv1.h"
16 #include "AliMUONHit.h"
17 #include "AliMUONHitMapA1.h"
18 #include "AliMUONPadHit.h"
19 #include "AliMUONTransientDigit.h"
20 #include "AliRun.h"
21 #include "AliRunDigitizer.h"
22 #include "AliRunLoader.h"
23 #include "AliLoader.h"
24
25 ClassImp(AliMUONDigitizerv1)
26
27 //___________________________________________
28 AliMUONDigitizerv1::AliMUONDigitizerv1() :AliDigitizer()
29 {
30 // Default ctor - don't use it
31   fHitMap = 0;
32   fTDList = 0;
33   if (GetDebug()>2) 
34     cerr<<"AliMUONDigitizerv1::AliMUONDigitizerv1"
35         <<"(AliRunDigitizer* manager) was processed"<<endl;
36 }
37
38 //___________________________________________
39 AliMUONDigitizerv1::AliMUONDigitizerv1(AliRunDigitizer* manager) 
40     :AliDigitizer(manager)
41 {
42 // ctor which should be used
43   fHitMap  = 0;
44   fTDList  = 0;
45   fDebug   = 0; 
46   fHits = new TClonesArray("AliMUONHit",1000);
47   if (GetDebug()>2) 
48     cerr<<"AliMUONDigitizerv1::AliMUONDigitizerv1"
49         <<"(AliRunDigitizer* manager) was processed"<<endl;
50 }
51
52 //------------------------------------------------------------------------
53 AliMUONDigitizerv1::~AliMUONDigitizerv1()
54 {
55 // Destructor
56   delete fHits;
57 }
58
59 //------------------------------------------------------------------------
60 void AliMUONDigitizerv1::AddTransientDigit(AliMUONTransientDigit * mTD)
61 {
62   // Choosing the maping of the cathode plane of the chamber:
63   Int_t iNchCpl= mTD->Chamber() + (mTD->Cathode()-1) * AliMUONConstants::NCh();
64   fTDList->AddAtAndExpand(mTD, fTDCounter);
65   fHitMap[iNchCpl]->SetHit( mTD->PadX(), mTD->PadY(), fTDCounter);
66   fTDCounter++;
67 }
68
69 //------------------------------------------------------------------------
70 Bool_t AliMUONDigitizerv1::ExistTransientDigit(AliMUONTransientDigit * mTD)
71 {
72   // Choosing the maping of the cathode plane of the chamber:
73   Int_t iNchCpl= mTD->Chamber() + (mTD->Cathode()-1) * AliMUONConstants::NCh();
74   return( fHitMap[iNchCpl]->TestHit(mTD->PadX(), mTD->PadY()) );
75 }
76
77 //------------------------------------------------------------------------
78 Bool_t AliMUONDigitizerv1::Init()
79 {
80 // Initialization 
81   return kTRUE;
82 }
83
84 //------------------------------------------------------------------------
85 void AliMUONDigitizerv1::UpdateTransientDigit(Int_t track, AliMUONTransientDigit * mTD)
86 {
87   // Choosing the maping of the cathode plane of the chamber:
88   Int_t iNchCpl= mTD->Chamber() + (mTD->Cathode()-1) * AliMUONConstants::NCh();
89   AliMUONTransientDigit *pdigit = 
90     static_cast<AliMUONTransientDigit*>(fHitMap[iNchCpl]->GetHit(mTD->PadX(),mTD->PadY()));
91   // update charge
92   //
93   Int_t iqpad  = mTD->Signal();        // charge per pad
94   pdigit->AddSignal(iqpad);
95   pdigit->AddPhysicsSignal(iqpad);              
96   // update list of tracks
97   //
98   Int_t charge;    
99   track=+ fMask;
100   if (fSignal) charge = iqpad;
101   else         charge = kBgTag;
102   pdigit->UpdateTrackList(track,charge);
103 }
104
105
106 //--------------------------------------------------------------------------
107 void AliMUONDigitizerv1::MakeTransientDigit(Int_t track, Int_t iHit, AliMUONHit * mHit)
108 {
109   AliMUON *pMUON  = (AliMUON *) gAlice->GetModule("MUON");
110   if (!pMUON) { 
111     cerr<<"AliMUONDigitizerv1::Digitize Error:"
112         <<" module MUON not found in the input file"<<endl;
113   }
114   if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::MakeTransientDigit starts"<<endl;
115   Int_t   ichamber = mHit->Chamber()-1;
116   AliMUONChamber & chamber = pMUON->Chamber(ichamber);
117   Float_t xhit = mHit->X();
118   Float_t yhit = mHit->Y();
119   Float_t zhit = mHit->Z();
120   Float_t eloss= mHit->Eloss(); 
121   Float_t tof  = mHit->Age();
122   // Variables for chamber response from AliMUONChamber::DisIntegration
123   Float_t newdigit[6][500];  // Pad information
124   Int_t nnew=0;              // Number of touched Pads per hit
125   Int_t digits[6];
126   
127   //
128   // Calls the charge disintegration method of the current chamber 
129   if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::MakeTransientDigit calling AliMUONChamber::DisIngtegration starts"<<endl;
130   chamber.DisIntegration(eloss, tof, xhit, yhit, zhit, nnew, newdigit);
131   // Creating a new TransientDigits from hit
132   for(Int_t iTD=0; iTD<nnew; iTD++) {
133     digits[0] = Int_t(newdigit[1][iTD]);  // Padx of the Digit
134     digits[1] = Int_t(newdigit[2][iTD]);  // Pady of the Digit
135     digits[2] = Int_t(newdigit[5][iTD]);  // Cathode plane
136     digits[3] = Int_t(newdigit[3][iTD]);  // Induced charge in the Pad
137     if (fSignal) digits[4] = Int_t(newdigit[3][iTD]);
138     else         digits[4] = 0;
139     digits[5] = iHit+fMask;    // Hit number in the list
140     if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::MakeTransientDigit " <<
141                         "PadX "<< digits[0] << " " <<
142                         "PadY "<< digits[1] << " " <<
143                         "Plane " << digits[2] << " " <<
144                         "Charge " << digits[3] <<" " <<
145                         "Hit " << digits[5] << endl;
146     // list of tracks
147     Int_t charge;    
148     track += fMask;
149     if (fSignal) charge = digits[3];
150     else         charge = kBgTag;
151     if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::MakeTransientDigit Creating AliMUONTransientDigit"<<endl;
152     AliMUONTransientDigit * mTD = new AliMUONTransientDigit(ichamber, digits);
153     mTD->AddToTrackList(track,charge);
154
155     if (!ExistTransientDigit(mTD)) {
156       AddTransientDigit(mTD);
157       if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::MakeTransientDigit Adding TransientDigit"<<endl;
158     }
159     else {
160       if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::MakeTransientDigit updating TransientDigit"<<endl;
161       UpdateTransientDigit(track, mTD);
162       delete mTD;
163     }
164   }
165 }
166 //-----------------------------------------------------------------------
167 void AliMUONDigitizerv1::Exec(Option_t* option)
168 {
169   TString optionString = option;
170   if (optionString.Data() == "deb") {
171     cout<<"AliMUONDigitizerv1::Exec: called with option deb "<<endl;
172     fDebug = 3;
173   }
174   AliMUONChamber*   chamber;
175   AliSegmentation*  c1Segmentation; //Cathode plane c1 of the chamber
176   AliSegmentation*  c2Segmentation; //Cathode place c2 of the chamber
177
178   if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::Digitize() starts"<<endl;
179   fTDList = new TObjArray;
180
181   //Loaders (We assume input0 to be the output too)
182   AliRunLoader * runloader;    // Input   loader
183   AliLoader    * gime; 
184  
185   runloader    = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
186   if (runloader == 0x0) {
187     cerr<<"AliMUONDigitizerv1::Digitize() opening file "<<fManager->GetInputFileName(0,0)<<endl;
188     return;  // RunDigitizer is not working.  
189   }
190   gime         = runloader->GetLoader("MUONLoader");
191   if (gime->TreeH()==0x0) {
192      Info("Digitize","TreeH is not loaded yet. Loading...");
193      gime->LoadHits("READ");
194      Info("Digitize","Now treeH is %#x. MUONLoader is %#x",gime->TreeH(),gime);
195   }
196
197   if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::Digitize() loaders"<<endl;
198
199    if (runloader->GetAliRun() == 0x0) runloader->LoadgAlice();
200   gAlice = runloader->GetAliRun();
201
202   // Getting Module MUON  
203   AliMUON *pMUON  = (AliMUON *) gAlice->GetDetector("MUON");
204   if (!pMUON) {
205     cerr<<"AliMUONDigitizerv1::Digitize Error:"
206         <<" module MUON not found in the input file"<<endl;
207     return;
208   } 
209
210   // Loading Event
211   Int_t currentevent = fManager->GetOutputEventNr();
212   
213   if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::Digitize() Event Number is "<<currentevent <<endl;
214   if ( (currentevent<10)                                                 || 
215        (Int_t(TMath::Log10(currentevent)) == TMath::Log10(currentevent) ) )
216     cout <<"ALiMUONDigitizerv1::Digitize() Event Number is "<< currentevent <<endl;
217
218   // Output file for digits 
219   AliRunLoader *  runloaderout  = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
220   AliLoader * gimeout         = runloaderout->GetLoader("MUONLoader"); 
221   // New branch per chamber for MUON digit in the tree of digits
222   if (gime->TreeD() == 0x0) gimeout->MakeDigitsContainer();
223   TTree* treeD = gimeout->TreeD();
224   pMUON->MakeBranchInTreeD(treeD);
225
226   // Array of pointer of the AliMUONHitMapA1:
227   //  two HitMaps per chamber, or one HitMap per cahtode plane
228   fHitMap= new AliMUONHitMapA1* [2*AliMUONConstants::NCh()];
229
230   //Loop over chambers for the definition AliMUONHitMap
231   for (Int_t i=0; i<AliMUONConstants::NCh(); i++) {
232     chamber = &(pMUON->Chamber(i));
233     c1Segmentation = chamber->SegmentationModel(1); // Cathode plane 1
234     fHitMap[i] = new AliMUONHitMapA1(c1Segmentation, fTDList);
235     c2Segmentation = chamber->SegmentationModel(2); // Cathode plane 2
236     fHitMap[i+AliMUONConstants::NCh()] = new AliMUONHitMapA1(c2Segmentation, fTDList);
237   }
238
239 // Loop over files to merge and to digitize
240     fSignal = kTRUE;
241     for (Int_t inputFile=0; inputFile<fManager->GetNinputs(); inputFile++) {
242       if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::Digitize() Input File is "<<inputFile<<endl;
243
244       // Hits Branch
245       TBranch *branchHits = 0;
246
247       // Connect MUON Hit branch
248       if (inputFile > 0 ) {
249         fSignal = kFALSE;
250          runloader    = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
251          if (runloader == 0x0) {
252            cerr<<"AliMUONDigitizerv1::Digitize() RunLoader for inputFile "<<inputFile<< " not found !!! "<<endl;
253          }
254          gime  = runloader->GetLoader("MUONLoader");
255          if (gime->TreeH() == 0x0) gime->LoadHits("READ");      
256       }
257
258       // Setting the address of fHits list
259       TTree *treeH = gime->TreeH();
260       if (treeH == 0x0) {
261         Error("Digitize","Can not get TreeH from input %d",inputFile);
262         Info("Digitize","Now treeH is %#x. MUONLoader is %#x",gime->TreeH(),gime);
263         return;
264       }
265       if (GetDebug()>2) {
266         cerr<<"AliMUONDigitizerv1::Exec inputFile is "<<inputFile<<" "<<endl;
267         cerr<<"AliMUONDigitizerv1::Exec treeH, fHits "<<treeH<<" "<<fHits<<endl;
268       }
269       if (treeH && fHits) {
270         branchHits = treeH->GetBranch("MUON");
271         if (branchHits) {
272           fHits->Clear();
273           branchHits->SetAddress(&fHits);
274         }
275         else
276           Error("Exec","branch MUON was not found");
277       }
278       if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::Exec branchHits = "<<branchHits<<endl;
279
280
281       fMask = fManager->GetMask(inputFile);
282       //
283       // Loop over tracks
284       Int_t itrack;
285       Int_t ntracks = (Int_t) treeH->GetEntries();
286       for (itrack = 0; itrack < ntracks; itrack++) {
287         if (GetDebug()>2) cerr<<"AliMUONDigitizerv1::Exec itrack = "<<itrack<<endl;
288         fHits->Clear();
289         treeH->GetEvent(itrack);
290         //
291         //  Loop over hits
292         Int_t ihit, ichamber;
293         AliMUONHit* mHit;
294         for(ihit = 0; ihit < fHits->GetEntriesFast(); ihit++) {
295           mHit = static_cast<AliMUONHit*>(fHits->At(ihit));
296           ichamber = mHit->Chamber()-1;  // chamber number
297           if (ichamber > AliMUONConstants::NCh()-1) {
298             cerr<<"AliMUONDigitizer: ERROR: "
299                 <<"fNch > AliMUONConstants::NCh()-1, fNch, NCh(): "
300                 <<ichamber<<", "<< AliMUONConstants::NCh()<<endl;
301             return;
302           }
303           chamber = &(pMUON->Chamber(ichamber));
304           //
305           //Dumping Hit content:
306           if (GetDebug()>2) {
307             cerr<<"AliMuonDigitizerv1::Exec ihit, ichamber, x, y, z, eloss " <<
308               ihit << " " << 
309               mHit->Chamber() << " " <<
310               mHit->X() << " " <<
311               mHit->Y() << " " <<
312               mHit->Z() << " " <<
313               mHit->Eloss() << " " << endl;
314           }
315           // 
316           // Inititializing Correlation
317           chamber->ChargeCorrelationInit();
318           if (ichamber < AliMUONConstants::NTrackingCh()) {
319             // Tracking Chamber
320             // Initialize hit position (cursor) in the segmentation model 
321             chamber->SigGenInit(mHit->X(), mHit->Y(), mHit->Z());
322           } else {
323             // Trigger Chamber 
324           }
325           MakeTransientDigit(itrack, ihit, mHit);
326         } // hit loop
327       } // track loop
328     } // end file loop
329     if (GetDebug()>2) cerr<<"AliMUONDigitizer::Exec End of hits, track and file loops"<<endl;
330
331     // Loop on cathodes
332     Int_t icat;
333     for(icat=0; icat<2; icat++) {
334       //
335       // Filling Digit List
336       Int_t tracks[kMAXTRACKS];
337       Int_t charges[kMAXTRACKS];
338       Int_t nentries = fTDList->GetEntriesFast();
339       Int_t digits[6];
340       for (Int_t nent = 0; nent < nentries; nent++) {
341         AliMUONTransientDigit *address = (AliMUONTransientDigit*)fTDList->At(nent);
342         if (address == 0) continue; 
343         Int_t ich = address->Chamber();
344         Int_t   q = address->Signal(); 
345         chamber = &(pMUON->Chamber(ich));
346         //
347         //  Digit Response (noise, threshold, saturation, ...)
348         AliMUONResponse * response = chamber->ResponseModel();
349         q = response->DigitResponse(q,address);
350         
351         if (!q) continue;
352         
353         digits[0] = address->PadX();
354         digits[1] = address->PadY();
355         digits[2] = address->Cathode()-1;
356         digits[3] = q;
357         digits[4] = address->Physics();
358         digits[5] = address->Hit();
359         
360         Int_t nptracks = address->GetNTracks();
361         
362         if (nptracks > kMAXTRACKS) {
363           if (GetDebug() >0) {
364             cerr<<"AliMUONDigitizer:Exec  nptracks > 10 "<<nptracks;
365             cerr<<"reset to max value "<<kMAXTRACKS<<endl;
366           }
367           nptracks = kMAXTRACKS;
368         }
369         if (nptracks > 2 && GetDebug() >2) {
370           cerr<<"AliMUONDigitizer::Exec  nptracks > 2 "<<nptracks<<endl;
371           //    printf("cat,ich,ix,iy,q %d %d %d %d %d \n",icat,ich,digits[0],digits[1],q);
372         }
373         for (Int_t tr = 0; tr < nptracks; tr++) {
374           tracks[tr]   = address->GetTrack(tr);
375           charges[tr]  = address->GetCharge(tr);
376         }      //end loop over list of tracks for one pad
377         // Sort list of tracks according to charge
378         if (nptracks > 1) {
379           SortTracks(tracks,charges,nptracks);
380         }
381         if (nptracks < kMAXTRACKS ) {
382           for (Int_t i = nptracks; i < kMAXTRACKS; i++) {
383             tracks[i]  = -1;
384             charges[i] = 0;
385           }
386         }
387         
388         // fill digits
389         if (GetDebug()>3) cerr<<"AliMUONDigitzerv1::Exex TransientDigit to Digit"<<endl;
390         if ( digits[2] == icat ) pMUON->AddDigits(ich,tracks,charges,digits);
391       }
392       // Filling list of digits per chamber for a given cathode.
393       treeD->Fill();
394       pMUON->ResetDigits();    
395     } // end loop cathode
396     fTDList->Delete();  
397     
398     for(Int_t ii = 0; ii < 2*AliMUONConstants::NCh(); ++ii) {
399       if (fHitMap[ii]) {
400         delete fHitMap[ii];
401         fHitMap[ii] = 0;
402       }
403     }
404     
405     if (GetDebug()>2) 
406       cerr<<"AliMUONDigitizer::Exec: writing the TreeD: "
407           <<treeD->GetName()<<endl;
408
409     gimeout->WriteDigits("OVERWRITE");
410     delete [] fHitMap;
411     delete fTDList;
412     
413     if (fHits)    fHits->Clear();
414 }
415
416 //------------------------------------------------------------------------
417 void AliMUONDigitizerv1::SortTracks(Int_t *tracks,Int_t *charges,Int_t ntr)
418 {
419   //
420   // Sort the list of tracks contributing to a given digit
421   // Only the 3 most significant tracks are acctually sorted
422   //
423   
424   //
425   //  Loop over signals, only 3 times
426   //
427   
428   Int_t qmax;
429   Int_t jmax;
430   Int_t idx[3] = {-2,-2,-2};
431   Int_t jch[3] = {-2,-2,-2};
432   Int_t jtr[3] = {-2,-2,-2};
433   Int_t i,j,imax;
434   
435   if (ntr<3) imax=ntr;
436   else imax=3;
437   for(i=0;i<imax;i++){
438     qmax=0;
439     jmax=0;
440     
441     for(j=0;j<ntr;j++){
442       
443       if((i == 1 && j == idx[i-1]) 
444          ||(i == 2 && (j == idx[i-1] || j == idx[i-2]))) continue;
445       
446       if(charges[j] > qmax) {
447         qmax = charges[j];
448         jmax=j;
449       }       
450     } 
451     
452     if(qmax > 0) {
453       idx[i]=jmax;
454       jch[i]=charges[jmax]; 
455       jtr[i]=tracks[jmax]; 
456     }
457     
458   } 
459   
460   for(i=0;i<3;i++){
461     if (jtr[i] == -2) {
462       charges[i]=0;
463       tracks[i]=0;
464     } else {
465       charges[i]=jch[i];
466       tracks[i]=jtr[i];
467     }
468   }
469 }