]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FIT/AliFITDigitizer.cxx
new geometryy from WT
[u/mrichter/AliRoot.git] / FIT / AliFITDigitizer.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  *    Produde digits from hits
19  *       digits is TObject and includes
20  *      We are writing array if C & A  TDC
21  *      C & A  ADC (will need for slow simulation)
22  *      TOF first particle C & A
23  *      mean time and time difference (vertex position)
24  *
25  *      Alla.Maevskaya@cern.ch 
26  ****************************************************************/
27
28
29 #include <TArrayI.h>
30 #include <TFile.h>
31 #include <TGraph.h>
32 #include <TH1F.h>
33 #include <TMath.h>
34 #include <TRandom.h>
35 #include <TTree.h> 
36
37 #include "AliLog.h"
38 #include "AliFITDigitizer.h"
39 #include "AliFIT.h"
40 #include "AliFITHits.h"
41 #include "AliFITDigit.h"
42 #include "AliDigitizationInput.h"
43 #include "AliRun.h"
44 #include <AliLoader.h>
45 #include <AliRunLoader.h>
46 #include <stdlib.h>
47
48 ClassImp(AliFITDigitizer)
49
50 //___________________________________________
51   AliFITDigitizer::AliFITDigitizer()  :AliDigitizer(),
52                                      fFIT(0),
53                                      fHits(0),
54                                      fDigits(0),
55                                      fNdigits(0)
56 {
57 // Default ctor - don't use it
58
59 }
60
61 //___________________________________________
62 AliFITDigitizer::AliFITDigitizer(AliDigitizationInput* digInput) 
63   :AliDigitizer(digInput),
64    fFIT(0),
65    fHits(0),
66    fDigits(0),
67    fNdigits(0)
68 {
69 // ctor which should be used
70  
71 }
72
73
74 //------------------------------------------------------------------------
75 AliFITDigitizer::~AliFITDigitizer()
76 {
77 // Destructor
78
79   AliDebug(1,"FIT");
80
81  }
82
83 //------------------------------------------------------------------------
84 Bool_t AliFITDigitizer::Init()
85 {
86 // Initialization
87   AliDebug(1," Init");
88  return kTRUE;
89
90 }
91  
92 //---------------------------------------------------------------------
93 void AliFITDigitizer::Digitize(Option_t* /*option*/)
94 {
95
96   /*
97     Produde digits from hits
98     digits is TObject and includes
99     We are writing array if C & A for each channel CFD, LED, QT0 and QT1
100     C & A  ADC (will need for slow simulation)
101   */
102   
103   
104   
105   //output loader 
106   AliDebug(1,"start...");
107   //input loader
108   //
109   // From hits to digits
110   //
111   
112   AliRunLoader *outRL = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
113   AliLoader * outFitLoader = outRL->GetLoader("FITLoader");
114   
115   fFIT  = static_cast<AliFIT*>(gAlice->GetDetector("FIT"));
116   if (!fFIT) {
117     AliError("Can not get FIT from gAlice");
118     return;
119   }  
120   fFIT->ResetDigits();
121   
122   DigitizeHits();
123   
124   //load digits    
125   outFitLoader->LoadDigits("UPDATE");
126   TTree *treeD  = outFitLoader->TreeD();
127   if (treeD == 0x0) {
128     outFitLoader->MakeTree("D");
129     treeD = outFitLoader->TreeD();
130   }
131   treeD->Reset();
132   fFIT  = (AliFIT*)outRL ->GetAliRun()->GetDetector("FIT");
133   // Make a branch in the tree 
134   fFIT->MakeBranch("D");
135   treeD->Fill();
136   
137   outFitLoader->WriteDigits("OVERWRITE");
138   fFIT->ResetDigits();
139   outFitLoader->UnloadDigits();
140 }
141 //____________________________________________________________________________
142 void AliFITDigitizer::DigitizeHits()
143 {
144
145    Int_t hit, nhits;
146   Float_t time[160], besttime[160];
147   Int_t countE[160];
148   Int_t timeCFD, timeLED, timeQT1, timeQT0;
149
150   Int_t threshold = 0; //photoelectrons
151   Float_t channelWidth = 24.4 ; 
152   Int_t pmt, mcp, volume, qt; 
153   //eqailized distance from IP
154   Float_t zdetC = 82;
155   Float_t zdetA = 375.;
156   Float_t c = 0.0299792458; // cm/ps
157   Float_t eqdistance = (zdetA - zdetC) /c;
158   
159   AliFITHits  *startHit;
160   TBranch *brHits=0;
161  
162   Int_t nFiles=fDigInput->GetNinputs();
163   for (Int_t inputFile=0; inputFile<nFiles;  inputFile++) {
164     if (inputFile < nFiles-1) {
165       AliWarning(Form("ignoring input stream %d", inputFile));
166       continue;
167      }
168
169   for (Int_t i0=0; i0<160; i0++)
170       {
171         time[i0]=besttime[i0]=999999; countE[i0]=0;
172       }
173   AliRunLoader * inRL = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inputFile));
174   AliLoader * fitLoader = inRL->GetLoader("FITLoader");
175   if (!inRL->GetAliRun()) inRL->LoadgAlice();
176   Int_t numpmt;
177   //read Hits 
178   fitLoader->LoadHits("READ");//probably it is necessary to load them before
179   fHits = fFIT->Hits ();
180   TTree *th = fitLoader->TreeH();
181   brHits = th->GetBranch("FIT");
182   if (brHits) {
183     fFIT->SetHitsAddressBranch(brHits);
184   }else{
185     AliWarning("Branch FIT hit not found for this event");
186     continue;      
187   } 
188   Int_t ntracks    = (Int_t) th->GetEntries();
189   if (ntracks<=0) return;
190   // Start loop on tracks in the hits containers
191   for (Int_t track=0; track<ntracks;track++) {
192     brHits->GetEntry(track);
193     nhits = fHits->GetEntriesFast();
194     for (hit=0;hit<nhits;hit++) 
195       {
196         startHit   = (AliFITHits*) fHits->UncheckedAt(hit);
197         if (!startHit) {
198           AliError("The unchecked hit doesn't exist");
199           break;
200         }
201         pmt = startHit->Pmt();
202         mcp = startHit->MCP();
203         volume = startHit->Volume();
204         numpmt= volume*80 +4*(mcp) + (pmt-1);
205         besttime[numpmt] = startHit->Time();
206         if(besttime[numpmt]<time[numpmt]) time[numpmt]=besttime[numpmt];
207         countE[numpmt]++;
208       } //hits loop
209   } //track loop
210   
211   for (Int_t ipmt=0; ipmt<160; ipmt++)
212     {
213       if (countE[ipmt]>threshold && time[ipmt]<50000 && time[ipmt]>0 ) {
214         //fill ADC
215         // QTC procedure:
216         // 1MIP ->200phe  ;
217           qt= countE[ipmt] /*/ph2Mip*/;  // 50mv/Mip amp in mV 
218           //  fill TDC
219           if (ipmt<80) time[ipmt] = time[ipmt] + eqdistance;
220           timeCFD = Int_t (gRandom->Gaus(time[ipmt], 50)/channelWidth ); 
221           timeLED =  Int_t (time[ipmt]/channelWidth );
222           timeQT0 = 0;
223           timeQT1 = Int_t (qt ) ;
224           AliDebug(1,Form("Digits:::::  numpmt %i  time CFD %i QTC  %i :: counts %i \n ",  ipmt, timeCFD, timeQT1,  countE[ipmt]) );
225           fFIT-> AddDigit(ipmt,   timeCFD, timeLED, timeQT0,  timeQT1, 0);
226       } //hitted PMTs
227     } //pmt loop
228   fitLoader->UnloadHits();
229   
230   }
231 }
232
233 //____________________________________________________________________________
234 void AliFITDigitizer::AddSDigit(Int_t npmt,  
235                                 Int_t timeCFD, Int_t timeLED, Int_t timeQT0, 
236                                 Int_t timeQT1) 
237  { 
238  
239 // Adds Digit 
240  
241   TClonesArray &ldigits = *fDigits;  
242          
243   new(ldigits[fNdigits++]) AliFITDigit(npmt,  
244                                        timeCFD, timeLED, timeQT0, timeQT1);
245          
246 }
247 //____________________________________________________________________________
248 void AliFITDigitizer::ResetDigits()
249 {
250
251 // Clears Digits
252
253   fNdigits = 0;
254   if (fDigits) fDigits->Clear();
255 }