]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMFTDigitizer.cxx
analysis code updated
[u/mrichter/AliRoot.git] / MFT / AliMFTDigitizer.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 //====================================================================================================================================================
17 //
18 //      Digitizer class for the ALICE Muon Forward Tracker
19 //
20 //      Contact author: antonio.uras@cern.ch
21 //
22 //====================================================================================================================================================
23
24 #include "AliRun.h"
25 #include "AliRunLoader.h"
26 #include "AliDigitizationInput.h"
27 #include "AliLoader.h"
28 #include "AliLog.h"
29 #include "AliMFTDigitizer.h"
30 #include "AliMFTDigit.h"
31 #include "AliMFT.h"
32 #include "AliMFTSegmentation.h"
33 #include "TObjArray.h"
34 #include "TClonesArray.h"
35 #include "AliDigitizer.h"
36
37 ClassImp(AliMFTDigitizer)
38
39 //====================================================================================================================================================
40     
41 AliMFTDigitizer::AliMFTDigitizer():
42   AliDigitizer(),
43   fNPlanes(0),
44   fSegmentation(0)
45 {
46
47   // default constructor 
48
49 }
50
51 //====================================================================================================================================================
52
53 AliMFTDigitizer::AliMFTDigitizer(AliDigitizationInput *digInp):
54   AliDigitizer(digInp),
55   fNPlanes(0),
56   fSegmentation(0)
57 {
58
59 }
60
61 //====================================================================================================================================================
62
63 void AliMFTDigitizer::Digitize(Option_t*) {
64
65   // This method is responsible for merging sdigits to a list of digits
66
67   AliDebug(1, "************************************************************************");
68   AliDebug(1, "************************ AliMFTDigitizer::Digitize *********************");
69   AliDebug(1, "************************************************************************");
70   
71   if (!fSegmentation) {
72     fSegmentation = new AliMFTSegmentation("AliMFTGeometry.root");
73     fNPlanes = fSegmentation -> GetNPlanes();
74   }
75
76   AliDebug(1, Form("nPlanes = %d",fNPlanes));
77
78   AliDebug(1,Form("Start with %i input(s) for event %i", fDigInput->GetNinputs(), fDigInput->GetOutputEventNr()));
79     
80   AliRunLoader *pInRunLoader=0;
81   AliLoader    *pInMFTLoader=0;
82   
83   TClonesArray sDigits[fNMaxPlanes];
84   for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) sDigits[iPlane].SetClass("AliMFTDigit");  // tmp storage for sdigits sum up from all input files
85   
86   // filling the arrays of sdigits...
87
88   for (Int_t iFile=0; iFile<fDigInput->GetNinputs(); iFile++) {
89     
90     pInRunLoader = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iFile));
91     pInMFTLoader = pInRunLoader->GetLoader("MFTLoader"); 
92     if (!pInMFTLoader) {
93       AliDebug(1,"no MFT lodader, checking in the other input \n"); 
94       continue;
95     }
96     
97     if (!pInRunLoader->GetAliRun()) pInRunLoader->LoadgAlice();
98     AliMFT *pInMFT = (AliMFT*) pInRunLoader->GetAliRun()->GetDetector("MFT"); 
99     
100     AliDebug(1, "pInMFTLoader->LoadSDigits()...");
101     pInMFTLoader->LoadSDigits();
102     AliDebug(1, "... done!");
103     AliDebug(1, "    pInMFTLoader->TreeS()->GetEntry(0);");
104     pInMFTLoader->TreeS()->GetEntry(0);
105     AliDebug(1, "... done!");
106     
107     for (Int_t iPlane=0; iPlane<pInMFT->GetSDigitsList()->GetEntries(); iPlane++) {      
108       for(Int_t iSDig=0; iSDig<((TClonesArray*)pInMFT->GetSDigitsList()->At(iPlane))->GetEntries(); iSDig++) {
109         AliDebug(2, Form("Reading digit %03d of plane %02d (A)", iSDig, iPlane));
110         AliMFTDigit *pSDig = (AliMFTDigit*) ((TClonesArray*)pInMFT->GetSDigitsList()->At(iPlane))->At(iSDig);
111         AliDebug(2, Form("Reading digit %03d of plane %02d (B)", iSDig, iPlane));
112         pSDig->AddOffset2TrackID(fDigInput->GetMask(iFile));   // -> To be introduced for merging (since all inputs count tracks independently from 0)
113         AliDebug(2, Form("Reading digit %03d of plane %02d (C)", iSDig, iPlane));
114         new ((sDigits[iPlane])[sDigits[iPlane].GetEntries()]) AliMFTDigit(*pSDig);  
115         AliDebug(2, Form("Reading digit %03d of plane %02d (D)", iSDig, iPlane));
116       }
117     }
118     
119     pInMFTLoader->UnloadSDigits();   
120     pInMFT->ResetSDigits();
121
122   }
123   
124   AliRunLoader *pOutRunLoader = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());  // open output stream (only 1 possible)
125   AliLoader    *pOutMFTLoader = pOutRunLoader->GetLoader("MFTLoader");                        
126   AliRun       *pAliRun       = pOutRunLoader->GetAliRun();
127   AliMFT       *pOutMFT       = (AliMFT*) pAliRun->GetDetector("MFT");      
128   pOutMFTLoader->MakeTree("D");   
129   pOutMFT->MakeBranch("D");
130   pOutMFT->SetTreeAddress();
131   
132   SDigits2Digits(sDigits, pOutMFT->GetDigitsList());   // here the sdigits are merged into digits
133   
134   pOutMFTLoader->TreeD()->Fill();              // fill the output tree with the list of digits
135   pOutMFTLoader->WriteDigits("OVERWRITE");     // serialize them to file
136   
137   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) sDigits[iPlane].Clear();   // remove all tmp sdigits
138   pOutMFTLoader->UnloadDigits();   
139   pOutMFT->ResetDigits(); 
140
141 }
142
143 //====================================================================================================================================================
144
145 void AliMFTDigitizer::SDigits2Digits(TClonesArray *pSDigitList, TObjArray *pDigitList) {   
146
147   TClonesArray *myDigitList[fNMaxPlanes] = {0};
148  
149   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) { 
150     myDigitList[iPlane] = (TClonesArray*)(*pDigitList)[iPlane];
151     if (myDigitList[iPlane]->GetEntries()!=0) AliErrorClass("Some of digits lists is not empty");   // in principle those lists should be empty 
152   }
153    
154   AliDebug(1,"starting loop over planes");
155    
156   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
157      
158     AliMFTDigit *newDig=NULL;
159     AliMFTDigit *oldDig=NULL;
160     pSDigitList[iPlane].Sort();
161
162     AliDebug(1,"starting loop over sdigits to create digits");
163     AliDebug(1, Form("MFT plane #%02d has %d SDigits", iPlane, Int_t(pSDigitList[iPlane].GetEntries())));
164
165     for (Int_t iSDig=0; iSDig<pSDigitList[iPlane].GetEntries(); iSDig++) {
166
167       newDig = (AliMFTDigit*) (pSDigitList[iPlane].At(iSDig));
168       Bool_t digitExists = kFALSE;
169       Int_t nDigits = myDigitList[iPlane]->GetEntries();
170       
171       for (Int_t iDig=0; iDig<nDigits; iDig++) {
172         oldDig = (AliMFTDigit*) (myDigitList[iPlane]->At(iDig));
173         if (newDig->GetDetElemID()==oldDig->GetDetElemID() &&
174             newDig->GetPixelX()==oldDig->GetPixelX() &&
175             newDig->GetPixelY()==oldDig->GetPixelY() &&
176             newDig->GetPixelZ()==oldDig->GetPixelZ()) {
177           digitExists = kTRUE;
178           MergeDigits(oldDig, newDig);
179           break;
180         }
181       }
182
183       if (!digitExists) new ((*myDigitList[iPlane])[myDigitList[iPlane]->GetEntries()]) AliMFTDigit(*newDig);
184
185     }
186     
187     AliDebug(1, Form("MFT plane #%02d has %d Digits", iPlane, Int_t(myDigitList[iPlane]->GetEntries())));
188     
189 //     for (Int_t iDigit=0; iDigit<myDigitList[iPlane]->GetEntries(); iDigit++) {
190 //       AliDebug(1, Form("Digit %03d of MFT plane %02d has pixel coordinates (%05d, %05d)", 
191 //                     iDigit, iPlane, ((AliMFTDigit*) myDigitList[iPlane]->At(iDigit))->GetPixelX(), ((AliMFTDigit*) myDigitList[iPlane]->At(iDigit))->GetPixelY()) );
192 //     }
193
194     AliDebug(1, "ending loop over sdigits to create digits");
195
196   }
197
198   AliDebug(1,"ending loop over layers");  
199
200 }
201
202 //====================================================================================================================================================
203
204 void AliMFTDigitizer::MergeDigits(AliMFTDigit *mainDig, AliMFTDigit *digToSum) {
205   
206   mainDig -> SetEloss(mainDig->GetEloss() + digToSum->GetEloss());
207   
208   Bool_t trackExists = kFALSE;
209   for (Int_t iTrack=0; iTrack<mainDig->GetNMCTracks(); iTrack++) {
210     if (digToSum->GetMCLabel(0) == mainDig->GetMCLabel(iTrack)) {
211       trackExists = kTRUE;
212       break;
213     }
214   }
215   
216   if (!trackExists) mainDig->AddMCLabel(digToSum->GetMCLabel(0));
217   
218 }
219
220 //====================================================================================================================================================
221