]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDSDigitizer.cxx
173f1c8bd42233e082764c0eabc7bf0cf354db69
[u/mrichter/AliRoot.git] / FMD / AliFMDSDigitizer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 // This is a TTask that constructs SDigits out of Hits
18 // A Summable Digits is the sum of all hits in a cell
19 // A threshold is applied 
20 //
21 //-- Author: Alla Maevskaia(INR)
22 //////////////////////////////////////////////////////////////////////////////
23
24 // --- ROOT system ---
25 #include "TTask.h"
26 #include "TTree.h"
27 #include "TSystem.h"
28 #include "TFile.h"
29 // --- Standard library ---
30
31 // --- AliRoot header files ---
32
33 #include "AliFMDdigit.h"
34 #include "AliFMDhit.h"
35 #include "AliFMD.h"
36 #include "AliFMDv1.h"
37 #include "AliFMDSDigitizer.h"
38 #include "AliRun.h"
39 #include "AliDetector.h"
40 #include "AliMC.h"
41
42 #include "TROOT.h"
43 #include "TFolder.h"
44 #include <stdlib.h>
45 #include <iostream.h>
46 #include <fstream.h>
47
48 ClassImp(AliFMDSDigitizer)
49
50            
51 //____________________________________________________________________________ 
52   AliFMDSDigitizer::AliFMDSDigitizer():TTask("AliFMDSDigitizer","") 
53 {
54   fNevents = 0 ;     // Number of events to digitize, 0 means all evens in current file
55   // add Task to //root/Tasks folder
56   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
57   roottasks->Add(this) ; 
58 }
59 //____________________________________________________________________________ 
60 AliFMDSDigitizer::AliFMDSDigitizer(char* HeaderFile, char *SDigitsFile):TTask("AliFMDSDigitizer","")
61 {
62   // ctor
63   fNevents = 0 ;    // Number of events to digitize, 0 means all events in current file
64   fSDigitsFile = SDigitsFile ;
65   fHeadersFile = HeaderFile ;
66   //add Task to //root/Tasks folder
67   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
68   roottasks->Add(this) ; 
69     
70 }
71
72 //____________________________________________________________________________ 
73   AliFMDSDigitizer::~AliFMDSDigitizer()
74 {
75   // dtor
76 }
77
78
79 //____________________________________________________________________________
80 void AliFMDSDigitizer::Exec(Option_t *option) { 
81   //Collects all hits in the same active volume into digit
82   TClonesArray * sdigits = new TClonesArray("AliFMDdigit",1000) ;
83   TFile * file = 0;
84
85   AliFMD * FMD = (AliFMD *) gAlice->GetDetector("FMD") ;
86   
87   if(fNevents == 0) 
88     fNevents = (Int_t) gAlice->TreeE()->GetEntries() ; 
89
90   cout<<"AliFMDSDigitizer-> Nevents"<<fNevents<<endl;
91   for(Int_t ievent = 0; ievent < fNevents; ievent++){
92     gAlice->GetEvent(ievent) ;
93     if(gAlice->TreeH()==0) return ;
94
95     if(gAlice->TreeS() == 0)            
96       gAlice->MakeTree("S") ;
97     
98     TClonesArray * FMDhits = FMD->Hits() ;
99
100     
101
102     Int_t nSdigits = 0 ;
103     
104     //Make branches
105     char branchname[20];
106     sprintf(branchname,"%s",FMD->GetName());  
107     
108     Int_t bufferSize = 16000 ;
109     char * file =0;
110     if(!fSDigitsFile.IsNull())
111       file = (char*) fSDigitsFile.Data() ; //ievent ;
112     else
113       if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
114         file = new char[30] ;
115         sprintf(file,"FMD.SDigits.root") ;
116         cout<<"CONFIG_SPLIT_FILE "<<file<<endl; 
117       }
118       else{
119         file = 0 ;
120         cout<<" FILE =0 "<<endl;}
121         cout<<"After CONFIG_SPLIT_FILE "<<file<<endl; 
122     
123     gAlice->MakeBranchInTree(gAlice->TreeS(),branchname,&sdigits,bufferSize,file);  
124     /*    
125     Int_t splitlevel = 0 ;
126     sprintf(branchname,"AliFMDSDigitizer");   
127     AliFMDSDigitizer * sd = this ;
128     gAlice->MakeBranchInTree(gAlice->TreeS(),branchname,"AliFMDSDigitizer",&sd, bufferSize, splitlevel,file); 
129     */
130
131     //Now made SDigits from hits, for PHOS it is the same
132   Int_t volume,sector,ring,charge;
133   Float_t e;
134   Float_t de[10][20][150];
135   Int_t ivol,isec,iring;
136   Int_t hit,nbytes;
137   TParticle *particle;
138   AliFMDhit  *fmdHit;
139
140  
141  // Event ------------------------- LOOP  
142  
143   for (ivol=1; ivol<=6; ivol++)
144     for(isec=0; isec<16; isec++)
145       for(iring=0; iring<128; iring++)
146         de[ivol][isec][iring]=0;
147       
148   if (FMD){
149     FMDhits   = FMD->Hits();
150     TTree *TH = gAlice->TreeH();
151     Stat_t ntracks    = TH->GetEntries();
152  
153    for (Int_t track=0; track<ntracks;track++) {
154      gAlice->ResetHits();
155      nbytes += TH->GetEvent(track);
156      particle=gAlice->Particle(track);
157      Int_t nhits =FMDhits->GetEntriesFast();
158      
159      for (hit=0;hit<nhits;hit++) {
160        fmdHit   = (AliFMDhit*)FMDhits->UncheckedAt(hit);
161        
162        volume = fmdHit->Volume();
163        sector = fmdHit->NumberOfSector();
164        ring   = fmdHit->NumberOfRing();
165        e = fmdHit->Edep();
166        de[volume][sector][ring]=de[volume][sector][ring]+e;
167
168      } //hit loop
169    } //track loop
170   }//if FMD
171   Int_t digit[5];
172   Float_t I=1.664*0.04*2.33/22400; // = 0.69e-6;
173   for(ivol=1; ivol<=6; ivol++){
174     for(isec=0; isec<16; isec++){
175       for(iring=0; iring<128; iring++){
176         if(de[ivol][isec][iring]>0.){
177           digit[0]=ivol;
178           digit[1]=isec;
179           digit[2]=iring;
180           charge=Int_t (de[ivol][isec][iring]/I);
181           digit[3]=charge;
182
183           //dinamic diapason from MIP(0.155MeV) to 30MIP(4.65MeV)
184           //1024 ADC channels 
185           Float_t channelWidth=(22400*30)/1024;
186           digit[4]=Int_t(digit[3]/channelWidth);
187
188           new((*sdigits)[nSdigits++]) AliFMDdigit(digit) ;
189           
190         } //de >threshold
191       }// iring loop
192     }//sector loop
193   } // volume loop
194
195
196   
197   
198     gAlice->TreeS()->Fill() ;
199     TFile *f1 = new TFile(file,"RECREATE");
200     f1->cd();
201     gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
202   }
203
204   delete sdigits ;
205   if(file)
206     file->Close() ;
207
208 }
209 //__________________________________________________________________
210 void AliFMDSDigitizer::SetSDigitsFile(char * file ){
211   if(!fSDigitsFile.IsNull())
212     cout << "Changing SDigits file from " <<(char *)fSDigitsFile.Data() << " to " << file << endl ;
213   fSDigitsFile=file ;
214 }
215 //__________________________________________________________________
216 void AliFMDSDigitizer::Print(Option_t* option)const{
217   cout << "------------------- "<< GetName() << " -------------" << endl ;
218   if(fSDigitsFile.IsNull())
219     cout << " Writing SDigitis to file galice.root "<< endl ;
220   else
221     cout << "    Writing SDigitis to file  " << (char*) fSDigitsFile.Data() << endl ;
222
223 }