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