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