]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSSDigitizer.cxx
new volume names
[u/mrichter/AliRoot.git] / PHOS / AliPHOSSDigitizer.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 /* $Id$ */
17
18 //_________________________________________________________________________
19 // This is a TTask that constructs SDigits out of Hits
20 // A Summable Digits is the sum of all hits in a cell
21 // A threshold is applied 
22 //
23 //*-- Author :  Dmitri Peressounko (SUBATECH & KI) 
24 //////////////////////////////////////////////////////////////////////////////
25
26 // --- ROOT system ---
27 #include "TTask.h"
28 #include "TTree.h"
29 #include "TSystem.h"
30 // --- Standard library ---
31
32 // --- AliRoot header files ---
33
34 #include "AliPHOSDigit.h"
35 #include "AliPHOSHit.h"
36 #include "AliPHOSv1.h"
37 #include "AliPHOSSDigitizer.h"
38
39 #include "TROOT.h"
40 #include "TFolder.h"
41
42 ClassImp(AliPHOSSDigitizer)
43
44            
45 //____________________________________________________________________________ 
46   AliPHOSSDigitizer::AliPHOSSDigitizer():TTask("AliPHOSSDigitizer","") 
47 {
48   // ctor
49   fA = 0;
50   fB = 10000000. ;
51   fPrimThreshold = 0.01 ;
52   fNevents = 0 ;     // Number of events to digitize, 0 means all evens in current file
53   // add Task to //root/Tasks folder
54   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
55   roottasks->Add(this) ; 
56
57 }
58 //____________________________________________________________________________ 
59 AliPHOSSDigitizer::AliPHOSSDigitizer(char* HeaderFile, char *SDigitsFile):TTask("AliPHOSSDigitizer","")
60 {
61   // ctor
62   fA = 0;
63   fB = 10000000.;
64   fPrimThreshold = 0.01 ;
65   fNevents = 0 ;         // Number of events to digitize, 0 means all events in current file
66   fSDigitsFile = SDigitsFile ;
67   fHeadersFile = HeaderFile ;
68   //add Task to //root/Tasks folder
69   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
70   roottasks->Add(this) ; 
71     
72 }
73
74 //____________________________________________________________________________ 
75   AliPHOSSDigitizer::~AliPHOSSDigitizer()
76 {
77   // dtor
78 }
79
80
81 //____________________________________________________________________________
82 void AliPHOSSDigitizer::Exec(Option_t *option) { 
83   //Collects all hits in the same active volume into digit
84
85   TFile * file = 0;
86
87   // if(gAlice->TreeE()==0)        //If gAlice not yet read/constructed
88   if(fHeadersFile.IsNull())
89     file = new TFile("galice.root", "update") ;
90   else
91     file = new TFile(fHeadersFile.Data(),"update") ;
92
93   gAlice = (AliRun *) file->Get("gAlice") ;
94   
95   
96   
97   TClonesArray * sdigits = new TClonesArray("AliPHOSDigit",1000) ;
98
99   AliPHOS * PHOS = (AliPHOS *) gAlice->GetDetector("PHOS") ;
100   
101     
102   if(fNevents == 0) 
103     fNevents = (Int_t) gAlice->TreeE()->GetEntries() ; 
104   
105   Int_t ievent ;
106   for(ievent = 0; ievent < fNevents; ievent++){
107     gAlice->GetEvent(ievent) ;
108     gAlice->SetEvent(ievent) ;
109     
110     if(gAlice->TreeH()==0){
111       cout << "AliPHOSSDigitizer: There is no Hit Tree" << endl;
112       return ;
113     }
114     
115     if(gAlice->TreeS() == 0)
116       gAlice->MakeTree("S") ;
117     
118     TClonesArray * hits = PHOS->Hits() ;
119     
120     sdigits->Clear();
121     Int_t nSdigits = 0 ;
122     
123     //Make branches
124     char branchname[20];
125     sprintf(branchname,"%s",PHOS->GetName());  
126     
127     Int_t bufferSize = 16000 ;
128     char * file =0;
129     if(!fSDigitsFile.IsNull())
130       file = (char*) fSDigitsFile.Data() ; //ievent ;
131     else
132       if(gSystem->Getenv("CONFIG_SPLIT_FILE")){ //generating file name
133         file = new char[30] ;
134         //      sprintf(file,"PHOS.SDigits%d.root",ievent) ;
135         sprintf(file,"PHOS.SDigits.root") ;
136       }
137       else
138         file = 0 ;
139     
140     gAlice->MakeBranchInTree(gAlice->TreeS(),branchname,&sdigits,bufferSize,file);  
141
142
143     Int_t splitlevel = 0 ;
144     sprintf(branchname,"AliPHOSSDigitizer");   
145     AliPHOSSDigitizer * sd = this ;
146     gAlice->MakeBranchInTree(gAlice->TreeS(),branchname,"AliPHOSSDigitizer",&sd, bufferSize, splitlevel,file); 
147     
148
149     //Now made SDigits from hits, for PHOS it is the same
150
151     Int_t itrack ;
152     for (itrack=0; itrack<gAlice->GetNtrack(); itrack++){
153       
154       //=========== Get the Hits Tree for the Primary track itrack
155       gAlice->ResetHits();    
156       gAlice->TreeH()->GetEvent(itrack);
157       
158       Int_t i;
159       for ( i = 0 ; i < hits->GetEntries() ; i++ ) {
160         AliPHOSHit * hit = (AliPHOSHit*)hits->At(i) ;
161         AliPHOSDigit * newdigit ;
162
163         // Assign primary number only if contribution is significant
164         if( hit->GetEnergy() > fPrimThreshold)
165           newdigit = new AliPHOSDigit( hit->GetPrimary(), hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
166         else
167           newdigit = new AliPHOSDigit( -1               , hit->GetId(), Digitize( hit->GetEnergy() ) ) ;
168         
169         new((*sdigits)[nSdigits]) AliPHOSDigit(* newdigit) ;
170         nSdigits++ ;  
171         
172         delete newdigit ;    
173       } 
174       
175     } // loop over tracks
176     
177     sdigits->Sort() ;
178     
179     nSdigits = sdigits->GetEntries() ;
180     sdigits->Expand(nSdigits) ;
181     
182     Int_t i ;
183     for (i = 0 ; i < nSdigits ; i++) { 
184       AliPHOSDigit * digit = (AliPHOSDigit *) sdigits->At(i) ; 
185       digit->SetIndexInList(i) ;     
186     }
187     
188     gAlice->TreeS()->Fill() ;
189     gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
190   }
191
192   delete sdigits ;
193   if(file)
194     file->Close() ;
195
196 }
197 //__________________________________________________________________
198 void AliPHOSSDigitizer::SetSDigitsFile(char * file ){
199   if(!fSDigitsFile.IsNull())
200     cout << "Changing SDigits file from " <<(char *)fSDigitsFile.Data() << " to " << file << endl ;
201   fSDigitsFile=file ;
202 }
203 //__________________________________________________________________
204 void AliPHOSSDigitizer::Print(Option_t* option)const{
205   cout << "------------------- "<< GetName() << " -------------" << endl ;
206   if(fSDigitsFile.IsNull())
207     cout << " Writing SDigitis to file galice.root "<< endl ;
208   else
209     cout << "    Writing SDigitis to file  " << (char*) fSDigitsFile.Data() << endl ;
210   cout << "   with digitization parameters A = " << fA << endl ;
211   cout << "                                B = " << fB << endl ;
212   cout << "Threshold for Primary assignment  = " << fPrimThreshold << endl ; 
213   cout << "---------------------------------------------------"<<endl ;
214
215 }
216 //__________________________________________________________________
217 Bool_t AliPHOSSDigitizer::operator==( AliPHOSSDigitizer const &sd )const{
218   if( (fA==sd.fA)&&(fB==sd.fB)&&(fPrimThreshold==sd.fPrimThreshold))
219     return kTRUE ;
220   else
221     return kFALSE ;
222 }