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 "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
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
59 //____________________________________________________________________
62 ; // This is here to keep Emacs for indenting the next line
66 //____________________________________________________________________
67 AliFMDInput::AliFMDInput()
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.
100 //____________________________________________________________________
101 AliFMDInput::AliFMDInput(const char* gAliceFile)
102 : fGAliceFile(gAliceFile),
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.
132 //____________________________________________________________________
134 AliFMDInput::NEvents() const
136 // Get number of events
137 if (fTreeE) return fTreeE->GetEntries();
141 //____________________________________________________________________
145 // Initialize the object. Get the needed loaders, and such.
147 // Check if we have been initialized
149 AliWarning("Already initialized");
152 if (fGAliceFile.IsNull()) fGAliceFile = "galice.root";
154 fLoader = AliRunLoader::Open(fGAliceFile.Data(), "Alice", "read");
156 AliError(Form("Coulnd't read the file %s", fGAliceFile.Data()));
161 if (fLoader->LoadgAlice()) return kFALSE;
162 fRun = fLoader->GetAliRun();
165 fFMD = static_cast<AliFMD*>(fRun->GetDetector("FMD"));
167 AliError("Failed to get detector FMD from loader");
171 // Get the FMD loader
172 fFMDLoader = fLoader->GetLoader("FMDLoader");
174 AliError("Failed to get detector FMD loader from loader");
177 if (fLoader->LoadHeader()) {
178 AliError("Failed to get event header information from loader");
181 fTreeE = fLoader->TreeE();
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;
190 AliError("No files");
195 while ((file = static_cast<TSystemFile*>(next()))) {
196 TString fname(file->GetName());
197 if (fname.Contains("AliESDs")) fChainE->AddFile(fname.Data());
199 fChainE->SetBranchAddress("ESD", &fMainESD);
202 if (TESTBIT(fTreeMask, kRaw)) {
203 AliInfo("Getting FMD raw data digits");
204 fArrayA = new TClonesArray("AliFMDDigit");
205 fReader = new AliRawReaderFile(-1);
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");
216 fGeoManager = TGeoManager::Import(fname.Data());
218 Fatal("Init", "No geometry manager found");
221 AliCDBManager* cdb = AliCDBManager::Instance();
222 AliCDBEntry* align = cdb->Get("FMD/Align/Data");
224 AliInfo("Got alignment data from CDB");
225 TClonesArray* array = dynamic_cast<TClonesArray*>(align->GetObject());
227 AliWarning("Invalid align data from CDB");
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",
247 //____________________________________________________________________
249 AliFMDInput::Begin(Int_t event)
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.
256 // Check if we have been initialized
258 AliError("Not initialized");
262 if (fLoader->GetEvent(event)) return kFALSE;
263 AliInfo(Form("Now in event %d/%d", event, NEvents()));
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();
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();
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();
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();
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;
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);
320 //____________________________________________________________________
324 // Process one event. The default implementation one or more of
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
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;
348 //____________________________________________________________________
350 AliFMDInput::ProcessHits()
352 // Read the hit tree, and pass each hit to the member function
355 AliError("No hit tree defined");
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;
363 AliError("No hit array defined");
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));
371 TParticle* track = 0;
372 if (TESTBIT(fTreeMask, kKinematics) && fStack) {
373 Int_t trackno = hit->Track();
374 track = fStack->Particle(trackno);
376 if (!ProcessHit(hit, track)) return kFALSE;
382 //____________________________________________________________________
384 AliFMDInput::ProcessDigits()
386 // Read the digit tree, and pass each digit to the member function
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;
403 //____________________________________________________________________
405 AliFMDInput::ProcessSDigits()
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;
424 //____________________________________________________________________
426 AliFMDInput::ProcessRawDigits()
428 // Read the digit tree, and pass each digit to the member function
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;
440 //____________________________________________________________________
442 AliFMDInput::ProcessRecPoints()
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;
460 //____________________________________________________________________
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.
469 // Check if we have been initialized
471 AliError("Not initialized");
474 // Possibly unload global kinematics information
475 if (TESTBIT(fTreeMask, kKinematics)) {
476 fLoader->UnloadKinematics();
480 // Possibly unload FMD Hit information
481 if (TESTBIT(fTreeMask, kHits)) {
482 fFMDLoader->UnloadHits();
485 // Possibly unload FMD Digit information
486 if (TESTBIT(fTreeMask, kDigits)) {
487 fFMDLoader->UnloadDigits();
490 // Possibly unload FMD Sdigit information
491 if (TESTBIT(fTreeMask, kSDigits)) {
492 fFMDLoader->UnloadSDigits();
495 // Possibly unload FMD RecPoints information
496 if (TESTBIT(fTreeMask, kRecPoints)) {
497 fFMDLoader->UnloadRecPoints();
500 AliInfo("Now out event");
504 //____________________________________________________________________
508 // Run over all events and files references in galice.root
511 if (!(retval = Init())) return retval;
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;
519 if (!retval) return retval;
525 //____________________________________________________________________