]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSBeamTestDigitizer.cxx
Stand Alone tracker updated to use displaced primary vertices in the bending plane...
[u/mrichter/AliRoot.git] / ITS / AliITSBeamTestDigitizer.cxx
CommitLineData
38300302 1////////////////////////////////////////////////////
2// Class to manage the //
3// ITS beam test conversion from rawdata //
4// to digits. It executes the digitization for //
5// SPD, SDD and SSD. //
7d62fb64 6// //
7// //
38300302 8// Origin: E. Crescio crescio@to.infn.it //
9// J. Conrad Jan.Conrad@cern.ch //
10////////////////////////////////////////////////////
11#include "AliHeader.h"
12#include "AliRun.h"
13#include "AliRunLoader.h"
14#include "AliITSEventHeader.h"
15#include "AliITSLoader.h"
38300302 16#include "AliITSBeamTestDigSDD.h"
17#include "AliITSBeamTestDigSPD.h"
18#include "AliITSBeamTestDigSSD.h"
19#include "AliITSBeamTestDigitizer.h"
20#include "AliRawReaderDate.h"
8ace09b6 21#include "AliRawReaderRoot.h"
38300302 22
5ba31760 23const TString AliITSBeamTestDigitizer::fgkDefaultDigitsFileName="ITS.Digits.root";
38300302 24
25ClassImp(AliITSBeamTestDigitizer)
26
27
28//_____________________________________________________________
7537d03c 29AliITSBeamTestDigitizer::AliITSBeamTestDigitizer():TTask(),
30fEvIn(0),
31fEvFin(0),
32fRunNumber(-1),
33fDATEEvType(7),
34fFlagHeader(kTRUE),
35fFlagInit(kFALSE),
36fOptDate(kFALSE),
37fDigitsFileName(0),
38fRawdataFileName(0),
39fPeriod(kNov04),
40fRunLoader(0),
41fLoader(0),
42fHeader(0){
38300302 43 //
44 // Default constructor
45 //
38300302 46 SetFlagInit();
8ace09b6 47 SetOptDate();
38300302 48}
49
50//_____________________________________________________________
7537d03c 51AliITSBeamTestDigitizer::AliITSBeamTestDigitizer(const Text_t* name, const Text_t* title, Char_t* opt,const char* filename):TTask(name,title),
52fEvIn(0),
53fEvFin(0),
54fRunNumber(-1),
55fDATEEvType(7),
56fFlagHeader(kTRUE),
57fFlagInit(kFALSE),
58fOptDate(kFALSE),
59fDigitsFileName(0),
60fRawdataFileName(0),
61fPeriod(kNov04),
62fRunLoader(0),
63fLoader(0),
64fHeader(0)
38300302 65{
66 //
67 // Standard constructor
68 //
8ace09b6 69 SetOptDate();
5ba31760 70
71 TString choice(opt);
72 Bool_t aug04 = choice.Contains("Aug04");
73 Bool_t nov04 = choice.Contains("Nov04");
74 if(aug04) fPeriod=kAug04;
75 if(nov04) fPeriod=kNov04;
7d62fb64 76 Init(filename);
38300302 77 }
78
79//_____________________________________________________________
7537d03c 80AliITSBeamTestDigitizer::AliITSBeamTestDigitizer(const Text_t* name, const Text_t* title, Int_t run, Char_t* opt,const char* filename):TTask(name,title),
81fEvIn(0),
82fEvFin(0),
83fRunNumber(run),
84fDATEEvType(7),
85fFlagHeader(kTRUE),
86fFlagInit(kFALSE),
87fOptDate(kFALSE),
88fDigitsFileName(0),
89fRawdataFileName(0),
90fPeriod(kNov04),
91fRunLoader(0),
92fLoader(0),
93fHeader(0){
38300302 94 //
95 // Constructor
96 //
8ace09b6 97 SetOptDate();
5ba31760 98 TString choice(opt);
99 Bool_t aug04 = choice.Contains("Aug04");
100 Bool_t nov04 = choice.Contains("Nov04");
101 if(aug04) fPeriod=kAug04;
102 if(nov04) fPeriod=kNov04;
103
7d62fb64 104 Init(filename);
38300302 105 }
106
107//___________________________________________________________
7d62fb64 108void AliITSBeamTestDigitizer::Init(const char* filename){
38300302 109
110 //
111 //Initialization of run loader and its loader
112 //creation of galice.root
113 //
5ba31760 114
38300302 115
7d62fb64 116 fRunLoader = AliRunLoader::Open(filename,AliConfig::GetDefaultEventFolderName(),"update");
38300302 117 if (fRunLoader == 0x0)
118 {
119 Error("AliITSBeamTestDigitizer","Can not load the session",filename);
120 return;
121 }
122 fRunLoader->LoadgAlice();
123 gAlice = fRunLoader->GetAliRun();
124
125 if(!gAlice) {
126 Error("AliITSBeamTestDigitizer","gAlice not found on file. Aborting.");
127 return;
128 }
7d62fb64 129 fRunLoader->MakeTree("E");
38300302 130 fLoader = (AliITSLoader*)fRunLoader->GetLoader("ITSLoader");
38300302 131
132 fDigitsFileName=fgkDefaultDigitsFileName;
7d62fb64 133 this->Add(new AliITSBeamTestDigSPD("DigSPD","SPD Digitization"));
134 this->Add(new AliITSBeamTestDigSDD("DigSDD","SDD Digitization"));
135 this->Add(new AliITSBeamTestDigSSD("DigSSD","SSD Digitization"));
38300302 136
7d62fb64 137 SetFlagInit(kTRUE);
38300302 138}
139
7d62fb64 140
38300302 141//______________________________________________________________________
7537d03c 142AliITSBeamTestDigitizer::AliITSBeamTestDigitizer(const AliITSBeamTestDigitizer &bt):TTask(bt),
143fEvIn(bt.fEvIn),
144fEvFin(bt.fEvFin),
145fRunNumber(bt.fRunNumber),
146fDATEEvType(bt.fDATEEvType),
147fFlagHeader(bt.fFlagHeader),
148fFlagInit(bt.fFlagInit),
149fOptDate(bt.fOptDate),
150fDigitsFileName(bt.fDigitsFileName),
151fRawdataFileName(bt.fRawdataFileName),
152fPeriod(bt.fPeriod),
153fRunLoader(bt.fRunLoader),
154fLoader(bt.fLoader),
155fHeader(bt.fHeader){
38300302 156 // Copy constructor.
38300302 157}
158//______________________________________________________________________
5ba31760 159AliITSBeamTestDigitizer& AliITSBeamTestDigitizer::operator=(const AliITSBeamTestDigitizer &source){
7537d03c 160 // Assignment operator.
161 this->~AliITSBeamTestDigitizer();
162 new(this) AliITSBeamTestDigitizer(source);
163 return *this;
38300302 164}
165
166
167//______________________________________________________________
168AliITSBeamTestDigitizer::~AliITSBeamTestDigitizer(){
169
170 //Destructor
7d62fb64 171 // if(fBt) delete fBt;
38300302 172 if(fLoader) delete fLoader;
5ba31760 173 if(fHeader) delete fHeader;
38300302 174}
175
176
177//_____________________________________________________________
178void AliITSBeamTestDigitizer::SetNumberOfEventsPerFile(Int_t nev)
179{
180 //Sets number of events per file
181
182 if(fRunLoader) fRunLoader->SetNumberOfEventsPerFile(nev);
183 else Warning("SetNumberOfEventsPerFile","fRunLoader is 0");
184}
185
186
187//____________________________________________________
188void AliITSBeamTestDigitizer::ExecDigitization(){
189
190 // Execution of digitisation for SPD,SDD and SSD
191
192 if(!GetFlagInit()){
193 Warning("ExecDigitization()","Run Init() please..");
194 return;
195 }
196 fLoader->SetDigitsFileName(fDigitsFileName);
197 fLoader->LoadDigits("recreate");
198
8ace09b6 199 AliRawReader* rd;
200
201 if(GetOptDate()) rd = new AliRawReaderDate(fRawdataFileName,fEvIn);
202 else rd = new AliRawReaderRoot(fRawdataFileName,fEvIn);
203
38300302 204 AliHeader* header = fRunLoader->GetHeader();
38300302 205 Int_t iev=fEvIn-1;
206
207
208 AliITSBeamTestDigSDD* digSDD = (AliITSBeamTestDigSDD*)fTasks->FindObject("DigSDD");
209 AliITSBeamTestDigSPD* digSPD = (AliITSBeamTestDigSPD*)fTasks->FindObject("DigSPD");
210 AliITSBeamTestDigSSD* digSSD = (AliITSBeamTestDigSSD*)fTasks->FindObject("DigSSD");
211
212
213 do{
214 iev++;
215 if(fEvFin!=0){
216 if(iev>fEvFin) break;
217 }
218 AliITSEventHeader* itsh = new AliITSEventHeader("ITSHeader");
219 fRunLoader->SetEventNumber(iev);
220
8ace09b6 221 rd->RequireHeader(fFlagHeader);
222 rd->SelectEvents(fDATEEvType);
38300302 223
8ace09b6 224 digSDD->SetRawReader(rd);
225 digSPD->SetRawReader(rd);
226 digSSD->SetRawReader(rd);
38300302 227
228 if(fLoader->TreeD() == 0x0) fLoader->MakeTree("D");
229
230 TTree* treeD = (TTree*)fLoader->TreeD();
231
232 // Make branches outside the dig-classes
233
234 TClonesArray* digitsSPD = new TClonesArray("AliITSdigitSPD",1000);
5ba31760 235 treeD->Branch("ITSDigitsSPD",&digitsSPD);
38300302 236
237 TClonesArray* digitsSDD = new TClonesArray("AliITSdigitSDD",1000);
5ba31760 238 treeD->Branch("ITSDigitsSDD",&digitsSDD);
38300302 239
240 TClonesArray* digitsSSD = new TClonesArray("AliITSdigitSSD",1000);
5ba31760 241 treeD->Branch("ITSDigitsSSD",&digitsSSD);
38300302 242
243
244 digSSD->SetTree(treeD);
245 digSDD->SetTree(treeD);
246 digSPD->SetTree(treeD);
5ba31760 247
7d62fb64 248 AliITSgeom* geom = fLoader->GetITSgeom();
249
250 digSSD->SetITSgeom(geom);
251 digSDD->SetITSgeom(geom);
252 digSPD->SetITSgeom(geom);
5ba31760 253
38300302 254 digSSD->SetITSEventHeader(itsh);
255 digSDD->SetITSEventHeader(itsh);
256 digSPD->SetITSEventHeader(itsh);
257
258 digSDD->SetBtPeriod(GetBeamTestPeriod());
5ba31760 259 if(GetBeamTestPeriod()==1)digSDD->SetThreshold(16);
260 else digSDD->SetThreshold(0);
38300302 261 ExecuteTask(0);
262
263 header->SetEventNrInRun(iev);
264 header->SetEvent(iev);
265 header->SetRun(fRunNumber);
266 fRunLoader->GetHeader()->AddDetectorEventHeader(itsh);
267 fRunLoader->TreeE()->Fill();
268 header->Reset(fRunNumber,iev);
269
270 delete digitsSPD;
271 delete digitsSDD;
272 delete digitsSSD;
273
8ace09b6 274 }while(rd->NextEvent());
38300302 275
8ace09b6 276
38300302 277 fRunLoader->WriteHeader("OVERWRITE");
278 fRunLoader->WriteRunLoader("OVERWRITE");
5ba31760 279
8ace09b6 280 delete rd;
38300302 281 fLoader->UnloadDigits();
282 fLoader->UnloadRawClusters();
283 fRunLoader->UnloadHeader();
5ba31760 284
38300302 285
286}
287
288
289
290//_______________________________________________
291void AliITSBeamTestDigitizer:: SetActive(const TString& subdet,Bool_t value){
292
293 //Sets active sub-tasks (detectors)
294
295 Bool_t sdd = subdet.Contains("SDD");
296 Bool_t spd = subdet.Contains("SPD");
297 Bool_t ssd = subdet.Contains("SSD");
298
299 if(sdd){
300 AliITSBeamTestDigSDD* digSDD = (AliITSBeamTestDigSDD*)fTasks->FindObject("DigSDD");
301 digSDD->SetActive(value);
302 }
303
304 if(spd){
305 AliITSBeamTestDigSPD* digSPD = (AliITSBeamTestDigSPD*)fTasks->FindObject("DigSPD");
306 digSPD->SetActive(value);
307
308 }
309
310 if(ssd){
311 AliITSBeamTestDigSSD* digSSD = (AliITSBeamTestDigSSD*)fTasks->FindObject("DigSSD");
312 digSSD->SetActive(value);
313
314 }
315
316
317
318}
319