]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDInput.cxx
Correct raw data reconstruction in case of trigger
[u/mrichter/AliRoot.git] / FMD / AliFMDInput.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, 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 /* $Id$ */
16 /** @file    AliFMDInput.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Mon Mar 27 12:42:40 2006
19     @brief   FMD utility classes for reading FMD data
20 */
21 //___________________________________________________________________
22 //
23 // The classes defined here, are utility classes for reading in data
24 // for the FMD.  They are  put in a seperate library to not polute the
25 // normal libraries.  The classes are intended to be used as base
26 // classes for customized class that do some sort of analysis on the
27 // various types of data produced by the FMD. 
28 //
29 // Latest changes by Christian Holm Christensen
30 //
31 #include "AliFMDInput.h"        // ALIFMDHIT_H
32 #include "AliLog.h"             // ALILOG_H
33 #include "AliLoader.h"          // ALILOADER_H
34 #include "AliRunLoader.h"       // ALIRUNLOADER_H
35 #include "AliRun.h"             // ALIRUN_H
36 #include "AliStack.h"           // ALISTACK_H
37 #include "AliRawReaderFile.h"   // ALIRAWREADERFILE_H
38 #include "AliFMD.h"             // ALIFMD_H
39 #include "AliFMDHit.h"          // ALIFMDHIT_H
40 #include "AliFMDDigit.h"        // ALIFMDDigit_H
41 #include "AliFMDRecPoint.h"     // ALIFMDRECPOINT_H
42 #include "AliFMDRawReader.h"    // ALIFMDRAWREADER_H
43 #include <AliESD.h>
44 #include <AliESDFMD.h>
45 #include <AliCDBManager.h>
46 #include <AliCDBEntry.h>
47 #include <AliAlignObjAngles.h>
48 #include <TTree.h>              // ROOT_TTree
49 #include <TChain.h>             // ROOT_TChain
50 #include <TParticle.h>          // ROOT_TParticle
51 #include <TString.h>            // ROOT_TString
52 #include <TDatabasePDG.h>       // ROOT_TDatabasePDG
53 #include <TMath.h>              // ROOT_TMath
54 #include <TGeoManager.h>        // ROOT_TGeoManager 
55 #include <TSystemDirectory.h>   // ROOT_TSystemDirectory
56 #include <Riostream.h>          // ROOT_Riostream
57
58 //____________________________________________________________________
59 ClassImp(AliFMDInput)
60 #if 0
61   ; // This is here to keep Emacs for indenting the next line
62 #endif
63
64
65 //____________________________________________________________________
66 AliFMDInput::AliFMDInput()
67   : fGAliceFile(""), 
68     fLoader(0),
69     fRun(0), 
70     fStack(0),
71     fFMDLoader(0), 
72     fReader(0),
73     fFMD(0),
74     fMainESD(0),
75     fESD(0),
76     fTreeE(0),
77     fTreeH(0),
78     fTreeD(0),
79     fTreeS(0),
80     fTreeR(0), 
81     fChainE(0),
82     fArrayE(0),
83     fArrayH(0),
84     fArrayD(0),
85     fArrayS(0), 
86     fArrayR(0), 
87     fArrayA(0), 
88     fTreeMask(0), 
89     fIsInit(kFALSE)
90 {
91   // Constructor of an FMD input object.  Specify what data to read in
92   // using the AddLoad member function.   Sub-classes should at a
93   // minimum overload the member function Event.   A full job can be
94   // executed using the member function Run. 
95 }
96
97   
98
99 //____________________________________________________________________
100 AliFMDInput::AliFMDInput(const char* gAliceFile)
101   : fGAliceFile(gAliceFile),
102     fLoader(0),
103     fRun(0), 
104     fStack(0),
105     fFMDLoader(0), 
106     fReader(0),
107     fFMD(0),
108     fMainESD(0),
109     fESD(0),
110     fTreeE(0),
111     fTreeH(0),
112     fTreeD(0),
113     fTreeS(0),
114     fTreeR(0), 
115     fChainE(0),
116     fArrayE(0),
117     fArrayH(0),
118     fArrayD(0),
119     fArrayS(0), 
120     fArrayR(0), 
121     fArrayA(0), 
122     fTreeMask(0), 
123     fIsInit(kFALSE)
124 {
125   // Constructor of an FMD input object.  Specify what data to read in
126   // using the AddLoad member function.   Sub-classes should at a
127   // minimum overload the member function Event.   A full job can be
128   // executed using the member function Run. 
129 }
130
131 //____________________________________________________________________
132 Int_t
133 AliFMDInput::NEvents() const 
134 {
135   // Get number of events
136   if (fTreeE) return fTreeE->GetEntries();
137   return -1;
138 }
139
140 //____________________________________________________________________
141 Bool_t
142 AliFMDInput::Init()
143 {
144   // Initialize the object.  Get the needed loaders, and such. 
145
146   // Check if we have been initialized
147   if (fIsInit) { 
148     AliWarning("Already initialized");
149     return fIsInit;
150   }
151   if (fGAliceFile.IsNull()) fGAliceFile = "galice.root";
152   // Get the loader
153   fLoader = AliRunLoader::Open(fGAliceFile.Data(), "Alice", "read");
154   if (!fLoader) {
155     AliError(Form("Coulnd't read the file %s", fGAliceFile.Data()));
156     return kFALSE;
157   }
158   
159   // Get the run 
160   if  (fLoader->LoadgAlice()) return kFALSE;
161   fRun = fLoader->GetAliRun();
162   
163   // Get the FMD 
164   fFMD = static_cast<AliFMD*>(fRun->GetDetector("FMD"));
165   if (!fFMD) {
166     AliError("Failed to get detector FMD from loader");
167     return kFALSE;
168   }
169   
170   // Get the FMD loader
171   fFMDLoader = fLoader->GetLoader("FMDLoader");
172   if (!fFMDLoader) {
173     AliError("Failed to get detector FMD loader from loader");
174     return kFALSE;
175   }
176   if (fLoader->LoadHeader()) { 
177     AliError("Failed to get event header information from loader");
178     return kFALSE;
179   }
180   fTreeE = fLoader->TreeE();
181
182   // Optionally, get the ESD files
183   if (TESTBIT(fTreeMask, kESD)) {
184     fChainE = new TChain("esdTree");
185     TSystemDirectory dir(".",".");
186     TList*           files = dir.GetListOfFiles();
187     TSystemFile*     file = 0;
188     if (!files) {
189       AliError("No files");
190       return kFALSE;
191     }
192     files->Sort();
193     TIter            next(files);
194     while ((file = static_cast<TSystemFile*>(next()))) {
195       TString fname(file->GetName());
196       if (fname.Contains("AliESDs")) fChainE->AddFile(fname.Data());
197     }
198     fChainE->SetBranchAddress("ESD", &fMainESD);
199   }
200     
201   if (TESTBIT(fTreeMask, kRaw)) {
202     AliInfo("Getting FMD raw data digits");
203     fArrayA = new TClonesArray("AliFMDDigit");
204     fReader = new AliRawReaderFile(-1);
205   }
206   
207   // Optionally, get the geometry 
208   if (TESTBIT(fTreeMask, kGeometry)) {
209     TString fname(fRun->GetGeometryFileName());
210     if (fname.IsNull()) {
211       Warning("Init", "No file name for the geometry from AliRun");
212       fname = gSystem->DirName(fGAliceFile);
213       fname.Append("/geometry.root");
214     }
215     fGeoManager = TGeoManager::Import(fname.Data());
216     if (!fGeoManager) {
217       Fatal("Init", "No geometry manager found");
218       return kFALSE;
219     }
220     AliCDBManager* cdb   = AliCDBManager::Instance();
221     AliCDBEntry*   align = cdb->Get("FMD/Align/Data");
222     if (align) {
223       AliInfo("Got alignment data from CDB");
224       TClonesArray* array = dynamic_cast<TClonesArray*>(align->GetObject());
225       if (!array) {
226         AliWarning("Invalid align data from CDB");
227       }
228       else {
229         Int_t nAlign = array->GetEntries();
230         for (Int_t i = 0; i < nAlign; i++) {
231           AliAlignObjAngles* a = static_cast<AliAlignObjAngles*>(array->At(i));
232           if (!a->ApplyToGeometry()) {
233             AliWarning(Form("Failed to apply alignment to %s", 
234                             a->GetVolPath()));
235           }
236         }
237       }
238     }
239   }
240
241   
242   fIsInit = kTRUE;
243   return fIsInit;
244 }
245
246 //____________________________________________________________________
247 Bool_t
248 AliFMDInput::Begin(Int_t event)
249 {
250   // Called at the begining of each event.  Per default, it gets the
251   // data trees and gets pointers to the output arrays.   Users can
252   // overload this, but should call this member function in the
253   // overloaded member function of the derived class. 
254
255   // Check if we have been initialized
256   if (!fIsInit) { 
257     AliError("Not initialized");
258     return fIsInit;
259   }
260   // Get the event 
261   if (fLoader->GetEvent(event)) return kFALSE;
262   AliInfo(Form("Now in event %d/%d", event, NEvents()));
263
264   // Possibly load global kinematics information 
265   if (TESTBIT(fTreeMask, kKinematics)) {
266     AliInfo("Getting kinematics");
267     if (fLoader->LoadKinematics()) return kFALSE;
268     fStack = fLoader->Stack();
269   }
270   // Possibly load FMD Hit information 
271   if (TESTBIT(fTreeMask, kHits)) {
272     AliInfo("Getting FMD hits");
273     if (fFMDLoader->LoadHits()) return kFALSE;
274     fTreeH = fFMDLoader->TreeH();
275     if (!fArrayH) fArrayH = fFMD->Hits(); 
276   }
277   // Possibly load FMD Digit information 
278   if (TESTBIT(fTreeMask, kDigits)) {
279     AliInfo("Getting FMD digits");
280     if (fFMDLoader->LoadDigits()) return kFALSE;
281     fTreeD = fFMDLoader->TreeD();
282     if (!fArrayD) fArrayD = fFMD->Digits();
283   }
284   // Possibly load FMD Sdigit information 
285   if (TESTBIT(fTreeMask, kSDigits)) {
286     AliInfo("Getting FMD summable digits");
287     if (fFMDLoader->LoadSDigits()) return kFALSE;
288     fTreeS = fFMDLoader->TreeS();
289     if (!fArrayS) fArrayS = fFMD->SDigits();
290   }
291   // Possibly load FMD RecPoints information 
292   if (TESTBIT(fTreeMask, kRecPoints)) {
293     AliInfo("Getting FMD reconstructed points");
294     if (fFMDLoader->LoadRecPoints()) return kFALSE;
295     fTreeR = fFMDLoader->TreeR();
296     if (!fArrayR) fArrayR = new TClonesArray("AliFMDRecPoint");
297     fTreeR->SetBranchAddress("FMD",  &fArrayR);
298   }  // Possibly load FMD ESD information 
299   if (TESTBIT(fTreeMask, kESD)) {
300     AliInfo("Getting FMD event summary data");
301     Int_t read = fChainE->GetEntry(event);
302     if (read <= 0) return kFALSE;
303     fESD = fMainESD->GetFMDData();
304     if (!fESD) return kFALSE;
305   }
306   // Possibly load FMD Digit information 
307   if (TESTBIT(fTreeMask, kRaw)) {
308     AliInfo("Getting FMD raw data digits");
309     if (!fReader->NextEvent()) return kFALSE;
310     AliFMDRawReader r(fReader, 0);
311     fArrayA->Clear();
312     r.ReadAdcs(fArrayA);
313   }
314   
315   return kTRUE;
316 }
317
318
319 //____________________________________________________________________
320 Bool_t 
321 AliFMDInput::Event()
322 {
323   // Process one event.  The default implementation one or more of 
324   //
325   //   -  ProcessHits       if the hits are loaded. 
326   //   -  ProcessDigits     if the digits are loaded. 
327   //   -  ProcessSDigits    if the sumbable digits are loaded. 
328   //   -  ProcessRecPoints  if the reconstructed points are loaded. 
329   //   -  ProcessESD        if the event summary data is loaded
330   // 
331   if (TESTBIT(fTreeMask, kHits)) 
332     if (!ProcessHits()) return kFALSE; 
333   if (TESTBIT(fTreeMask, kDigits)) 
334     if (!ProcessDigits()) return kFALSE;
335   if (TESTBIT(fTreeMask, kSDigits)) 
336     if (!ProcessSDigits()) return kFALSE;
337   if (TESTBIT(fTreeMask, kRaw)) 
338     if (!ProcessRawDigits()) return kFALSE;
339   if (TESTBIT(fTreeMask, kRecPoints)) 
340     if (!ProcessRecPoints()) return kFALSE;
341   if (TESTBIT(fTreeMask, kESD))
342     if (!ProcessESD(fESD)) return kFALSE;
343   
344   return kTRUE;
345 }
346
347 //____________________________________________________________________
348 Bool_t 
349 AliFMDInput::ProcessHits()
350 {
351   if (!fTreeH) {
352     AliError("No hit tree defined");
353     return kFALSE;
354   }
355   Int_t nTracks = fTreeH->GetEntries();
356   for (Int_t i = 0; i < nTracks; i++) {
357     Int_t hitRead  = fTreeH->GetEntry(i);
358     if (hitRead <= 0) continue;
359     if (!fArrayH) {
360       AliError("No hit array defined");
361       return kFALSE;
362     }
363     Int_t nHit = fArrayH->GetEntries();
364     if (nHit <= 0) continue;
365     for (Int_t j = 0; j < nHit; j++) {
366       AliFMDHit* hit = static_cast<AliFMDHit*>(fArrayH->At(j));
367       if (!hit) continue;
368       TParticle* track = 0;
369       if (TESTBIT(fTreeMask, kKinematics) && fStack) {
370         Int_t trackno = hit->Track();
371         track = fStack->Particle(trackno);
372       }
373       if (!ProcessHit(hit, track)) return kFALSE;
374     }    
375   }
376   return kTRUE;
377 }
378
379 //____________________________________________________________________
380 Bool_t 
381 AliFMDInput::ProcessDigits()
382 {
383   // Read the digit tree, and pass each digit to the member function
384   // ProcessDigit.
385   Int_t nEv = fTreeD->GetEntries();
386   for (Int_t i = 0; i < nEv; i++) {
387     Int_t digitRead  = fTreeD->GetEntry(i);
388     if (digitRead <= 0) continue;
389     Int_t nDigit = fArrayD->GetEntries();
390     if (nDigit <= 0) continue;
391     for (Int_t j = 0; j < nDigit; j++) {
392       AliFMDDigit* digit = static_cast<AliFMDDigit*>(fArrayD->At(j));
393       if (!digit) continue;
394       if (!ProcessDigit(digit)) return kFALSE;
395     }    
396   }
397   return kTRUE;
398 }
399
400 //____________________________________________________________________
401 Bool_t 
402 AliFMDInput::ProcessSDigits()
403 {
404   // Read the summable digit tree, and pass each sumable digit to the
405   // member function ProcessSdigit.
406   Int_t nEv = fTreeD->GetEntries();
407   for (Int_t i = 0; i < nEv; i++) {
408     Int_t sdigitRead  = fTreeS->GetEntry(i);
409     if (sdigitRead <= 0) continue;
410     Int_t nSdigit = fArrayS->GetEntries();
411     if (nSdigit <= 0) continue;
412     for (Int_t j = 0; j < nSdigit; j++) {
413       AliFMDSDigit* sdigit = static_cast<AliFMDSDigit*>(fArrayS->At(j));
414       if (!sdigit) continue;
415       if (!ProcessSDigit(sdigit)) return kFALSE;
416     }    
417   }
418   return kTRUE;
419 }
420
421 //____________________________________________________________________
422 Bool_t 
423 AliFMDInput::ProcessRawDigits()
424 {
425   // Read the digit tree, and pass each digit to the member function
426   // ProcessDigit.
427   Int_t nDigit = fArrayA->GetEntries();
428   if (nDigit <= 0) return kTRUE;
429   for (Int_t j = 0; j < nDigit; j++) {
430     AliFMDDigit* digit = static_cast<AliFMDDigit*>(fArrayA->At(j));
431     if (!digit) continue;
432     if (!ProcessRawDigit(digit)) return kFALSE;
433   }    
434   return kTRUE;
435 }
436
437 //____________________________________________________________________
438 Bool_t 
439 AliFMDInput::ProcessRecPoints()
440 {
441   // Read the reconstrcted points tree, and pass each reconstruction
442   // object (AliFMDRecPoint) to either ProcessRecPoint.
443   Int_t nEv = fTreeR->GetEntries();
444   for (Int_t i = 0; i < nEv; i++) {
445     Int_t recRead  = fTreeR->GetEntry(i);
446     if (recRead <= 0) continue;
447     Int_t nRecPoint = fArrayR->GetEntries();
448     for (Int_t j = 0; j < nRecPoint; j++) {
449       AliFMDRecPoint* recPoint = static_cast<AliFMDRecPoint*>(fArrayR->At(j));
450       if (!recPoint) continue;
451       if (!ProcessRecPoint(recPoint)) return kFALSE;
452     }    
453   }
454   return kTRUE;
455 }
456
457 //____________________________________________________________________
458 Bool_t
459 AliFMDInput::End()
460 {
461   // Called at the end of each event.  Per default, it unloads the
462   // data trees and resets the pointers to the output arrays.   Users
463   // can overload this, but should call this member function in the
464   // overloaded member function of the derived class. 
465
466   // Check if we have been initialized
467   if (!fIsInit) { 
468     AliError("Not initialized");
469     return fIsInit;
470   }
471   // Possibly unload global kinematics information 
472   if (TESTBIT(fTreeMask, kKinematics)) {
473     fLoader->UnloadKinematics();
474     // fTreeK = 0;
475     fStack = 0;
476   }
477   // Possibly unload FMD Hit information 
478   if (TESTBIT(fTreeMask, kHits)) {
479     fFMDLoader->UnloadHits();
480     fTreeH = 0;
481   }
482   // Possibly unload FMD Digit information 
483   if (TESTBIT(fTreeMask, kDigits)) {
484     fFMDLoader->UnloadDigits();
485     fTreeD = 0;
486   }
487   // Possibly unload FMD Sdigit information 
488   if (TESTBIT(fTreeMask, kSDigits)) {
489     fFMDLoader->UnloadSDigits();
490     fTreeS = 0;
491   }
492   // Possibly unload FMD RecPoints information 
493   if (TESTBIT(fTreeMask, kRecPoints)) {
494     fFMDLoader->UnloadRecPoints();
495     fTreeR = 0;
496   }
497   AliInfo("Now out event");
498   return kTRUE;
499 }
500
501 //____________________________________________________________________
502 Bool_t
503 AliFMDInput::Run()
504 {
505   // Run over all events and files references in galice.root 
506
507   Bool_t retval;
508   if (!(retval = Init())) return retval;
509
510   Int_t nEvents = NEvents();
511   for (Int_t event = 0; event < nEvents; event++) {
512     if (!(retval = Begin(event))) break;
513     if (!(retval = Event())) break;
514     if (!(retval = End())) break;
515   }
516   if (!retval) return retval;
517   retval = Finish();
518   return retval;
519 }
520
521 //====================================================================
522 ClassImp(AliFMDInputHits)
523 #if 0
524   ;
525 #endif
526
527
528 //====================================================================
529 ClassImp(AliFMDInputDigits)
530 #if 0
531   ;
532 #endif
533
534 //====================================================================
535 ClassImp(AliFMDInputSDigits)
536 #if 0
537   ;
538 #endif
539
540 //====================================================================
541 ClassImp(AliFMDInputRecPoints)
542 #if 0
543   ;
544 #endif
545
546 //____________________________________________________________________
547 //
548 // EOF
549 //