1 /**************************************************************************
2 * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
21 //___________________________________________________________________
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.
29 // Latest changes by Christian Holm Christensen
31 #include "AliFMDInput.h" // ALIFMDHIT_H
32 #include "AliFMDDebug.h" // ALIFMDDEBUG_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 "AliRawReaderRoot.h" // ALIRAWREADERROOT_H
39 #include "AliRawReaderDate.h" // ALIRAWREADERDATE_H
40 #include "AliFMD.h" // ALIFMD_H
41 #include "AliFMDHit.h" // ALIFMDHIT_H
42 #include "AliFMDDigit.h" // ALIFMDDigit_H
43 #include "AliFMDSDigit.h" // ALIFMDDigit_H
44 #include "AliFMDRecPoint.h" // ALIFMDRECPOINT_H
45 #include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
47 #include <AliESDFMD.h>
48 #include <AliESDEvent.h>
49 #include <AliCDBManager.h>
50 #include <AliCDBEntry.h>
51 #include <AliAlignObjParams.h>
52 #include <TTree.h> // ROOT_TTree
53 #include <TChain.h> // ROOT_TChain
54 #include <TParticle.h> // ROOT_TParticle
55 #include <TString.h> // ROOT_TString
56 #include <TDatabasePDG.h> // ROOT_TDatabasePDG
57 #include <TMath.h> // ROOT_TMath
58 #include <TGeoManager.h> // ROOT_TGeoManager
59 #include <TSystemDirectory.h> // ROOT_TSystemDirectory
60 #include <Riostream.h> // ROOT_Riostream
61 #include <TFile.h> // ROOT_TFile
62 #include <TStreamerInfo.h>
65 //____________________________________________________________________
68 ; // This is here to keep Emacs for indenting the next line
72 //____________________________________________________________________
73 AliFMDInput::AliFMDInput()
74 : TNamed("AliFMDInput", "Input handler for various FMD data"),
105 // Constructor of an FMD input object. Specify what data to read in
106 // using the AddLoad member function. Sub-classes should at a
107 // minimum overload the member function Event. A full job can be
108 // executed using the member function Run.
113 //____________________________________________________________________
114 AliFMDInput::AliFMDInput(const char* gAliceFile)
115 : TNamed("AliFMDInput", "Input handler for various FMD data"),
116 fGAliceFile(gAliceFile),
146 // Constructor of an FMD input object. Specify what data to read in
147 // using the AddLoad member function. Sub-classes should at a
148 // minimum overload the member function Event. A full job can be
149 // executed using the member function Run.
152 //____________________________________________________________________
154 AliFMDInput::NEvents() const
156 // Get number of events
157 if (fTreeE) return fTreeE->GetEntries();
161 //____________________________________________________________________
165 // Initialize the object. Get the needed loaders, and such.
167 // Check if we have been initialized
169 AliWarning("Already initialized");
172 if (fGAliceFile.IsNull()) fGAliceFile = "galice.root";
174 fLoader = AliRunLoader::Open(fGAliceFile.Data(), "Alice", "read");
176 AliError(Form("Coulnd't read the file %s", fGAliceFile.Data()));
181 if (fLoader->LoadgAlice()) return kFALSE;
182 fRun = fLoader->GetAliRun();
185 fFMD = static_cast<AliFMD*>(fRun->GetDetector("FMD"));
187 AliError("Failed to get detector FMD from loader");
191 // Get the FMD loader
192 fFMDLoader = fLoader->GetLoader("FMDLoader");
194 AliError("Failed to get detector FMD loader from loader");
197 if (fLoader->LoadHeader()) {
198 AliError("Failed to get event header information from loader");
201 fTreeE = fLoader->TreeE();
203 // Optionally, get the ESD files
204 if (TESTBIT(fTreeMask, kESD)) {
205 fChainE = new TChain("esdTree");
206 TSystemDirectory dir(".",".");
207 TList* files = dir.GetListOfFiles();
208 TSystemFile* file = 0;
210 AliError("No files");
215 while ((file = static_cast<TSystemFile*>(next()))) {
216 TString fname(file->GetName());
217 if (fname.Contains("AliESDs")) fChainE->AddFile(fname.Data());
219 fESDEvent = new AliESDEvent();
220 fESDEvent->ReadFromTree(fChainE);
221 // fChainE->SetBranchAddress("ESD", &fMainESD);
225 if (TESTBIT(fTreeMask, kRaw)) {
226 AliInfo("Getting FMD raw data digits");
227 fArrayA = new TClonesArray("AliFMDDigit");
228 if (!fRawFile.IsNull() && fRawFile.EndsWith(".root"))
229 fReader = new AliRawReaderRoot(fRawFile.Data());
230 else if (!fRawFile.IsNull() && fRawFile.EndsWith(".raw"))
231 fReader = new AliRawReaderDate(fRawFile.Data());
233 fReader = new AliRawReaderFile(-1);
237 // Optionally, get the geometry
238 if (TESTBIT(fTreeMask, kGeometry)) {
239 TString fname(fRun->GetGeometryFileName());
240 if (fname.IsNull()) {
241 Warning("Init", "No file name for the geometry from AliRun");
242 fname = gSystem->DirName(fGAliceFile);
243 fname.Append("/geometry.root");
245 fGeoManager = TGeoManager::Import(fname.Data());
247 Fatal("Init", "No geometry manager found");
250 AliCDBManager* cdb = AliCDBManager::Instance();
251 AliCDBEntry* align = cdb->Get("FMD/Align/Data");
253 AliInfo("Got alignment data from CDB");
254 TClonesArray* array = dynamic_cast<TClonesArray*>(align->GetObject());
256 AliWarning("Invalid align data from CDB");
259 Int_t nAlign = array->GetEntries();
260 for (Int_t i = 0; i < nAlign; i++) {
261 AliAlignObjParams* a = static_cast<AliAlignObjParams*>(array->At(i));
262 if (!a->ApplyToGeometry()) {
263 AliWarning(Form("Failed to apply alignment to %s",
276 //____________________________________________________________________
278 AliFMDInput::Begin(Int_t event)
280 // Called at the begining of each event. Per default, it gets the
281 // data trees and gets pointers to the output arrays. Users can
282 // overload this, but should call this member function in the
283 // overloaded member function of the derived class.
285 // Check if we have been initialized
287 AliError("Not initialized");
292 if (fLoader->GetEvent(event)) return kFALSE;
293 AliInfo(Form("Now in event %8d/%8d", event, NEvents()));
295 // Possibly load global kinematics information
296 if (TESTBIT(fTreeMask, kKinematics) || TESTBIT(fTreeMask, kTracks)) {
297 // AliInfo("Getting kinematics");
298 if (fLoader->LoadKinematics("READ")) return kFALSE;
299 fStack = fLoader->Stack();
302 // Possibly load FMD Hit information
303 if (TESTBIT(fTreeMask, kHits) || TESTBIT(fTreeMask, kTracks)) {
304 // AliInfo("Getting FMD hits");
305 if (!fFMDLoader || fFMDLoader->LoadHits("READ")) return kFALSE;
306 fTreeH = fFMDLoader->TreeH();
307 if (!fArrayH) fArrayH = fFMD->Hits();
310 // Possibly load heaedr information
311 if (TESTBIT(fTreeMask, kHeader)) {
312 // AliInfo("Getting FMD hits");
313 if (!fLoader /* || fLoader->LoadHeader()*/) return kFALSE;
314 fHeader = fLoader->GetHeader();
317 // Possibly load FMD Digit information
318 if (TESTBIT(fTreeMask, kDigits)) {
319 // AliInfo("Getting FMD digits");
320 if (!fFMDLoader || fFMDLoader->LoadDigits("READ")) return kFALSE;
321 fTreeD = fFMDLoader->TreeD();
323 if (!fArrayD) fArrayD = fFMD->Digits();
327 AliWarning(Form("Failed to load FMD Digits"));
331 // Possibly load FMD Sdigit information
332 if (TESTBIT(fTreeMask, kSDigits)) {
333 // AliInfo("Getting FMD summable digits");
334 if (!fFMDLoader || fFMDLoader->LoadSDigits("READ")) return kFALSE;
335 fTreeS = fFMDLoader->TreeS();
336 if (!fArrayS) fArrayS = fFMD->SDigits();
339 // Possibly load FMD RecPoints information
340 if (TESTBIT(fTreeMask, kRecPoints)) {
341 // AliInfo("Getting FMD reconstructed points");
342 if (!fFMDLoader || fFMDLoader->LoadRecPoints("READ")) return kFALSE;
343 fTreeR = fFMDLoader->TreeR();
344 if (!fArrayR) fArrayR = new TClonesArray("AliFMDRecPoint");
345 fTreeR->SetBranchAddress("FMD", &fArrayR);
348 // Possibly load FMD ESD information
349 if (TESTBIT(fTreeMask, kESD)) {
350 // AliInfo("Getting FMD event summary data");
351 Int_t read = fChainE->GetEntry(event);
352 if (read <= 0) return kFALSE;
353 fESD = fESDEvent->GetFMDData();
354 if (!fESD) return kFALSE;
356 TFile* f = fChainE->GetFile();
358 TObject* o = f->GetStreamerInfoList()->FindObject("AliFMDMap");
360 TStreamerInfo* info = static_cast<TStreamerInfo*>(o);
361 std::cout << "AliFMDMap class version read is "
362 << info->GetClassVersion() << std::endl;
365 // fESD->CheckNeedUShort(fChainE->GetFile());
369 // Possibly load FMD Digit information
370 if (TESTBIT(fTreeMask, kRaw)) {
371 // AliInfo("Getting FMD raw data digits");
372 if (!fReader->NextEvent()) return kFALSE;
373 AliFMDRawReader r(fReader, 0);
382 //____________________________________________________________________
386 // Process one event. The default implementation one or more of
388 // - ProcessHits if the hits are loaded.
389 // - ProcessDigits if the digits are loaded.
390 // - ProcessSDigits if the sumbable digits are loaded.
391 // - ProcessRecPoints if the reconstructed points are loaded.
392 // - ProcessESD if the event summary data is loaded
394 if (TESTBIT(fTreeMask, kHits))
395 if (!ProcessHits()) return kFALSE;
396 if (TESTBIT(fTreeMask, kTracks))
397 if (!ProcessTracks()) return kFALSE;
398 if (TESTBIT(fTreeMask, kDigits))
399 if (!ProcessDigits()) return kFALSE;
400 if (TESTBIT(fTreeMask, kSDigits))
401 if (!ProcessSDigits()) return kFALSE;
402 if (TESTBIT(fTreeMask, kRaw))
403 if (!ProcessRawDigits()) return kFALSE;
404 if (TESTBIT(fTreeMask, kRecPoints))
405 if (!ProcessRecPoints()) return kFALSE;
406 if (TESTBIT(fTreeMask, kESD))
407 if (!ProcessESDs()) return kFALSE;
412 //____________________________________________________________________
414 AliFMDInput::ProcessHits()
416 // Read the hit tree, and pass each hit to the member function
419 AliError("No hit tree defined");
423 AliError("No hit array defined");
427 Int_t nTracks = fTreeH->GetEntries();
428 for (Int_t i = 0; i < nTracks; i++) {
429 Int_t hitRead = fTreeH->GetEntry(i);
430 if (hitRead <= 0) continue;
432 Int_t nHit = fArrayH->GetEntries();
433 if (nHit <= 0) continue;
435 for (Int_t j = 0; j < nHit; j++) {
436 AliFMDHit* hit = static_cast<AliFMDHit*>(fArrayH->At(j));
439 TParticle* track = 0;
440 if (TESTBIT(fTreeMask, kKinematics) && fStack) {
441 Int_t trackno = hit->Track();
442 track = fStack->Particle(trackno);
444 if (!ProcessHit(hit, track)) return kFALSE;
450 //____________________________________________________________________
452 AliFMDInput::ProcessTracks()
454 // Read the hit tree, and pass each hit to the member function
457 AliError("No track tree defined");
461 AliError("No hit tree defined");
465 AliError("No hit array defined");
469 // Int_t nTracks = fStack->GetNtrack();
470 Int_t nTracks = fTreeH->GetEntries();
471 for (Int_t i = 0; i < nTracks; i++) {
472 Int_t trackno = nTracks - i - 1;
473 TParticle* track = fStack->Particle(trackno);
474 if (!track) continue;
476 // Get the hits for this track.
477 Int_t hitRead = fTreeH->GetEntry(i);
478 Int_t nHit = fArrayH->GetEntries();
479 if (nHit == 0 || hitRead <= 0) {
480 // Let user code see the track, even if there's no hits.
481 if (!ProcessTrack(trackno, track, 0)) return kFALSE;
485 // Loop over the hits corresponding to this track.
486 for (Int_t j = 0; j < nHit; j++) {
487 AliFMDHit* hit = static_cast<AliFMDHit*>(fArrayH->At(j));
488 if (!ProcessTrack(trackno, track, hit)) return kFALSE;
494 //____________________________________________________________________
496 AliFMDInput::ProcessDigits()
498 // Read the digit tree, and pass each digit to the member function
501 AliError("No digit tree defined");
505 AliError("No digit array defined");
509 Int_t nEv = fTreeD->GetEntries();
510 for (Int_t i = 0; i < nEv; i++) {
511 Int_t digitRead = fTreeD->GetEntry(i);
512 if (digitRead <= 0) continue;
513 Int_t nDigit = fArrayD->GetEntries();
514 AliFMDDebug(0, ("Got %5d digits for this event", nDigit));
515 if (nDigit <= 0) continue;
516 for (Int_t j = 0; j < nDigit; j++) {
517 AliFMDDigit* digit = static_cast<AliFMDDigit*>(fArrayD->At(j));
518 if (!digit) continue;
519 if (!ProcessDigit(digit)) return kFALSE;
525 //____________________________________________________________________
527 AliFMDInput::ProcessSDigits()
529 // Read the summable digit tree, and pass each sumable digit to the
530 // member function ProcessSdigit.
532 AliWarning("No sdigit tree defined");
533 return kTRUE; // Empty SDigits is fine
536 AliWarning("No sdigit array defined");
537 return kTRUE; // Empty SDigits is fine
540 Int_t nEv = fTreeS->GetEntries();
541 for (Int_t i = 0; i < nEv; i++) {
542 Int_t sdigitRead = fTreeS->GetEntry(i);
543 if (sdigitRead <= 0) continue;
544 Int_t nSdigit = fArrayS->GetEntries();
545 AliFMDDebug(0, ("Got %5d digits for this event", nSdigit));
546 if (nSdigit <= 0) continue;
547 for (Int_t j = 0; j < nSdigit; j++) {
548 AliFMDSDigit* sdigit = static_cast<AliFMDSDigit*>(fArrayS->At(j));
549 if (!sdigit) continue;
550 if (!ProcessSDigit(sdigit)) return kFALSE;
556 //____________________________________________________________________
558 AliFMDInput::ProcessRawDigits()
560 // Read the digit tree, and pass each digit to the member function
563 AliError("No raw digit array defined");
567 Int_t nDigit = fArrayA->GetEntries();
568 if (nDigit <= 0) return kTRUE;
569 for (Int_t j = 0; j < nDigit; j++) {
570 AliFMDDigit* digit = static_cast<AliFMDDigit*>(fArrayA->At(j));
571 if (!digit) continue;
572 if (!ProcessRawDigit(digit)) return kFALSE;
577 //____________________________________________________________________
579 AliFMDInput::ProcessRecPoints()
581 // Read the reconstrcted points tree, and pass each reconstruction
582 // object (AliFMDRecPoint) to either ProcessRecPoint.
584 AliError("No recpoint tree defined");
588 AliError("No recpoints array defined");
592 Int_t nEv = fTreeR->GetEntries();
593 for (Int_t i = 0; i < nEv; i++) {
594 Int_t recRead = fTreeR->GetEntry(i);
595 if (recRead <= 0) continue;
596 Int_t nRecPoint = fArrayR->GetEntries();
597 for (Int_t j = 0; j < nRecPoint; j++) {
598 AliFMDRecPoint* recPoint = static_cast<AliFMDRecPoint*>(fArrayR->At(j));
599 if (!recPoint) continue;
600 if (!ProcessRecPoint(recPoint)) return kFALSE;
606 //____________________________________________________________________
608 AliFMDInput::ProcessESDs()
610 // Process event summary data
611 if (!fESD) return kFALSE;
612 for (UShort_t det = 1; det <= 3; det++) {
613 Char_t rings[] = { 'I', (det == 1 ? '\0' : 'O'), '\0' };
614 for (Char_t* rng = rings; *rng != '\0'; rng++) {
615 UShort_t nsec = (*rng == 'I' ? 20 : 40);
616 UShort_t nstr = (*rng == 'I' ? 512 : 256);
617 for (UShort_t sec = 0; sec < nsec; sec++) {
618 for (UShort_t str = 0; str < nstr; str++) {
619 Float_t eta = fESD->Eta(det,*rng,sec,str);
620 Float_t mult = fESD->Multiplicity(det,*rng,sec,str);
621 if (!fESD->IsAngleCorrected())
622 mult *= TMath::Abs(TMath::Cos(2.*TMath::ATan(TMath::Exp(-eta))));
623 if (!ProcessESD(det, *rng, sec, str, eta, mult)) continue;
631 //____________________________________________________________________
635 // Called at the end of each event. Per default, it unloads the
636 // data trees and resets the pointers to the output arrays. Users
637 // can overload this, but should call this member function in the
638 // overloaded member function of the derived class.
640 // Check if we have been initialized
642 AliError("Not initialized");
645 // Possibly unload global kinematics information
646 if (TESTBIT(fTreeMask, kKinematics) || TESTBIT(fTreeMask, kTracks)) {
647 fLoader->UnloadKinematics();
651 // Possibly unload FMD Hit information
652 if (TESTBIT(fTreeMask, kHits) || TESTBIT(fTreeMask, kTracks)) {
653 fFMDLoader->UnloadHits();
656 // Possibly unload FMD Digit information
657 if (TESTBIT(fTreeMask, kDigits)) {
658 fFMDLoader->UnloadDigits();
661 // Possibly unload FMD Sdigit information
662 if (TESTBIT(fTreeMask, kSDigits)) {
663 fFMDLoader->UnloadSDigits();
666 // Possibly unload FMD RecPoints information
667 if (TESTBIT(fTreeMask, kRecPoints)) {
668 fFMDLoader->UnloadRecPoints();
671 // AliInfo("Now out event");
675 //____________________________________________________________________
679 // Run over all events and files references in galice.root
682 if (!(retval = Init())) return retval;
684 Int_t nEvents = NEvents();
685 for (Int_t event = 0; event < nEvents; event++) {
686 if (!(retval = Begin(event))) break;
687 if (!(retval = Event())) break;
688 if (!(retval = End())) break;
690 if (!retval) return retval;
695 //__________________________________________________________________
697 AliFMDInput::MakeLogScale(Int_t n, Double_t min, Double_t max)
699 // Service function to define a logarithmic axis.
702 // min Minimum of axis
703 // max Maximum of axis
707 for (Int_t i = 1; i < n+1; i++) bins[i] = bins[i-1] + (max-min)/n;
710 Float_t dp = n / TMath::Log10(max / min);
711 Float_t pmin = TMath::Log10(min);
712 for (Int_t i = 1; i < n+1; i++) {
713 Float_t p = pmin + i / dp;
714 bins[i] = TMath::Power(10, p);
721 //____________________________________________________________________