]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDSDigitizer.cxx
Check if the provided path is a rootfile or a directory by using the
[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 "TFile.h"
43 #include "TTask.h"
44 #include "TTree.h"
45 #include "TSystem.h"
46 #include "TROOT.h"
47 #include "TFolder.h"
48 #include <stdlib.h>
49 #include <iostream.h>
50 #include <fstream.h>
51
52 ClassImp(AliFMDSDigitizer)
53
54 //____________________________________________________________________________ 
55   AliFMDSDigitizer::AliFMDSDigitizer():TTask("AliFMDSDigitizer","") 
56 {
57   // ctor
58   fNevents = 0 ;     
59   fSDigits = 0 ;
60   fHits = 0 ;
61
62 }
63            
64 //____________________________________________________________________________ 
65   AliFMDSDigitizer::AliFMDSDigitizer(char* HeaderFile,char *SdigitsFile ):TTask("AliFMDSDigitizer","") 
66 {
67   fNevents = 0 ;     // Number of events to digitize, 0 means all evens in current file
68   // add Task to //root/Tasks folder
69   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
70   roottasks->Add(this) ; 
71 }
72
73 //____________________________________________________________________________ 
74   AliFMDSDigitizer::~AliFMDSDigitizer()
75 {
76   // dtor
77 }
78
79 //---------------------------------------------------------------------
80 void AliFMDSDigitizer::SetRingsSi1(Int_t ringsSi1)
81 {
82   fRingsSi1=ringsSi1;
83 }
84 void AliFMDSDigitizer::SetSectorsSi1(Int_t sectorsSi1)
85 {
86   fSectorsSi1=sectorsSi1;
87 }
88 void AliFMDSDigitizer::SetRingsSi2(Int_t ringsSi2)
89 {
90   fRingsSi2=ringsSi2;
91 }
92 void AliFMDSDigitizer::SetSectorsSi2(Int_t sectorsSi2)
93 {
94   fSectorsSi2=sectorsSi2;
95 }
96
97 //____________________________________________________________________________
98 void AliFMDSDigitizer::Exec(Option_t *option) { 
99
100
101
102
103   Int_t NumberOfRings[5]=
104   {fRingsSi1,fRingsSi2,fRingsSi1,fRingsSi2,fRingsSi1};
105   Int_t NumberOfSectors[5]=
106   {fSectorsSi1,fSectorsSi2,fSectorsSi1,fSectorsSi2,fSectorsSi1};
107
108   // Initialise Hit array
109   fHits = new TClonesArray ("AliFMDhit", 1000);
110   fSDigits = new TClonesArray ("AliFMDdigit", 1000);
111
112   AliFMD *FMD = (AliFMD *) gAlice->GetDetector ("FMD");
113
114   if (fNevents == 0)
115     fNevents = (Int_t) gAlice->TreeE ()->GetEntries ();
116
117   for (Int_t ievent = 0; ievent < fNevents; ievent++)
118     {
119       gAlice->GetEvent (ievent);
120       if (gAlice->TreeH () == 0)
121         return;
122       if (gAlice->TreeS () == 0)
123         gAlice->MakeTree ("S");
124
125
126
127       
128             //Make branches
129       char branchname[20];
130        sprintf (branchname, "%s", FMD->GetName ());
131       //Make branch for digits
132         FMD->MakeBranch ("S");
133     
134        //Now made SDigits from hits, for PHOS it is the same
135       Int_t volume, sector, ring, charge;
136       Float_t e;
137       Float_t de[10][50][150];
138       Int_t ivol, isec, iring;
139       Int_t hit, nbytes;
140       TParticle *particle;
141       AliFMDhit *fmdHit;
142       TClonesArray *FMDhits = FMD->Hits ();
143
144       // Event ------------------------- LOOP  
145
146       for (ivol = 0; ivol < 10; ivol++)
147         for (isec = 0; isec < 50; isec++)
148           for (iring = 0; iring < 150; iring++)
149             de[ivol][isec][iring] = 0;
150
151       if (FMD)
152         {
153           FMDhits = FMD->Hits ();
154           TTree *TH = gAlice->TreeH ();
155           Stat_t ntracks = TH->GetEntries ();
156           for (Int_t track = 0; track < ntracks; track++)
157             {
158               gAlice->ResetHits ();
159               nbytes += TH->GetEvent (track);
160               particle = gAlice->Particle (track);
161               Int_t nhits = FMDhits->GetEntriesFast ();
162
163               for (hit = 0; hit < nhits; hit++)
164                 {
165                   fmdHit = (AliFMDhit *) FMDhits->UncheckedAt (hit);
166
167                   volume = fmdHit->Volume ();
168                   sector = fmdHit->NumberOfSector ();
169                   ring = fmdHit->NumberOfRing ();
170                   e = fmdHit->Edep ();
171                   de[volume][sector][ring] += e;
172                 }               //hit loop
173             }                   //track loop
174         }                       //if FMD
175
176
177       Int_t digit[5];
178       Float_t I = 1.664 * 0.04 * 2.33 / 22400;  // = 0.69e-6;
179       for (ivol = 1; ivol < 6; ivol++)
180         { 
181           for (isec = 1; isec <= NumberOfSectors[ivol-1]; isec++)
182             { 
183               for (iring = 1; iring <= NumberOfRings[ivol-1]; iring++)
184                 {
185                       digit[0] = ivol;
186                       digit[1] = isec;
187                       digit[2] = iring;
188                       charge = Int_t (de[ivol][isec][iring] / I);
189
190                       digit[3] = charge;
191                       //dinamic diapason from MIP(0.155MeV) to 30MIP(4.65MeV)
192                       //1024 ADC channels 
193                       Float_t channelWidth = (22400 * 30) / 1024;
194
195                       digit[4] = Int_t (digit[3] / channelWidth);
196                       FMD->AddSDigit(digit);
197
198                 }               // iring loop
199             }                   //sector loop
200         }                       // volume loop
201       
202       gAlice->TreeS()->Reset();
203       gAlice->TreeS()->Fill();
204       gAlice->TreeS()->Write(0,TObject::kOverwrite) ;
205     }                           //event loop
206
207
208 }
209  
210 //__________________________________________________________________
211 void AliFMDSDigitizer::SetSDigitsFile(char * file ){
212   if(!fSDigitsFile.IsNull())
213     cout << "Changing SDigits file from " <<(char *)fSDigitsFile.Data() << " to " << file << endl ;
214   fSDigitsFile=file ;
215 }
216 //__________________________________________________________________
217 void AliFMDSDigitizer::Print(Option_t* option)const
218 {
219   cout << "------------------- "<< GetName() << " -------------" << endl ;
220   if(fSDigitsFile.IsNull())
221     cout << " Writing SDigitis to file galice.root "<< endl ;
222   else
223     cout << "    Writing SDigitis to file  " << (char*) fSDigitsFile.Data() << endl ;
224
225 }