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 "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
45 #include <AliESDFMD.h>
46 #include <AliESDEvent.h>
47 #include <AliCDBManager.h>
48 #include <AliCDBEntry.h>
49 #include <AliAlignObjParams.h>
50 #include <TTree.h> // ROOT_TTree
51 #include <TChain.h> // ROOT_TChain
52 #include <TParticle.h> // ROOT_TParticle
53 #include <TString.h> // ROOT_TString
54 #include <TDatabasePDG.h> // ROOT_TDatabasePDG
55 #include <TMath.h> // ROOT_TMath
56 #include <TGeoManager.h> // ROOT_TGeoManager
57 #include <TSystemDirectory.h> // ROOT_TSystemDirectory
58 #include <Riostream.h> // ROOT_Riostream
59 #include <TFile.h> // ROOT_TFile
60 #include <TStreamerInfo.h>
63 //____________________________________________________________________
66 ; // This is here to keep Emacs for indenting the next line
70 //____________________________________________________________________
71 AliFMDInput::AliFMDInput()
99 // Constructor of an FMD input object. Specify what data to read in
100 // using the AddLoad member function. Sub-classes should at a
101 // minimum overload the member function Event. A full job can be
102 // executed using the member function Run.
107 //____________________________________________________________________
108 AliFMDInput::AliFMDInput(const char* gAliceFile)
109 : fGAliceFile(gAliceFile),
136 // Constructor of an FMD input object. Specify what data to read in
137 // using the AddLoad member function. Sub-classes should at a
138 // minimum overload the member function Event. A full job can be
139 // executed using the member function Run.
142 //____________________________________________________________________
144 AliFMDInput::NEvents() const
146 // Get number of events
147 if (fTreeE) return fTreeE->GetEntries();
151 //____________________________________________________________________
155 // Initialize the object. Get the needed loaders, and such.
157 // Check if we have been initialized
159 AliWarning("Already initialized");
162 if (fGAliceFile.IsNull()) fGAliceFile = "galice.root";
164 fLoader = AliRunLoader::Open(fGAliceFile.Data(), "Alice", "read");
166 AliError(Form("Coulnd't read the file %s", fGAliceFile.Data()));
171 if (fLoader->LoadgAlice()) return kFALSE;
172 fRun = fLoader->GetAliRun();
175 fFMD = static_cast<AliFMD*>(fRun->GetDetector("FMD"));
177 AliError("Failed to get detector FMD from loader");
181 // Get the FMD loader
182 fFMDLoader = fLoader->GetLoader("FMDLoader");
184 AliError("Failed to get detector FMD loader from loader");
187 if (fLoader->LoadHeader()) {
188 AliError("Failed to get event header information from loader");
191 fTreeE = fLoader->TreeE();
193 // Optionally, get the ESD files
194 if (TESTBIT(fTreeMask, kESD)) {
195 fChainE = new TChain("esdTree");
196 TSystemDirectory dir(".",".");
197 TList* files = dir.GetListOfFiles();
198 TSystemFile* file = 0;
200 AliError("No files");
205 while ((file = static_cast<TSystemFile*>(next()))) {
206 TString fname(file->GetName());
207 if (fname.Contains("AliESDs")) fChainE->AddFile(fname.Data());
209 fESDEvent = new AliESDEvent();
210 fESDEvent->ReadFromTree(fChainE);
211 // fChainE->SetBranchAddress("ESD", &fMainESD);
215 if (TESTBIT(fTreeMask, kRaw)) {
216 AliInfo("Getting FMD raw data digits");
217 fArrayA = new TClonesArray("AliFMDDigit");
218 fReader = new AliRawReaderFile(-1);
221 // Optionally, get the geometry
222 if (TESTBIT(fTreeMask, kGeometry)) {
223 TString fname(fRun->GetGeometryFileName());
224 if (fname.IsNull()) {
225 Warning("Init", "No file name for the geometry from AliRun");
226 fname = gSystem->DirName(fGAliceFile);
227 fname.Append("/geometry.root");
229 fGeoManager = TGeoManager::Import(fname.Data());
231 Fatal("Init", "No geometry manager found");
234 AliCDBManager* cdb = AliCDBManager::Instance();
235 AliCDBEntry* align = cdb->Get("FMD/Align/Data");
237 AliInfo("Got alignment data from CDB");
238 TClonesArray* array = dynamic_cast<TClonesArray*>(align->GetObject());
240 AliWarning("Invalid align data from CDB");
243 Int_t nAlign = array->GetEntries();
244 for (Int_t i = 0; i < nAlign; i++) {
245 AliAlignObjParams* a = static_cast<AliAlignObjParams*>(array->At(i));
246 if (!a->ApplyToGeometry()) {
247 AliWarning(Form("Failed to apply alignment to %s",
260 //____________________________________________________________________
262 AliFMDInput::Begin(Int_t event)
264 // Called at the begining of each event. Per default, it gets the
265 // data trees and gets pointers to the output arrays. Users can
266 // overload this, but should call this member function in the
267 // overloaded member function of the derived class.
269 // Check if we have been initialized
271 AliError("Not initialized");
276 if (fLoader->GetEvent(event)) return kFALSE;
277 AliInfo(Form("Now in event %8d/%8d", event, NEvents()));
279 // Possibly load global kinematics information
280 if (TESTBIT(fTreeMask, kKinematics) || TESTBIT(fTreeMask, kTracks)) {
281 // AliInfo("Getting kinematics");
282 if (fLoader->LoadKinematics()) return kFALSE;
283 fStack = fLoader->Stack();
286 // Possibly load FMD Hit information
287 if (TESTBIT(fTreeMask, kHits) || TESTBIT(fTreeMask, kTracks)) {
288 // AliInfo("Getting FMD hits");
289 if (!fFMDLoader || fFMDLoader->LoadHits()) return kFALSE;
290 fTreeH = fFMDLoader->TreeH();
291 if (!fArrayH) fArrayH = fFMD->Hits();
294 // Possibly load heaedr information
295 if (TESTBIT(fTreeMask, kHeader)) {
296 // AliInfo("Getting FMD hits");
297 if (!fLoader /* || fLoader->LoadHeader()*/) return kFALSE;
298 fHeader = fLoader->GetHeader();
301 // Possibly load FMD Digit information
302 if (TESTBIT(fTreeMask, kDigits)) {
303 // AliInfo("Getting FMD digits");
304 if (!fFMDLoader || fFMDLoader->LoadDigits()) return kFALSE;
305 fTreeD = fFMDLoader->TreeD();
307 if (!fArrayD) fArrayD = fFMD->Digits();
311 AliWarning(Form("Failed to load FMD Digits"));
315 // Possibly load FMD Sdigit information
316 if (TESTBIT(fTreeMask, kSDigits)) {
317 // AliInfo("Getting FMD summable digits");
318 if (!fFMDLoader || fFMDLoader->LoadSDigits()) return kFALSE;
319 fTreeS = fFMDLoader->TreeS();
320 if (!fArrayS) fArrayS = fFMD->SDigits();
323 // Possibly load FMD RecPoints information
324 if (TESTBIT(fTreeMask, kRecPoints)) {
325 // AliInfo("Getting FMD reconstructed points");
326 if (!fFMDLoader || fFMDLoader->LoadRecPoints()) return kFALSE;
327 fTreeR = fFMDLoader->TreeR();
328 if (!fArrayR) fArrayR = new TClonesArray("AliFMDRecPoint");
329 fTreeR->SetBranchAddress("FMD", &fArrayR);
332 // Possibly load FMD ESD information
333 if (TESTBIT(fTreeMask, kESD)) {
334 // AliInfo("Getting FMD event summary data");
335 Int_t read = fChainE->GetEntry(event);
336 if (read <= 0) return kFALSE;
337 fESD = fESDEvent->GetFMDData();
338 if (!fESD) return kFALSE;
340 TFile* f = fChainE->GetFile();
342 TObject* o = f->GetStreamerInfoList()->FindObject("AliFMDMap");
344 TStreamerInfo* info = static_cast<TStreamerInfo*>(o);
345 std::cout << "AliFMDMap class version read is "
346 << info->GetClassVersion() << std::endl;
349 // fESD->CheckNeedUShort(fChainE->GetFile());
353 // Possibly load FMD Digit information
354 if (TESTBIT(fTreeMask, kRaw)) {
355 // AliInfo("Getting FMD raw data digits");
356 if (!fReader->NextEvent()) return kFALSE;
357 AliFMDRawReader r(fReader, 0);
366 //____________________________________________________________________
370 // Process one event. The default implementation one or more of
372 // - ProcessHits if the hits are loaded.
373 // - ProcessDigits if the digits are loaded.
374 // - ProcessSDigits if the sumbable digits are loaded.
375 // - ProcessRecPoints if the reconstructed points are loaded.
376 // - ProcessESD if the event summary data is loaded
378 if (TESTBIT(fTreeMask, kHits))
379 if (!ProcessHits()) return kFALSE;
380 if (TESTBIT(fTreeMask, kTracks))
381 if (!ProcessTracks()) return kFALSE;
382 if (TESTBIT(fTreeMask, kDigits))
383 if (!ProcessDigits()) return kFALSE;
384 if (TESTBIT(fTreeMask, kSDigits))
385 if (!ProcessSDigits()) return kFALSE;
386 if (TESTBIT(fTreeMask, kRaw))
387 if (!ProcessRawDigits()) return kFALSE;
388 if (TESTBIT(fTreeMask, kRecPoints))
389 if (!ProcessRecPoints()) return kFALSE;
390 if (TESTBIT(fTreeMask, kESD))
391 if (!ProcessESDs()) return kFALSE;
396 //____________________________________________________________________
398 AliFMDInput::ProcessHits()
400 // Read the hit tree, and pass each hit to the member function
403 AliError("No hit tree defined");
407 AliError("No hit array defined");
411 Int_t nTracks = fTreeH->GetEntries();
412 for (Int_t i = 0; i < nTracks; i++) {
413 Int_t hitRead = fTreeH->GetEntry(i);
414 if (hitRead <= 0) continue;
416 Int_t nHit = fArrayH->GetEntries();
417 if (nHit <= 0) continue;
419 for (Int_t j = 0; j < nHit; j++) {
420 AliFMDHit* hit = static_cast<AliFMDHit*>(fArrayH->At(j));
423 TParticle* track = 0;
424 if (TESTBIT(fTreeMask, kKinematics) && fStack) {
425 Int_t trackno = hit->Track();
426 track = fStack->Particle(trackno);
428 if (!ProcessHit(hit, track)) return kFALSE;
434 //____________________________________________________________________
436 AliFMDInput::ProcessTracks()
438 // Read the hit tree, and pass each hit to the member function
441 AliError("No track tree defined");
445 AliError("No hit tree defined");
449 AliError("No hit array defined");
453 // Int_t nTracks = fStack->GetNtrack();
454 Int_t nTracks = fTreeH->GetEntries();
455 for (Int_t i = 0; i < nTracks; i++) {
456 Int_t trackno = nTracks - i - 1;
457 TParticle* track = fStack->Particle(trackno);
458 if (!track) continue;
460 // Get the hits for this track.
461 Int_t hitRead = fTreeH->GetEntry(i);
462 Int_t nHit = fArrayH->GetEntries();
463 if (nHit == 0 || hitRead <= 0) {
464 // Let user code see the track, even if there's no hits.
465 if (!ProcessTrack(trackno, track, 0)) return kFALSE;
469 // Loop over the hits corresponding to this track.
470 for (Int_t j = 0; j < nHit; j++) {
471 AliFMDHit* hit = static_cast<AliFMDHit*>(fArrayH->At(j));
472 if (!ProcessTrack(trackno, track, hit)) return kFALSE;
478 //____________________________________________________________________
480 AliFMDInput::ProcessDigits()
482 // Read the digit tree, and pass each digit to the member function
484 Int_t nEv = fTreeD->GetEntries();
485 for (Int_t i = 0; i < nEv; i++) {
486 Int_t digitRead = fTreeD->GetEntry(i);
487 if (digitRead <= 0) continue;
488 Int_t nDigit = fArrayD->GetEntries();
489 if (nDigit <= 0) continue;
490 for (Int_t j = 0; j < nDigit; j++) {
491 AliFMDDigit* digit = static_cast<AliFMDDigit*>(fArrayD->At(j));
492 if (!digit) continue;
493 if (!ProcessDigit(digit)) return kFALSE;
499 //____________________________________________________________________
501 AliFMDInput::ProcessSDigits()
503 // Read the summable digit tree, and pass each sumable digit to the
504 // member function ProcessSdigit.
505 Int_t nEv = fTreeD->GetEntries();
506 for (Int_t i = 0; i < nEv; i++) {
507 Int_t sdigitRead = fTreeS->GetEntry(i);
508 if (sdigitRead <= 0) continue;
509 Int_t nSdigit = fArrayS->GetEntries();
510 if (nSdigit <= 0) continue;
511 for (Int_t j = 0; j < nSdigit; j++) {
512 AliFMDSDigit* sdigit = static_cast<AliFMDSDigit*>(fArrayS->At(j));
513 if (!sdigit) continue;
514 if (!ProcessSDigit(sdigit)) return kFALSE;
520 //____________________________________________________________________
522 AliFMDInput::ProcessRawDigits()
524 // Read the digit tree, and pass each digit to the member function
526 Int_t nDigit = fArrayA->GetEntries();
527 if (nDigit <= 0) return kTRUE;
528 for (Int_t j = 0; j < nDigit; j++) {
529 AliFMDDigit* digit = static_cast<AliFMDDigit*>(fArrayA->At(j));
530 if (!digit) continue;
531 if (!ProcessRawDigit(digit)) return kFALSE;
536 //____________________________________________________________________
538 AliFMDInput::ProcessRecPoints()
540 // Read the reconstrcted points tree, and pass each reconstruction
541 // object (AliFMDRecPoint) to either ProcessRecPoint.
542 Int_t nEv = fTreeR->GetEntries();
543 for (Int_t i = 0; i < nEv; i++) {
544 Int_t recRead = fTreeR->GetEntry(i);
545 if (recRead <= 0) continue;
546 Int_t nRecPoint = fArrayR->GetEntries();
547 for (Int_t j = 0; j < nRecPoint; j++) {
548 AliFMDRecPoint* recPoint = static_cast<AliFMDRecPoint*>(fArrayR->At(j));
549 if (!recPoint) continue;
550 if (!ProcessRecPoint(recPoint)) return kFALSE;
556 //____________________________________________________________________
558 AliFMDInput::ProcessESDs()
560 // Process event summary data
561 if (!fESD) return kFALSE;
562 for (UShort_t det = 1; det <= 3; det++) {
563 Char_t rings[] = { 'I', (det == 1 ? '\0' : 'O'), '\0' };
564 for (Char_t* rng = rings; *rng != '\0'; rng++) {
565 UShort_t nsec = (*rng == 'I' ? 20 : 40);
566 UShort_t nstr = (*rng == 'I' ? 512 : 256);
567 for (UShort_t sec = 0; sec < nsec; sec++) {
568 for (UShort_t str = 0; str < nstr; str++) {
569 Float_t eta = fESD->Eta(det,*rng,sec,str);
570 Float_t mult = fESD->Multiplicity(det,*rng,sec,str);
571 if (!fESD->IsAngleCorrected())
572 mult *= TMath::Abs(TMath::Cos(2.*TMath::ATan(TMath::Exp(-eta))));
573 if (!ProcessESD(det, *rng, sec, str, eta, mult)) continue;
581 //____________________________________________________________________
585 // Called at the end of each event. Per default, it unloads the
586 // data trees and resets the pointers to the output arrays. Users
587 // can overload this, but should call this member function in the
588 // overloaded member function of the derived class.
590 // Check if we have been initialized
592 AliError("Not initialized");
595 // Possibly unload global kinematics information
596 if (TESTBIT(fTreeMask, kKinematics) || TESTBIT(fTreeMask, kTracks)) {
597 fLoader->UnloadKinematics();
601 // Possibly unload FMD Hit information
602 if (TESTBIT(fTreeMask, kHits) || TESTBIT(fTreeMask, kTracks)) {
603 fFMDLoader->UnloadHits();
606 // Possibly unload FMD Digit information
607 if (TESTBIT(fTreeMask, kDigits)) {
608 fFMDLoader->UnloadDigits();
611 // Possibly unload FMD Sdigit information
612 if (TESTBIT(fTreeMask, kSDigits)) {
613 fFMDLoader->UnloadSDigits();
616 // Possibly unload FMD RecPoints information
617 if (TESTBIT(fTreeMask, kRecPoints)) {
618 fFMDLoader->UnloadRecPoints();
621 // AliInfo("Now out event");
625 //____________________________________________________________________
629 // Run over all events and files references in galice.root
632 if (!(retval = Init())) return retval;
634 Int_t nEvents = NEvents();
635 for (Int_t event = 0; event < nEvents; event++) {
636 if (!(retval = Begin(event))) break;
637 if (!(retval = Event())) break;
638 if (!(retval = End())) break;
640 if (!retval) return retval;
645 //__________________________________________________________________
647 AliFMDInput::MakeLogScale(Int_t n, Double_t min, Double_t max)
649 // Service function to define a logarithmic axis.
652 // min Minimum of axis
653 // max Maximum of axis
657 for (Int_t i = 1; i < n+1; i++) bins[i] = bins[i-1] + (max-min)/n;
660 Float_t dp = n / TMath::Log10(max / min);
661 Float_t pmin = TMath::Log10(min);
662 for (Int_t i = 1; i < n+1; i++) {
663 Float_t p = pmin + i / dp;
664 bins[i] = TMath::Power(10, p);
671 //____________________________________________________________________