1 /**************************************************************************
2 * Copyright(c) 1998-1999, 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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for running the reconstruction //
22 // Clusters and tracks are created for all detectors and all events by //
25 // AliReconstruction rec; //
28 // The Run method returns kTRUE in case of successful execution. //
30 // If the input to the reconstruction are not simulated digits but raw data, //
31 // this can be specified by an argument of the Run method or by the method //
33 // rec.SetInput("..."); //
35 // The input formats and the corresponding argument are: //
36 // - DDL raw data files: directory name, ends with "/" //
37 // - raw data root file: root file name, extension ".root" //
38 // - raw data DATE file: DATE file name, any other non-empty string //
39 // - MC root files : empty string, default //
41 // By default all events are reconstructed. The reconstruction can be //
42 // limited to a range of events by giving the index of the first and the //
43 // last event as an argument to the Run method or by calling //
45 // rec.SetEventRange(..., ...); //
47 // The index -1 (default) can be used for the last event to indicate no //
48 // upper limit of the event range. //
50 // In case of raw-data reconstruction the user can modify the default //
51 // number of events per digits/clusters/tracks file. In case the option //
52 // is not used the number is set 1. In case the user provides 0, than //
53 // the number of events is equal to the number of events inside the //
54 // raw-data file (i.e. one digits/clusters/tracks file): //
56 // rec.SetNumberOfEventsPerFile(...); //
59 // The name of the galice file can be changed from the default //
60 // "galice.root" by passing it as argument to the AliReconstruction //
61 // constructor or by //
63 // rec.SetGAliceFile("..."); //
65 // The local reconstruction can be switched on or off for individual //
68 // rec.SetRunLocalReconstruction("..."); //
70 // The argument is a (case sensitive) string with the names of the //
71 // detectors separated by a space. The special string "ALL" selects all //
72 // available detectors. This is the default. //
74 // The reconstruction of the primary vertex position can be switched off by //
76 // rec.SetRunVertexFinder(kFALSE); //
78 // The tracking and the creation of ESD tracks can be switched on for //
79 // selected detectors by //
81 // rec.SetRunTracking("..."); //
83 // Uniform/nonuniform field tracking switches (default: uniform field) //
85 // rec.SetUniformFieldTracking(); ( rec.SetUniformFieldTracking(kFALSE); ) //
87 // The filling of additional ESD information can be steered by //
89 // rec.SetFillESD("..."); //
91 // Again, for both methods the string specifies the list of detectors. //
92 // The default is "ALL". //
94 // The call of the shortcut method //
96 // rec.SetRunReconstruction("..."); //
98 // is equivalent to calling SetRunLocalReconstruction, SetRunTracking and //
99 // SetFillESD with the same detector selecting string as argument. //
101 // The reconstruction requires digits or raw data as input. For the creation //
102 // of digits and raw data have a look at the class AliSimulation. //
104 // The input data of a detector can be replaced by the corresponding HLT //
105 // data by calling (usual detector string) //
106 // SetUseHLTData("..."); //
109 ///////////////////////////////////////////////////////////////////////////////
116 #include <TGeoGlobalMagField.h>
117 #include <TGeoManager.h>
119 #include <TLorentzVector.h>
121 #include <TObjArray.h>
122 #include <TPRegexp.h>
123 #include <TParameter.h>
124 #include <TPluginManager.h>
126 #include <TProofOutputFile.h>
130 #include "AliAlignObj.h"
131 #include "AliCDBEntry.h"
132 #include "AliCDBManager.h"
133 #include "AliCDBStorage.h"
134 #include "AliCTPRawStream.h"
135 #include "AliCascadeVertexer.h"
136 #include "AliCentralTrigger.h"
137 #include "AliCodeTimer.h"
139 #include "AliDetectorRecoParam.h"
140 #include "AliESDCaloCells.h"
141 #include "AliESDCaloCluster.h"
142 #include "AliESDEvent.h"
143 #include "AliESDMuonTrack.h"
144 #include "AliESDPmdTrack.h"
145 #include "AliESDTagCreator.h"
146 #include "AliESDVertex.h"
147 #include "AliESDcascade.h"
148 #include "AliESDfriend.h"
149 #include "AliESDkink.h"
150 #include "AliESDpid.h"
151 #include "AliESDtrack.h"
152 #include "AliESDtrack.h"
153 #include "AliEventInfo.h"
154 #include "AliGRPObject.h"
155 #include "AliGRPRecoParam.h"
156 #include "AliGenEventHeader.h"
157 #include "AliGeomManager.h"
158 #include "AliGlobalQADataMaker.h"
159 #include "AliHeader.h"
162 #include "AliMultiplicity.h"
164 #include "AliPlaneEff.h"
166 #include "AliQADataMakerRec.h"
167 #include "AliQAManager.h"
168 #include "AliRawEvent.h"
169 #include "AliRawEventHeaderBase.h"
170 #include "AliRawHLTManager.h"
171 #include "AliRawReaderDate.h"
172 #include "AliRawReaderFile.h"
173 #include "AliRawReaderRoot.h"
174 #include "AliReconstruction.h"
175 #include "AliReconstructor.h"
177 #include "AliRunInfo.h"
178 #include "AliRunLoader.h"
179 #include "AliSysInfo.h" // memory snapshots
180 #include "AliTrackPointArray.h"
181 #include "AliTracker.h"
182 #include "AliTriggerClass.h"
183 #include "AliTriggerCluster.h"
184 #include "AliTriggerConfiguration.h"
185 #include "AliV0vertexer.h"
186 #include "AliVertexer.h"
187 #include "AliVertexerTracks.h"
189 ClassImp(AliReconstruction)
191 //_____________________________________________________________________________
192 const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
194 //_____________________________________________________________________________
195 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
197 fUniformField(kFALSE),
198 fRunVertexFinder(kTRUE),
199 fRunVertexFinderTracks(kTRUE),
200 fRunHLTTracking(kFALSE),
201 fRunMuonTracking(kFALSE),
203 fRunCascadeFinder(kTRUE),
204 fStopOnError(kFALSE),
205 fWriteAlignmentData(kFALSE),
206 fWriteESDfriend(kFALSE),
207 fFillTriggerESD(kTRUE),
215 fRunLocalReconstruction("ALL"),
219 fUseTrackingErrorsForAlignment(""),
220 fGAliceFileName(gAliceFilename),
225 fNumberOfEventsPerFile((UInt_t)-1),
227 fLoadAlignFromCDB(kTRUE),
228 fLoadAlignData("ALL"),
235 fParentRawReader(NULL),
239 fDiamondProfileSPD(NULL),
240 fDiamondProfile(NULL),
241 fDiamondProfileTPC(NULL),
245 fAlignObjArray(NULL),
248 fInitCDBCalled(kFALSE),
249 fSetRunNumberFromDataCalled(kFALSE),
255 fSameQACycle(kFALSE),
257 fRunPlaneEff(kFALSE),
266 fIsNewRunLoader(kFALSE),
270 // create reconstruction object with default parameters
273 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
274 fReconstructor[iDet] = NULL;
275 fLoader[iDet] = NULL;
276 fTracker[iDet] = NULL;
278 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
279 fQACycles[iDet] = 999999 ;
280 fQAWriteExpert[iDet] = kFALSE ;
286 //_____________________________________________________________________________
287 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
289 fUniformField(rec.fUniformField),
290 fRunVertexFinder(rec.fRunVertexFinder),
291 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
292 fRunHLTTracking(rec.fRunHLTTracking),
293 fRunMuonTracking(rec.fRunMuonTracking),
294 fRunV0Finder(rec.fRunV0Finder),
295 fRunCascadeFinder(rec.fRunCascadeFinder),
296 fStopOnError(rec.fStopOnError),
297 fWriteAlignmentData(rec.fWriteAlignmentData),
298 fWriteESDfriend(rec.fWriteESDfriend),
299 fFillTriggerESD(rec.fFillTriggerESD),
301 fCleanESD(rec.fCleanESD),
302 fV0DCAmax(rec.fV0DCAmax),
303 fV0CsPmin(rec.fV0CsPmin),
307 fRunLocalReconstruction(rec.fRunLocalReconstruction),
308 fRunTracking(rec.fRunTracking),
309 fFillESD(rec.fFillESD),
310 fLoadCDB(rec.fLoadCDB),
311 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
312 fGAliceFileName(rec.fGAliceFileName),
313 fRawInput(rec.fRawInput),
314 fEquipIdMap(rec.fEquipIdMap),
315 fFirstEvent(rec.fFirstEvent),
316 fLastEvent(rec.fLastEvent),
317 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
319 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
320 fLoadAlignData(rec.fLoadAlignData),
321 fUseHLTData(rec.fUseHLTData),
327 fParentRawReader(NULL),
329 fRecoParam(rec.fRecoParam),
331 fDiamondProfileSPD(rec.fDiamondProfileSPD),
332 fDiamondProfile(rec.fDiamondProfile),
333 fDiamondProfileTPC(rec.fDiamondProfileTPC),
337 fAlignObjArray(rec.fAlignObjArray),
338 fCDBUri(rec.fCDBUri),
340 fInitCDBCalled(rec.fInitCDBCalled),
341 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
342 fQADetectors(rec.fQADetectors),
344 fQATasks(rec.fQATasks),
346 fRunGlobalQA(rec.fRunGlobalQA),
347 fSameQACycle(rec.fSameQACycle),
348 fRunPlaneEff(rec.fRunPlaneEff),
357 fIsNewRunLoader(rec.fIsNewRunLoader),
363 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
364 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
366 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
367 fReconstructor[iDet] = NULL;
368 fLoader[iDet] = NULL;
369 fTracker[iDet] = NULL;
372 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
373 fQACycles[iDet] = rec.fQACycles[iDet];
374 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
377 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
378 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
382 //_____________________________________________________________________________
383 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
385 // assignment operator
386 // Used in PROOF mode
387 // Be very careful while modifing it!
388 // Simple rules to follow:
389 // for persistent data members - use their assignment operators
390 // for non-persistent ones - do nothing or take the default values from constructor
391 // TSelector members should not be touched
392 if(&rec == this) return *this;
394 fUniformField = rec.fUniformField;
395 fRunVertexFinder = rec.fRunVertexFinder;
396 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
397 fRunHLTTracking = rec.fRunHLTTracking;
398 fRunMuonTracking = rec.fRunMuonTracking;
399 fRunV0Finder = rec.fRunV0Finder;
400 fRunCascadeFinder = rec.fRunCascadeFinder;
401 fStopOnError = rec.fStopOnError;
402 fWriteAlignmentData = rec.fWriteAlignmentData;
403 fWriteESDfriend = rec.fWriteESDfriend;
404 fFillTriggerESD = rec.fFillTriggerESD;
406 fCleanESD = rec.fCleanESD;
407 fV0DCAmax = rec.fV0DCAmax;
408 fV0CsPmin = rec.fV0CsPmin;
412 fRunLocalReconstruction = rec.fRunLocalReconstruction;
413 fRunTracking = rec.fRunTracking;
414 fFillESD = rec.fFillESD;
415 fLoadCDB = rec.fLoadCDB;
416 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
417 fGAliceFileName = rec.fGAliceFileName;
418 fRawInput = rec.fRawInput;
419 fEquipIdMap = rec.fEquipIdMap;
420 fFirstEvent = rec.fFirstEvent;
421 fLastEvent = rec.fLastEvent;
422 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
424 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
425 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
428 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
429 fLoadAlignData = rec.fLoadAlignData;
430 fUseHLTData = rec.fUseHLTData;
432 delete fRunInfo; fRunInfo = NULL;
433 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
435 fEventInfo = rec.fEventInfo;
439 fParentRawReader = NULL;
441 fRecoParam = rec.fRecoParam;
443 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
444 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
445 delete fLoader[iDet]; fLoader[iDet] = NULL;
446 delete fTracker[iDet]; fTracker[iDet] = NULL;
449 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
450 fQACycles[iDet] = rec.fQACycles[iDet];
451 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
454 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
455 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
456 delete fDiamondProfile; fDiamondProfile = NULL;
457 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
458 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
459 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
461 delete fGRPData; fGRPData = NULL;
462 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
463 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
465 delete fAlignObjArray; fAlignObjArray = NULL;
468 fSpecCDBUri.Delete();
469 fInitCDBCalled = rec.fInitCDBCalled;
470 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
471 fQADetectors = rec.fQADetectors;
473 fQATasks = rec.fQATasks;
475 fRunGlobalQA = rec.fRunGlobalQA;
476 fSameQACycle = rec.fSameQACycle;
477 fRunPlaneEff = rec.fRunPlaneEff;
486 fIsNewRunLoader = rec.fIsNewRunLoader;
493 //_____________________________________________________________________________
494 AliReconstruction::~AliReconstruction()
501 if (fAlignObjArray) {
502 fAlignObjArray->Delete();
503 delete fAlignObjArray;
505 fSpecCDBUri.Delete();
507 AliCodeTimer::Instance()->Print();
510 //_____________________________________________________________________________
511 void AliReconstruction::InitCDB()
513 // activate a default CDB storage
514 // First check if we have any CDB storage set, because it is used
515 // to retrieve the calibration and alignment constants
516 AliCodeTimerAuto("");
518 if (fInitCDBCalled) return;
519 fInitCDBCalled = kTRUE;
521 AliCDBManager* man = AliCDBManager::Instance();
522 if (man->IsDefaultStorageSet())
524 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
525 AliWarning("Default CDB storage has been already set !");
526 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
527 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
528 fCDBUri = man->GetDefaultStorage()->GetURI();
531 if (fCDBUri.Length() > 0)
533 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
534 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
535 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
537 fCDBUri="local://$ALICE_ROOT/OCDB";
538 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
539 AliWarning("Default CDB storage not yet set !!!!");
540 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
541 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
544 man->SetDefaultStorage(fCDBUri);
547 // Now activate the detector specific CDB storage locations
548 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
549 TObject* obj = fSpecCDBUri[i];
551 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
552 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
553 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
554 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
556 AliSysInfo::AddStamp("InitCDB");
559 //_____________________________________________________________________________
560 void AliReconstruction::SetDefaultStorage(const char* uri) {
561 // Store the desired default CDB storage location
562 // Activate it later within the Run() method
568 //_____________________________________________________________________________
569 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
570 // Store a detector-specific CDB storage location
571 // Activate it later within the Run() method
573 AliCDBPath aPath(calibType);
574 if(!aPath.IsValid()){
575 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
576 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
577 if(!strcmp(calibType, fgkDetectorName[iDet])) {
578 aPath.SetPath(Form("%s/*", calibType));
579 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
583 if(!aPath.IsValid()){
584 AliError(Form("Not a valid path or detector: %s", calibType));
589 // // check that calibType refers to a "valid" detector name
590 // Bool_t isDetector = kFALSE;
591 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
592 // TString detName = fgkDetectorName[iDet];
593 // if(aPath.GetLevel0() == detName) {
594 // isDetector = kTRUE;
600 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
604 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
605 if (obj) fSpecCDBUri.Remove(obj);
606 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
610 //_____________________________________________________________________________
611 Bool_t AliReconstruction::SetRunNumberFromData()
613 // The method is called in Run() in order
614 // to set a correct run number.
615 // In case of raw data reconstruction the
616 // run number is taken from the raw data header
618 if (fSetRunNumberFromDataCalled) return kTRUE;
619 fSetRunNumberFromDataCalled = kTRUE;
621 AliCDBManager* man = AliCDBManager::Instance();
624 if(fRawReader->NextEvent()) {
625 if(man->GetRun() > 0) {
626 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
628 man->SetRun(fRawReader->GetRunNumber());
629 fRawReader->RewindEvents();
632 if(man->GetRun() > 0) {
633 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
636 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
642 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
644 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
649 // read run number from gAlice
650 if(rl->GetHeader()) {
651 man->SetRun(rl->GetHeader()->GetRun());
656 AliError("Neither run-loader header nor RawReader objects are found !");
668 //_____________________________________________________________________________
669 void AliReconstruction::SetCDBLock() {
670 // Set CDB lock: from now on it is forbidden to reset the run number
671 // or the default storage or to activate any further storage!
673 AliCDBManager::Instance()->SetLock(1);
676 //_____________________________________________________________________________
677 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
679 // Read the alignment objects from CDB.
680 // Each detector is supposed to have the
681 // alignment objects in DET/Align/Data CDB path.
682 // All the detector objects are then collected,
683 // sorted by geometry level (starting from ALIC) and
684 // then applied to the TGeo geometry.
685 // Finally an overlaps check is performed.
687 // Load alignment data from CDB and fill fAlignObjArray
688 if(fLoadAlignFromCDB){
690 TString detStr = detectors;
691 TString loadAlObjsListOfDets = "";
693 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
694 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
695 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
697 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
699 loadAlObjsListOfDets += fgkDetectorName[iDet];
700 loadAlObjsListOfDets += " ";
702 } // end loop over detectors
704 if(AliGeomManager::GetNalignable("GRP") != 0)
705 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
706 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
707 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
709 // Check if the array with alignment objects was
710 // provided by the user. If yes, apply the objects
711 // to the present TGeo geometry
712 if (fAlignObjArray) {
713 if (gGeoManager && gGeoManager->IsClosed()) {
714 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
715 AliError("The misalignment of one or more volumes failed!"
716 "Compare the list of simulated detectors and the list of detector alignment data!");
721 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
727 if (fAlignObjArray) {
728 fAlignObjArray->Delete();
729 delete fAlignObjArray; fAlignObjArray=NULL;
735 //_____________________________________________________________________________
736 void AliReconstruction::SetGAliceFile(const char* fileName)
738 // set the name of the galice file
740 fGAliceFileName = fileName;
743 //_____________________________________________________________________________
744 void AliReconstruction::SetInput(const char* input)
746 // In case the input string starts with 'mem://', we run in an online mode
747 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
748 // file is assumed. One can give as an input:
749 // mem://: - events taken from DAQ monitoring libs online
751 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
752 if (input) fRawInput = input;
755 //_____________________________________________________________________________
756 void AliReconstruction::SetOption(const char* detector, const char* option)
758 // set options for the reconstruction of a detector
760 TObject* obj = fOptions.FindObject(detector);
761 if (obj) fOptions.Remove(obj);
762 fOptions.Add(new TNamed(detector, option));
765 //_____________________________________________________________________________
766 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
768 // Set custom reconstruction parameters for a given detector
769 // Single set of parameters for all the events
771 // First check if the reco-params are global
772 if(!strcmp(detector, "GRP")) {
774 fRecoParam.AddDetRecoParam(kNDetectors,par);
778 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
779 if(!strcmp(detector, fgkDetectorName[iDet])) {
781 fRecoParam.AddDetRecoParam(iDet,par);
788 //_____________________________________________________________________________
789 Bool_t AliReconstruction::SetFieldMap(Float_t l3Cur, Float_t diCur, Float_t l3Pol,
790 Float_t diPol, Float_t beamenergy,
791 const Char_t *beamtype, const Char_t *path)
793 //------------------------------------------------
794 // The magnetic field map, defined externally...
795 // L3 current 30000 A -> 0.5 T
796 // L3 current 12000 A -> 0.2 T
797 // dipole current 6000 A
798 // The polarities must be the same
799 //------------------------------------------------
800 const Float_t l3NominalCurrent1=30000.; // (A)
801 const Float_t l3NominalCurrent2=12000.; // (A)
802 const Float_t diNominalCurrent =6000. ; // (A)
804 const Float_t tolerance=0.03; // relative current tolerance
805 const Float_t zero=77.; // "zero" current (A)
807 TString s=(l3Pol < 0) ? "L3: -" : "L3: +";
809 AliMagF::BMap_t map = AliMagF::k5kG;
813 l3Cur = TMath::Abs(l3Cur);
814 if (TMath::Abs(l3Cur-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
815 fcL3 = l3Cur/l3NominalCurrent1;
818 } else if (TMath::Abs(l3Cur-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
819 fcL3 = l3Cur/l3NominalCurrent2;
822 } else if (l3Cur <= zero) {
824 map = AliMagF::k5kGUniform;
826 fUniformField=kTRUE; // track with the uniform (zero) B field
828 AliError(Form("Wrong L3 current (%f A)!",l3Cur));
832 diCur = TMath::Abs(diCur);
833 if (TMath::Abs(diCur-diNominalCurrent)/diNominalCurrent < tolerance) {
834 // 3% current tolerance...
835 fcDip = diCur/diNominalCurrent;
837 } else if (diCur <= zero) { // some small current..
841 AliError(Form("Wrong dipole current (%f A)!",diCur));
845 if (l3Pol!=diPol && (map==AliMagF::k5kG || map==AliMagF::k2kG) && fcDip!=0) {
846 AliError("L3 and Dipole polarities must be the same");
850 if (l3Pol<0) fcL3 = -fcL3;
851 if (diPol<0) fcDip = -fcDip;
853 AliMagF::BeamType_t btype = AliMagF::kNoBeamField;
854 TString btypestr = beamtype;
856 TPRegexp protonBeam("(proton|p)\\s*-?\\s*\\1");
857 TPRegexp ionBeam("(lead|pb|ion|a)\\s*-?\\s*\\1");
858 if (btypestr.Contains(ionBeam)) btype = AliMagF::kBeamTypeAA;
859 else if (btypestr.Contains(protonBeam)) btype = AliMagF::kBeamTypepp;
861 AliInfo(Form("Cannot determine the beam type from %s, assume no LHC magnet field",beamtype));
864 AliMagF* fld = new AliMagF("MagneticFieldMap", s.Data(), 2, fcL3, fcDip, 10., map, path,
866 TGeoGlobalMagField::Instance()->SetField( fld );
867 TGeoGlobalMagField::Instance()->Lock();
873 Bool_t AliReconstruction::InitGRP() {
874 //------------------------------------
875 // Initialization of the GRP entry
876 //------------------------------------
877 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
881 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
884 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
886 fGRPData = new AliGRPObject();
887 fGRPData->ReadValuesFromMap(m);
891 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
892 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
896 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
900 AliError("No GRP entry found in OCDB!");
904 TString lhcState = fGRPData->GetLHCState();
905 if (lhcState==AliGRPObject::GetInvalidString()) {
906 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
907 lhcState = "UNKNOWN";
910 TString beamType = fGRPData->GetBeamType();
911 if (beamType==AliGRPObject::GetInvalidString()) {
912 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
913 beamType = "UNKNOWN";
916 Float_t beamEnergy = fGRPData->GetBeamEnergy();
917 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
918 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
921 // energy is provided in MeV*120
924 TString runType = fGRPData->GetRunType();
925 if (runType==AliGRPObject::GetInvalidString()) {
926 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
930 Int_t activeDetectors = fGRPData->GetDetectorMask();
931 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
932 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
933 activeDetectors = 1074790399;
936 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
940 // Process the list of active detectors
941 if (activeDetectors) {
942 UInt_t detMask = activeDetectors;
943 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
944 fRunTracking = MatchDetectorList(fRunTracking,detMask);
945 fFillESD = MatchDetectorList(fFillESD,detMask);
946 fQADetectors = MatchDetectorList(fQADetectors,detMask);
947 fLoadCDB.Form("%s %s %s %s",
948 fRunLocalReconstruction.Data(),
951 fQADetectors.Data());
952 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
953 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1)) {
954 // switch off the vertexer
955 AliInfo("SPD is not in the list of active detectors. Vertexer switched off.");
956 fRunVertexFinder = kFALSE;
958 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
959 // switch off the reading of CTP raw-data payload
960 if (fFillTriggerESD) {
961 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
962 fFillTriggerESD = kFALSE;
967 AliInfo("===================================================================================");
968 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
969 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
970 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
971 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
972 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
973 AliInfo("===================================================================================");
975 //*** Dealing with the magnetic field map
976 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {AliInfo("Running with the externally locked B field !");}
978 // Construct the field map out of the information retrieved from GRP.
981 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
982 if (l3Current == AliGRPObject::GetInvalidFloat()) {
983 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
987 Char_t l3Polarity = fGRPData->GetL3Polarity();
988 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
989 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
994 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
995 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
996 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1000 Char_t diPolarity = fGRPData->GetDipolePolarity();
1001 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1002 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1007 TObjString *l3Current=
1008 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
1010 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1013 TObjString *l3Polarity=
1014 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
1016 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1021 TObjString *diCurrent=
1022 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
1024 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1027 TObjString *diPolarity=
1028 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
1030 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1036 if ( !SetFieldMap(l3Current, diCurrent, l3Polarity ? -1:1, diPolarity ? -1:1) )
1037 AliFatal("Failed to creat a B field map ! Exiting...");
1038 AliInfo("Running with the B field constructed out of GRP !");
1040 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1044 //*** Get the diamond profiles from OCDB
1045 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1047 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1049 AliError("No SPD diamond profile found in OCDB!");
1052 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1054 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1056 AliError("No diamond profile found in OCDB!");
1059 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1061 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1063 AliError("No TPC diamond profile found in OCDB!");
1069 //_____________________________________________________________________________
1070 Bool_t AliReconstruction::LoadCDB()
1072 AliCodeTimerAuto("");
1074 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1076 TString detStr = fLoadCDB;
1077 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1078 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1079 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1084 //_____________________________________________________________________________
1085 Bool_t AliReconstruction::Run(const char* input)
1088 AliCodeTimerAuto("");
1091 if (GetAbort() != TSelector::kContinue) return kFALSE;
1093 TChain *chain = NULL;
1094 if (fRawReader && (chain = fRawReader->GetChain())) {
1097 gProof->AddInput(this);
1099 outputFile.SetProtocol("root",kTRUE);
1100 outputFile.SetHost(gSystem->HostName());
1101 outputFile.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1102 AliInfo(Form("Output file with ESDs is %s",outputFile.GetUrl()));
1103 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.GetUrl()));
1105 chain->Process("AliReconstruction");
1108 chain->Process(this);
1113 if (GetAbort() != TSelector::kContinue) return kFALSE;
1115 if (GetAbort() != TSelector::kContinue) return kFALSE;
1116 //******* The loop over events
1117 AliInfo("Starting looping over events");
1119 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1120 (fRawReader && fRawReader->NextEvent())) {
1121 if (!ProcessEvent(iEvent)) {
1122 Abort("ProcessEvent",TSelector::kAbortFile);
1128 if (GetAbort() != TSelector::kContinue) return kFALSE;
1130 if (GetAbort() != TSelector::kContinue) return kFALSE;
1136 //_____________________________________________________________________________
1137 void AliReconstruction::InitRawReader(const char* input)
1139 AliCodeTimerAuto("");
1141 // Init raw-reader and
1142 // set the input in case of raw data
1143 if (input) fRawInput = input;
1144 fRawReader = AliRawReader::Create(fRawInput.Data());
1146 AliInfo("Reconstruction will run over digits");
1148 if (!fEquipIdMap.IsNull() && fRawReader)
1149 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1151 if (!fUseHLTData.IsNull()) {
1152 // create the RawReaderHLT which performs redirection of HLT input data for
1153 // the specified detectors
1154 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1156 fParentRawReader=fRawReader;
1157 fRawReader=pRawReader;
1159 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1162 AliSysInfo::AddStamp("CreateRawReader");
1165 //_____________________________________________________________________________
1166 void AliReconstruction::InitRun(const char* input)
1168 // Initialization of raw-reader,
1169 // run number, CDB etc.
1170 AliCodeTimerAuto("");
1171 AliSysInfo::AddStamp("Start");
1173 // Initialize raw-reader if any
1174 InitRawReader(input);
1176 // Initialize the CDB storage
1179 // Set run number in CDBManager (if it is not already set by the user)
1180 if (!SetRunNumberFromData()) {
1181 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1185 // Set CDB lock: from now on it is forbidden to reset the run number
1186 // or the default storage or to activate any further storage!
1191 //_____________________________________________________________________________
1192 void AliReconstruction::Begin(TTree *)
1194 // Initialize AlReconstruction before
1195 // going into the event loop
1196 // Should follow the TSelector convention
1197 // i.e. initialize only the object on the client side
1198 AliCodeTimerAuto("");
1200 AliReconstruction *reco = NULL;
1202 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1205 AliSysInfo::AddStamp("ReadInputInBegin");
1208 // Import ideal TGeo geometry and apply misalignment
1210 TString geom(gSystem->DirName(fGAliceFileName));
1211 geom += "/geometry.root";
1212 AliGeomManager::LoadGeometry(geom.Data());
1214 Abort("LoadGeometry", TSelector::kAbortProcess);
1217 AliSysInfo::AddStamp("LoadGeom");
1218 TString detsToCheck=fRunLocalReconstruction;
1219 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1220 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1223 AliSysInfo::AddStamp("CheckGeom");
1226 if (!MisalignGeometry(fLoadAlignData)) {
1227 Abort("MisalignGeometry", TSelector::kAbortProcess);
1230 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1231 AliSysInfo::AddStamp("MisalignGeom");
1234 Abort("InitGRP", TSelector::kAbortProcess);
1237 AliSysInfo::AddStamp("InitGRP");
1240 Abort("LoadCDB", TSelector::kAbortProcess);
1243 AliSysInfo::AddStamp("LoadCDB");
1245 // Read the reconstruction parameters from OCDB
1246 if (!InitRecoParams()) {
1247 AliWarning("Not all detectors have correct RecoParam objects initialized");
1249 AliSysInfo::AddStamp("InitRecoParams");
1251 if (fInput && gProof) {
1252 if (reco) *reco = *this;
1254 gProof->AddInputData(gGeoManager,kTRUE);
1256 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1257 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1258 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1259 magFieldMap->SetName("MagneticFieldMap");
1260 gProof->AddInputData(magFieldMap,kTRUE);
1265 //_____________________________________________________________________________
1266 void AliReconstruction::SlaveBegin(TTree*)
1268 // Initialization related to run-loader,
1269 // vertexer, trackers, recontructors
1270 // In proof mode it is executed on the slave
1271 AliCodeTimerAuto("");
1273 TProofOutputFile *outProofFile = NULL;
1275 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1278 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1280 AliGeomManager::SetGeometry(tgeo);
1282 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1283 Int_t runNumber = -1;
1284 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1285 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1286 man->SetCacheFlag(kTRUE);
1287 man->SetLock(kTRUE);
1291 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1292 TGeoGlobalMagField::Instance()->SetField(map);
1294 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1295 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1296 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1297 fOutput->Add(outProofFile);
1299 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1302 // get the run loader
1303 if (!InitRunLoader()) {
1304 Abort("InitRunLoader", TSelector::kAbortProcess);
1307 AliSysInfo::AddStamp("LoadLoader");
1309 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1312 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1313 Abort("CreateTrackers", TSelector::kAbortProcess);
1316 AliSysInfo::AddStamp("CreateTrackers");
1318 // create the ESD output file and tree
1319 if (!outProofFile) {
1320 ffile = TFile::Open("AliESDs.root", "RECREATE");
1321 ffile->SetCompressionLevel(2);
1322 if (!ffile->IsOpen()) {
1323 Abort("OpenESDFile", TSelector::kAbortProcess);
1328 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1329 Abort(Form("Problems opening output PROOF file: %s/%s",
1330 outProofFile->GetDir(), outProofFile->GetFileName()),
1331 TSelector::kAbortProcess);
1336 ftree = new TTree("esdTree", "Tree with ESD objects");
1337 fesd = new AliESDEvent();
1338 fesd->CreateStdContent();
1340 fesd->WriteToTree(ftree);
1341 if (fWriteESDfriend) {
1343 // Since we add the branch manually we must
1344 // book and add it after WriteToTree
1345 // otherwise it is created twice,
1346 // once via writetotree and once here.
1347 // The case for AliESDfriend is now
1348 // caught also in AlIESDEvent::WriteToTree but
1349 // be careful when changing the name (AliESDfriend is not
1350 // a TNamed so we had to hardwire it)
1351 fesdf = new AliESDfriend();
1352 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1353 br->SetFile("AliESDfriends.root");
1354 fesd->AddObject(fesdf);
1356 ftree->GetUserInfo()->Add(fesd);
1358 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1359 fhltesd = new AliESDEvent();
1360 fhltesd->CreateStdContent();
1362 // read the ESD template from CDB
1363 // HLT is allowed to put non-std content to its ESD, the non-std
1364 // objects need to be created before invocation of WriteToTree in
1365 // order to create all branches. Initialization is done from an
1366 // ESD layout template in CDB
1367 AliCDBManager* man = AliCDBManager::Instance();
1368 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1369 AliCDBEntry* hltESDConfig=NULL;
1370 if (man->GetId(hltESDConfigPath)!=NULL &&
1371 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1372 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1374 // init all internal variables from the list of objects
1375 pESDLayout->GetStdContent();
1377 // copy content and create non-std objects
1378 *fhltesd=*pESDLayout;
1381 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1382 hltESDConfigPath.GetPath().Data()));
1386 fhltesd->WriteToTree(fhlttree);
1387 fhlttree->GetUserInfo()->Add(fhltesd);
1389 ProcInfo_t procInfo;
1390 gSystem->GetProcInfo(&procInfo);
1391 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1394 //Initialize the QA and start of cycle
1396 fQAManager = AliQAManager::QAManager("rec") ;
1397 fQAManager->SetActiveDetectors(fQADetectors) ;
1398 for (Int_t det = 0 ; det < AliQA::kNDET ; det++) {
1399 fQAManager->SetCycleLength(AliQA::DETECTORINDEX_t(det), fQACycles[det]) ;
1400 fQAManager->SetWriteExpert(AliQA::DETECTORINDEX_t(det)) ;
1402 if (!fRawReader && fQATasks.Contains(AliQA::kRAWS))
1403 fQATasks.ReplaceAll(Form("%d",AliQA::kRAWS), "") ;
1404 fQAManager->SetTasks(fQATasks) ;
1405 fQAManager->InitQADataMaker(AliCDBManager::Instance()->GetRun()) ;
1409 Bool_t sameCycle = kFALSE ;
1411 fQAManager = AliQAManager::QAManager("rec") ;
1412 AliQADataMaker *qadm = fQAManager->GetQADataMaker(AliQA::kGLOBAL);
1413 AliInfo(Form("Initializing the global QA data maker"));
1414 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1415 qadm->StartOfCycle(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1416 TObjArray **arr=qadm->Init(AliQA::kRECPOINTS);
1417 AliTracker::SetResidualsArray(arr);
1420 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1421 qadm->StartOfCycle(AliQA::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1422 qadm->Init(AliQA::kESDS);
1426 //Initialize the Plane Efficiency framework
1427 if (fRunPlaneEff && !InitPlaneEff()) {
1428 Abort("InitPlaneEff", TSelector::kAbortProcess);
1432 if (strcmp(gProgName,"alieve") == 0)
1433 fRunAliEVE = InitAliEVE();
1438 //_____________________________________________________________________________
1439 Bool_t AliReconstruction::Process(Long64_t entry)
1441 // run the reconstruction over a single entry
1442 // from the chain with raw data
1443 AliCodeTimerAuto("");
1445 TTree *currTree = fChain->GetTree();
1446 AliRawEvent *event = new AliRawEvent;
1447 currTree->SetBranchAddress("rawevent",&event);
1448 currTree->GetEntry(entry);
1449 fRawReader = new AliRawReaderRoot(event);
1450 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1458 //_____________________________________________________________________________
1459 void AliReconstruction::Init(TTree *tree)
1462 AliError("The input tree is not found!");
1468 //_____________________________________________________________________________
1469 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1471 // run the reconstruction over a single event
1472 // The event loop is steered in Run method
1474 AliCodeTimerAuto("");
1476 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1477 fRunLoader->SetEventNumber(iEvent);
1478 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1480 fRunLoader->TreeE()->Fill();
1481 if (fRawReader && fRawReader->UseAutoSaveESD())
1482 fRunLoader->TreeE()->AutoSave("SaveSelf");
1485 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1489 AliInfo(Form("processing event %d", iEvent));
1491 fRunLoader->GetEvent(iEvent);
1493 // Fill Event-info object
1495 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo);
1496 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1498 // Set the reco-params
1500 TString detStr = fLoadCDB;
1501 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1502 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1503 AliReconstructor *reconstructor = GetReconstructor(iDet);
1504 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1505 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1506 reconstructor->SetRecoParam(par);
1508 fQAManager->SetRecoParam(iDet, par) ;
1516 fQAManager->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1517 fQAManager->RunOneEvent(fRawReader) ;
1519 // local single event reconstruction
1520 if (!fRunLocalReconstruction.IsNull()) {
1521 TString detectors=fRunLocalReconstruction;
1522 // run HLT event reconstruction first
1523 // ;-( IsSelected changes the string
1524 if (IsSelected("HLT", detectors) &&
1525 !RunLocalEventReconstruction("HLT")) {
1526 if (fStopOnError) {CleanUp(); return kFALSE;}
1528 detectors=fRunLocalReconstruction;
1529 detectors.ReplaceAll("HLT", "");
1530 if (!RunLocalEventReconstruction(detectors)) {
1531 if (fStopOnError) {CleanUp(); return kFALSE;}
1535 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1536 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1537 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1538 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1540 // Set magnetic field from the tracker
1541 fesd->SetMagneticField(AliTracker::GetBz());
1542 fhltesd->SetMagneticField(AliTracker::GetBz());
1544 // Set most probable pt, for B=0 tracking
1545 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1546 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1547 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1549 // Fill raw-data error log into the ESD
1550 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1553 if (fRunVertexFinder) {
1554 if (!RunVertexFinder(fesd)) {
1555 if (fStopOnError) {CleanUp(); return kFALSE;}
1560 if (!fRunTracking.IsNull()) {
1561 if (fRunMuonTracking) {
1562 if (!RunMuonTracking(fesd)) {
1563 if (fStopOnError) {CleanUp(); return kFALSE;}
1569 if (!fRunTracking.IsNull()) {
1570 if (!RunTracking(fesd)) {
1571 if (fStopOnError) {CleanUp(); return kFALSE;}
1576 if (!fFillESD.IsNull()) {
1577 TString detectors=fFillESD;
1578 // run HLT first and on hltesd
1579 // ;-( IsSelected changes the string
1580 if (IsSelected("HLT", detectors) &&
1581 !FillESD(fhltesd, "HLT")) {
1582 if (fStopOnError) {CleanUp(); return kFALSE;}
1585 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1586 if (detectors.Contains("ALL")) {
1588 for (Int_t idet=0; idet<kNDetectors; ++idet){
1589 detectors += fgkDetectorName[idet];
1593 detectors.ReplaceAll("HLT", "");
1594 if (!FillESD(fesd, detectors)) {
1595 if (fStopOnError) {CleanUp(); return kFALSE;}
1599 // fill Event header information from the RawEventHeader
1600 if (fRawReader){FillRawEventHeaderESD(fesd);}
1603 AliESDpid::MakePID(fesd);
1605 if (fFillTriggerESD) {
1606 if (!FillTriggerESD(fesd)) {
1607 if (fStopOnError) {CleanUp(); return kFALSE;}
1614 // Propagate track to the beam pipe (if not already done by ITS)
1616 const Int_t ntracks = fesd->GetNumberOfTracks();
1617 const Double_t kBz = fesd->GetMagneticField();
1618 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1621 UShort_t *selectedIdx=new UShort_t[ntracks];
1623 for (Int_t itrack=0; itrack<ntracks; itrack++){
1624 const Double_t kMaxStep = 1; //max step over the material
1627 AliESDtrack *track = fesd->GetTrack(itrack);
1628 if (!track) continue;
1630 AliExternalTrackParam *tpcTrack =
1631 (AliExternalTrackParam *)track->GetTPCInnerParam();
1635 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1638 Int_t n=trkArray.GetEntriesFast();
1639 selectedIdx[n]=track->GetID();
1640 trkArray.AddLast(tpcTrack);
1643 //Tracks refitted by ITS should already be at the SPD vertex
1644 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1647 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1648 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1653 // Improve the reconstructed primary vertex position using the tracks
1655 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1656 if(fesd->GetPrimaryVertexSPD()) {
1657 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1658 if(vtitle.Contains("cosmics")) {
1659 runVertexFinderTracks=kFALSE;
1663 if (runVertexFinderTracks) {
1664 // TPC + ITS primary vertex
1665 ftVertexer->SetITSMode();
1666 ftVertexer->SetConstraintOff();
1667 // get cuts for vertexer from AliGRPRecoParam
1669 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1670 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1671 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1672 ftVertexer->SetCuts(cutsVertexer);
1673 delete [] cutsVertexer; cutsVertexer = NULL;
1674 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1675 ftVertexer->SetVtxStart(fDiamondProfile);
1677 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1679 if (pvtx->GetStatus()) {
1680 fesd->SetPrimaryVertexTracks(pvtx);
1681 for (Int_t i=0; i<ntracks; i++) {
1682 AliESDtrack *t = fesd->GetTrack(i);
1683 t->RelateToVertex(pvtx, kBz, kVeryBig);
1688 // TPC-only primary vertex
1689 ftVertexer->SetTPCMode();
1690 ftVertexer->SetConstraintOff();
1691 // get cuts for vertexer from AliGRPRecoParam
1693 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1694 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1695 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1696 ftVertexer->SetCuts(cutsVertexer);
1697 delete [] cutsVertexer; cutsVertexer = NULL;
1698 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1699 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1701 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1703 if (pvtx->GetStatus()) {
1704 fesd->SetPrimaryVertexTPC(pvtx);
1705 for (Int_t i=0; i<ntracks; i++) {
1706 AliESDtrack *t = fesd->GetTrack(i);
1707 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1713 delete[] selectedIdx;
1715 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1720 AliV0vertexer vtxer;
1721 vtxer.Tracks2V0vertices(fesd);
1723 if (fRunCascadeFinder) {
1725 AliCascadeVertexer cvtxer;
1726 cvtxer.V0sTracks2CascadeVertices(fesd);
1731 if (fCleanESD) CleanESD(fesd);
1734 fQAManager->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1735 fQAManager->RunOneEvent(fesd) ;
1738 AliQADataMaker *qadm = fQAManager->GetQADataMaker(AliQA::kGLOBAL);
1739 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1740 if (qadm && fQATasks.Contains(Form("%d", AliQA::kESDS)))
1741 qadm->Exec(AliQA::kESDS, fesd);
1744 if (fWriteESDfriend) {
1745 // fesdf->~AliESDfriend();
1746 // new (fesdf) AliESDfriend(); // Reset...
1747 fesd->GetESDfriend(fesdf);
1751 // Auto-save the ESD tree in case of prompt reco @P2
1752 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1753 ftree->AutoSave("SaveSelf");
1754 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1755 if (friendfile) friendfile->Save();
1762 if (fRunAliEVE) RunAliEVE();
1766 if (fWriteESDfriend) {
1767 fesdf->~AliESDfriend();
1768 new (fesdf) AliESDfriend(); // Reset...
1771 ProcInfo_t procInfo;
1772 gSystem->GetProcInfo(&procInfo);
1773 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1776 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1777 if (fReconstructor[iDet])
1778 fReconstructor[iDet]->SetRecoParam(NULL);
1781 if (fRunQA || fRunGlobalQA)
1782 fQAManager->Increment() ;
1787 //_____________________________________________________________________________
1788 void AliReconstruction::SlaveTerminate()
1790 // Finalize the run on the slave side
1791 // Called after the exit
1792 // from the event loop
1793 AliCodeTimerAuto("");
1795 if (fIsNewRunLoader) { // galice.root didn't exist
1796 fRunLoader->WriteHeader("OVERWRITE");
1797 fRunLoader->CdGAFile();
1798 fRunLoader->Write(0, TObject::kOverwrite);
1801 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1802 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1804 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1805 cdbMapCopy->SetOwner(1);
1806 cdbMapCopy->SetName("cdbMap");
1807 TIter iter(cdbMap->GetTable());
1810 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1811 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1812 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1813 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1816 TList *cdbListCopy = new TList();
1817 cdbListCopy->SetOwner(1);
1818 cdbListCopy->SetName("cdbList");
1820 TIter iter2(cdbList);
1823 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1824 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1827 ftree->GetUserInfo()->Add(cdbMapCopy);
1828 ftree->GetUserInfo()->Add(cdbListCopy);
1833 if (fWriteESDfriend)
1834 ftree->SetBranchStatus("ESDfriend*",0);
1835 // we want to have only one tree version number
1836 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1839 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1840 if (fRunPlaneEff && !FinishPlaneEff()) {
1841 AliWarning("Finish PlaneEff evaluation failed");
1844 // End of cycle for the in-loop
1846 fQAManager->EndOfCycle() ;
1849 AliQADataMaker *qadm = fQAManager->GetQADataMaker(AliQA::kGLOBAL);
1851 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1852 qadm->EndOfCycle(AliQA::kRECPOINTS);
1853 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1854 qadm->EndOfCycle(AliQA::kESDS);
1862 //_____________________________________________________________________________
1863 void AliReconstruction::Terminate()
1865 // Create tags for the events in the ESD tree (the ESD tree is always present)
1866 // In case of empty events the tags will contain dummy values
1867 AliCodeTimerAuto("");
1869 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1870 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQA::Instance()->GetQA(), AliQA::Instance()->GetEventSpecies(), AliQA::kNDET, AliRecoParam::kNSpecies);
1872 // Cleanup of CDB manager: cache and active storages!
1873 AliCDBManager::Instance()->ClearCache();
1876 //_____________________________________________________________________________
1877 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1879 // run the local reconstruction
1881 static Int_t eventNr=0;
1882 AliCodeTimerAuto("")
1884 TString detStr = detectors;
1885 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1886 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1887 AliReconstructor* reconstructor = GetReconstructor(iDet);
1888 if (!reconstructor) continue;
1889 AliLoader* loader = fLoader[iDet];
1890 // Matthias April 2008: temporary fix to run HLT reconstruction
1891 // although the HLT loader is missing
1892 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1894 reconstructor->Reconstruct(fRawReader, NULL);
1897 reconstructor->Reconstruct(dummy, NULL);
1902 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1905 // conversion of digits
1906 if (fRawReader && reconstructor->HasDigitConversion()) {
1907 AliInfo(Form("converting raw data digits into root objects for %s",
1908 fgkDetectorName[iDet]));
1909 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1910 // fgkDetectorName[iDet]));
1911 loader->LoadDigits("update");
1912 loader->CleanDigits();
1913 loader->MakeDigitsContainer();
1914 TTree* digitsTree = loader->TreeD();
1915 reconstructor->ConvertDigits(fRawReader, digitsTree);
1916 loader->WriteDigits("OVERWRITE");
1917 loader->UnloadDigits();
1919 // local reconstruction
1920 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1921 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1922 loader->LoadRecPoints("update");
1923 loader->CleanRecPoints();
1924 loader->MakeRecPointsContainer();
1925 TTree* clustersTree = loader->TreeR();
1926 if (fRawReader && !reconstructor->HasDigitConversion()) {
1927 reconstructor->Reconstruct(fRawReader, clustersTree);
1929 loader->LoadDigits("read");
1930 TTree* digitsTree = loader->TreeD();
1932 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1933 if (fStopOnError) return kFALSE;
1935 reconstructor->Reconstruct(digitsTree, clustersTree);
1937 loader->UnloadDigits();
1940 TString detQAStr(fQADetectors) ;
1942 fQAManager->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1943 fQAManager->RunOneEventInOneDetector(iDet, clustersTree) ;
1945 loader->WriteRecPoints("OVERWRITE");
1946 loader->UnloadRecPoints();
1947 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1949 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1950 AliError(Form("the following detectors were not found: %s",
1952 if (fStopOnError) return kFALSE;
1958 //_____________________________________________________________________________
1959 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1961 // run the barrel tracking
1963 AliCodeTimerAuto("")
1965 AliVertexer *vertexer = CreateVertexer();
1966 if (!vertexer) return kFALSE;
1968 AliInfo("running the ITS vertex finder");
1969 AliESDVertex* vertex = NULL;
1971 fLoader[0]->LoadRecPoints();
1972 TTree* cltree = fLoader[0]->TreeR();
1974 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
1975 vertex = vertexer->FindVertexForCurrentEvent(cltree);
1978 AliError("Can't get the ITS cluster tree");
1980 fLoader[0]->UnloadRecPoints();
1983 AliError("Can't get the ITS loader");
1986 AliWarning("Vertex not found");
1987 vertex = new AliESDVertex();
1988 vertex->SetName("default");
1991 vertex->SetName("reconstructed");
1996 vertex->GetXYZ(vtxPos);
1997 vertex->GetSigmaXYZ(vtxErr);
1999 esd->SetPrimaryVertexSPD(vertex);
2000 // if SPD multiplicity has been determined, it is stored in the ESD
2001 AliMultiplicity *mult = vertexer->GetMultiplicity();
2002 if(mult)esd->SetMultiplicity(mult);
2004 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2005 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2014 //_____________________________________________________________________________
2015 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2017 // run the HLT barrel tracking
2019 AliCodeTimerAuto("")
2022 AliError("Missing runLoader!");
2026 AliInfo("running HLT tracking");
2028 // Get a pointer to the HLT reconstructor
2029 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2030 if (!reconstructor) return kFALSE;
2033 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2034 TString detName = fgkDetectorName[iDet];
2035 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2036 reconstructor->SetOption(detName.Data());
2037 AliTracker *tracker = reconstructor->CreateTracker();
2039 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2040 if (fStopOnError) return kFALSE;
2044 Double_t vtxErr[3]={0.005,0.005,0.010};
2045 const AliESDVertex *vertex = esd->GetVertex();
2046 vertex->GetXYZ(vtxPos);
2047 tracker->SetVertex(vtxPos,vtxErr);
2049 fLoader[iDet]->LoadRecPoints("read");
2050 TTree* tree = fLoader[iDet]->TreeR();
2052 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2055 tracker->LoadClusters(tree);
2057 if (tracker->Clusters2Tracks(esd) != 0) {
2058 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2062 tracker->UnloadClusters();
2070 //_____________________________________________________________________________
2071 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2073 // run the muon spectrometer tracking
2075 AliCodeTimerAuto("")
2078 AliError("Missing runLoader!");
2081 Int_t iDet = 7; // for MUON
2083 AliInfo("is running...");
2085 // Get a pointer to the MUON reconstructor
2086 AliReconstructor *reconstructor = GetReconstructor(iDet);
2087 if (!reconstructor) return kFALSE;
2090 TString detName = fgkDetectorName[iDet];
2091 AliDebug(1, Form("%s tracking", detName.Data()));
2092 AliTracker *tracker = reconstructor->CreateTracker();
2094 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2099 fLoader[iDet]->LoadRecPoints("read");
2101 tracker->LoadClusters(fLoader[iDet]->TreeR());
2103 Int_t rv = tracker->Clusters2Tracks(esd);
2107 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2111 fLoader[iDet]->UnloadRecPoints();
2113 tracker->UnloadClusters();
2121 //_____________________________________________________________________________
2122 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2124 // run the barrel tracking
2125 static Int_t eventNr=0;
2126 AliCodeTimerAuto("")
2128 AliInfo("running tracking");
2130 //Fill the ESD with the T0 info (will be used by the TOF)
2131 if (fReconstructor[11] && fLoader[11]) {
2132 fLoader[11]->LoadRecPoints("READ");
2133 TTree *treeR = fLoader[11]->TreeR();
2135 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2139 // pass 1: TPC + ITS inwards
2140 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2141 if (!fTracker[iDet]) continue;
2142 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2145 fLoader[iDet]->LoadRecPoints("read");
2146 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2147 TTree* tree = fLoader[iDet]->TreeR();
2149 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2152 fTracker[iDet]->LoadClusters(tree);
2153 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2155 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2156 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2159 // preliminary PID in TPC needed by the ITS tracker
2161 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2162 AliESDpid::MakePID(esd);
2164 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2167 // pass 2: ALL backwards
2169 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2170 if (!fTracker[iDet]) continue;
2171 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2174 if (iDet > 1) { // all except ITS, TPC
2176 fLoader[iDet]->LoadRecPoints("read");
2177 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2178 tree = fLoader[iDet]->TreeR();
2180 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2183 fTracker[iDet]->LoadClusters(tree);
2184 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2188 if (iDet>1) // start filling residuals for the "outer" detectors
2189 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2191 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2192 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2197 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2198 fTracker[iDet]->UnloadClusters();
2199 fLoader[iDet]->UnloadRecPoints();
2201 // updated PID in TPC needed by the ITS tracker -MI
2203 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2204 AliESDpid::MakePID(esd);
2206 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2208 //stop filling residuals for the "outer" detectors
2209 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2211 // pass 3: TRD + TPC + ITS refit inwards
2213 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2214 if (!fTracker[iDet]) continue;
2215 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2218 if (iDet<2) // start filling residuals for TPC and ITS
2219 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2221 if (fTracker[iDet]->RefitInward(esd) != 0) {
2222 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2225 // run postprocessing
2226 if (fTracker[iDet]->PostProcess(esd) != 0) {
2227 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2230 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2233 // write space-points to the ESD in case alignment data output
2235 if (fWriteAlignmentData)
2236 WriteAlignmentData(esd);
2238 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2239 if (!fTracker[iDet]) continue;
2241 fTracker[iDet]->UnloadClusters();
2242 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2243 fLoader[iDet]->UnloadRecPoints();
2244 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2246 // stop filling residuals for TPC and ITS
2247 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2253 //_____________________________________________________________________________
2254 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2256 // Remove the data which are not needed for the physics analysis.
2259 Int_t nTracks=esd->GetNumberOfTracks();
2260 Int_t nV0s=esd->GetNumberOfV0s();
2262 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2264 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2265 Bool_t rc=esd->Clean(cleanPars);
2267 nTracks=esd->GetNumberOfTracks();
2268 nV0s=esd->GetNumberOfV0s();
2270 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2275 //_____________________________________________________________________________
2276 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2278 // fill the event summary data
2280 AliCodeTimerAuto("")
2281 static Int_t eventNr=0;
2282 TString detStr = detectors;
2284 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2285 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2286 AliReconstructor* reconstructor = GetReconstructor(iDet);
2287 if (!reconstructor) continue;
2288 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2289 TTree* clustersTree = NULL;
2290 if (fLoader[iDet]) {
2291 fLoader[iDet]->LoadRecPoints("read");
2292 clustersTree = fLoader[iDet]->TreeR();
2293 if (!clustersTree) {
2294 AliError(Form("Can't get the %s clusters tree",
2295 fgkDetectorName[iDet]));
2296 if (fStopOnError) return kFALSE;
2299 if (fRawReader && !reconstructor->HasDigitConversion()) {
2300 reconstructor->FillESD(fRawReader, clustersTree, esd);
2302 TTree* digitsTree = NULL;
2303 if (fLoader[iDet]) {
2304 fLoader[iDet]->LoadDigits("read");
2305 digitsTree = fLoader[iDet]->TreeD();
2307 AliError(Form("Can't get the %s digits tree",
2308 fgkDetectorName[iDet]));
2309 if (fStopOnError) return kFALSE;
2312 reconstructor->FillESD(digitsTree, clustersTree, esd);
2313 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2315 if (fLoader[iDet]) {
2316 fLoader[iDet]->UnloadRecPoints();
2320 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2321 AliError(Form("the following detectors were not found: %s",
2323 if (fStopOnError) return kFALSE;
2325 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2330 //_____________________________________________________________________________
2331 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2333 // Reads the trigger decision which is
2334 // stored in Trigger.root file and fills
2335 // the corresponding esd entries
2337 AliCodeTimerAuto("")
2339 AliInfo("Filling trigger information into the ESD");
2342 AliCTPRawStream input(fRawReader);
2343 if (!input.Next()) {
2344 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2347 if (esd->GetTriggerMask() != input.GetClassMask())
2348 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2349 input.GetClassMask(),esd->GetTriggerMask()));
2350 if (esd->GetOrbitNumber() != input.GetOrbitID())
2351 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2352 input.GetOrbitID(),esd->GetOrbitNumber()));
2353 if (esd->GetBunchCrossNumber() != input.GetBCID())
2354 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2355 input.GetBCID(),esd->GetBunchCrossNumber()));
2358 // Here one has to add the filling of trigger inputs and
2359 // interaction records
2369 //_____________________________________________________________________________
2370 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2373 // Filling information from RawReader Header
2376 if (!fRawReader) return kFALSE;
2378 AliInfo("Filling information from RawReader Header");
2380 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2381 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2382 esd->SetPeriodNumber(fRawReader->GetPeriod());
2384 esd->SetTimeStamp(fRawReader->GetTimestamp());
2385 esd->SetEventType(fRawReader->GetType());
2391 //_____________________________________________________________________________
2392 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2394 // check whether detName is contained in detectors
2395 // if yes, it is removed from detectors
2397 // check if all detectors are selected
2398 if ((detectors.CompareTo("ALL") == 0) ||
2399 detectors.BeginsWith("ALL ") ||
2400 detectors.EndsWith(" ALL") ||
2401 detectors.Contains(" ALL ")) {
2406 // search for the given detector
2407 Bool_t result = kFALSE;
2408 if ((detectors.CompareTo(detName) == 0) ||
2409 detectors.BeginsWith(detName+" ") ||
2410 detectors.EndsWith(" "+detName) ||
2411 detectors.Contains(" "+detName+" ")) {
2412 detectors.ReplaceAll(detName, "");
2416 // clean up the detectors string
2417 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2418 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2419 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2424 //_____________________________________________________________________________
2425 Bool_t AliReconstruction::InitRunLoader()
2427 // get or create the run loader
2429 if (gAlice) delete gAlice;
2432 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2433 // load all base libraries to get the loader classes
2434 TString libs = gSystem->GetLibraries();
2435 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2436 TString detName = fgkDetectorName[iDet];
2437 if (detName == "HLT") continue;
2438 if (libs.Contains("lib" + detName + "base.so")) continue;
2439 gSystem->Load("lib" + detName + "base.so");
2441 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2443 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2448 fRunLoader->CdGAFile();
2449 fRunLoader->LoadgAlice();
2451 //PH This is a temporary fix to give access to the kinematics
2452 //PH that is needed for the labels of ITS clusters
2453 fRunLoader->LoadHeader();
2454 fRunLoader->LoadKinematics();
2456 } else { // galice.root does not exist
2458 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2460 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2461 AliConfig::GetDefaultEventFolderName(),
2464 AliError(Form("could not create run loader in file %s",
2465 fGAliceFileName.Data()));
2469 fIsNewRunLoader = kTRUE;
2470 fRunLoader->MakeTree("E");
2472 if (fNumberOfEventsPerFile > 0)
2473 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2475 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2481 //_____________________________________________________________________________
2482 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2484 // get the reconstructor object and the loader for a detector
2486 if (fReconstructor[iDet]) {
2487 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2488 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2489 fReconstructor[iDet]->SetRecoParam(par);
2491 return fReconstructor[iDet];
2494 // load the reconstructor object
2495 TPluginManager* pluginManager = gROOT->GetPluginManager();
2496 TString detName = fgkDetectorName[iDet];
2497 TString recName = "Ali" + detName + "Reconstructor";
2499 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2501 AliReconstructor* reconstructor = NULL;
2502 // first check if a plugin is defined for the reconstructor
2503 TPluginHandler* pluginHandler =
2504 pluginManager->FindHandler("AliReconstructor", detName);
2505 // if not, add a plugin for it
2506 if (!pluginHandler) {
2507 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2508 TString libs = gSystem->GetLibraries();
2509 if (libs.Contains("lib" + detName + "base.so") ||
2510 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2511 pluginManager->AddHandler("AliReconstructor", detName,
2512 recName, detName + "rec", recName + "()");
2514 pluginManager->AddHandler("AliReconstructor", detName,
2515 recName, detName, recName + "()");
2517 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2519 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2520 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2522 if (reconstructor) {
2523 TObject* obj = fOptions.FindObject(detName.Data());
2524 if (obj) reconstructor->SetOption(obj->GetTitle());
2525 reconstructor->Init();
2526 fReconstructor[iDet] = reconstructor;
2529 // get or create the loader
2530 if (detName != "HLT") {
2531 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2532 if (!fLoader[iDet]) {
2533 AliConfig::Instance()
2534 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2536 // first check if a plugin is defined for the loader
2538 pluginManager->FindHandler("AliLoader", detName);
2539 // if not, add a plugin for it
2540 if (!pluginHandler) {
2541 TString loaderName = "Ali" + detName + "Loader";
2542 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2543 pluginManager->AddHandler("AliLoader", detName,
2544 loaderName, detName + "base",
2545 loaderName + "(const char*, TFolder*)");
2546 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2548 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2550 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2551 fRunLoader->GetEventFolder());
2553 if (!fLoader[iDet]) { // use default loader
2554 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2556 if (!fLoader[iDet]) {
2557 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2558 if (fStopOnError) return NULL;
2560 fRunLoader->AddLoader(fLoader[iDet]);
2561 fRunLoader->CdGAFile();
2562 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2563 fRunLoader->Write(0, TObject::kOverwrite);
2568 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2569 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2570 reconstructor->SetRecoParam(par);
2572 return reconstructor;
2575 //_____________________________________________________________________________
2576 AliVertexer* AliReconstruction::CreateVertexer()
2578 // create the vertexer
2579 // Please note that the caller is the owner of the
2582 AliVertexer* vertexer = NULL;
2583 AliReconstructor* itsReconstructor = GetReconstructor(0);
2584 if (itsReconstructor) {
2585 vertexer = itsReconstructor->CreateVertexer();
2588 AliWarning("couldn't create a vertexer for ITS");
2594 //_____________________________________________________________________________
2595 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2597 // create the trackers
2598 AliInfo("Creating trackers");
2600 TString detStr = detectors;
2601 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2602 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2603 AliReconstructor* reconstructor = GetReconstructor(iDet);
2604 if (!reconstructor) continue;
2605 TString detName = fgkDetectorName[iDet];
2606 if (detName == "HLT") {
2607 fRunHLTTracking = kTRUE;
2610 if (detName == "MUON") {
2611 fRunMuonTracking = kTRUE;
2616 fTracker[iDet] = reconstructor->CreateTracker();
2617 if (!fTracker[iDet] && (iDet < 7)) {
2618 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2619 if (fStopOnError) return kFALSE;
2621 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2627 //_____________________________________________________________________________
2628 void AliReconstruction::CleanUp()
2630 // delete trackers and the run loader and close and delete the file
2632 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2633 delete fReconstructor[iDet];
2634 fReconstructor[iDet] = NULL;
2635 fLoader[iDet] = NULL;
2636 delete fTracker[iDet];
2637 fTracker[iDet] = NULL;
2645 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2646 delete fDiamondProfileSPD;
2647 fDiamondProfileSPD = NULL;
2648 delete fDiamondProfile;
2649 fDiamondProfile = NULL;
2650 delete fDiamondProfileTPC;
2651 fDiamondProfileTPC = NULL;
2657 delete fParentRawReader;
2658 fParentRawReader=NULL;
2667 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2669 // Write space-points which are then used in the alignment procedures
2670 // For the moment only ITS, TPC, TRD and TOF
2672 Int_t ntracks = esd->GetNumberOfTracks();
2673 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2675 AliESDtrack *track = esd->GetTrack(itrack);
2678 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2679 nsp += track->GetNcls(iDet);
2681 if (iDet==0) { // ITS "extra" clusters
2682 track->GetClusters(iDet,idx);
2683 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2688 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2689 track->SetTrackPointArray(sp);
2691 for (Int_t iDet = 5; iDet >= 0; iDet--) {
2692 AliTracker *tracker = fTracker[iDet];
2693 if (!tracker) continue;
2694 Int_t nspdet = track->GetClusters(iDet,idx);
2696 if (iDet==0) // ITS "extra" clusters
2697 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2699 if (nspdet <= 0) continue;
2703 while (isp2 < nspdet) {
2704 Bool_t isvalid=kTRUE;
2706 Int_t index=idx[isp++];
2707 if (index < 0) continue;
2709 TString dets = fgkDetectorName[iDet];
2710 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2711 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2712 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2713 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2714 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2716 isvalid = tracker->GetTrackPoint(index,p);
2719 if (!isvalid) continue;
2720 if (iDet==0 && (isp-1)>=6) p.SetExtra();
2721 sp->AddPoint(isptrack,&p); isptrack++;
2728 //_____________________________________________________________________________
2729 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2731 // The method reads the raw-data error log
2732 // accumulated within the rawReader.
2733 // It extracts the raw-data errors related to
2734 // the current event and stores them into
2735 // a TClonesArray inside the esd object.
2737 if (!fRawReader) return;
2739 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2741 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2743 if (iEvent != log->GetEventNumber()) continue;
2745 esd->AddRawDataErrorLog(log);
2750 //_____________________________________________________________________________
2751 void AliReconstruction::CheckQA()
2753 // check the QA of SIM for this run and remove the detectors
2754 // with status Fatal
2756 // TString newRunLocalReconstruction ;
2757 // TString newRunTracking ;
2758 // TString newFillESD ;
2760 // for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2761 // TString detName(AliQA::GetDetName(iDet)) ;
2762 // AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2763 // if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, specie, AliQA::kFATAL)) {
2764 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
2765 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
2767 // if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2768 // fRunLocalReconstruction.Contains("ALL") ) {
2769 // newRunLocalReconstruction += detName ;
2770 // newRunLocalReconstruction += " " ;
2772 // if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2773 // fRunTracking.Contains("ALL") ) {
2774 // newRunTracking += detName ;
2775 // newRunTracking += " " ;
2777 // if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2778 // fFillESD.Contains("ALL") ) {
2779 // newFillESD += detName ;
2780 // newFillESD += " " ;
2784 // fRunLocalReconstruction = newRunLocalReconstruction ;
2785 // fRunTracking = newRunTracking ;
2786 // fFillESD = newFillESD ;
2789 //_____________________________________________________________________________
2790 Int_t AliReconstruction::GetDetIndex(const char* detector)
2792 // return the detector index corresponding to detector
2794 for (index = 0; index < kNDetectors ; index++) {
2795 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2800 //_____________________________________________________________________________
2801 Bool_t AliReconstruction::FinishPlaneEff() {
2803 // Here execute all the necessary operationis, at the end of the tracking phase,
2804 // in case that evaluation of PlaneEfficiencies was required for some detector.
2805 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2807 // This Preliminary version works only FOR ITS !!!!!
2808 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2811 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2814 //for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2815 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2816 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2817 if(fTracker[iDet]) {
2818 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2819 TString name=planeeff->GetName();
2821 TFile* pefile = TFile::Open(name, "RECREATE");
2822 ret=(Bool_t)planeeff->Write();
2824 if(planeeff->GetCreateHistos()) {
2825 TString hname=planeeff->GetName();
2826 hname+="Histo.root";
2827 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
2833 //_____________________________________________________________________________
2834 Bool_t AliReconstruction::InitPlaneEff() {
2836 // Here execute all the necessary operations, before of the tracking phase,
2837 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2838 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2839 // which should be updated/recalculated.
2841 // This Preliminary version will work only FOR ITS !!!!!
2842 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2845 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2847 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2851 //_____________________________________________________________________________
2852 Bool_t AliReconstruction::InitAliEVE()
2854 // This method should be called only in case
2855 // AliReconstruction is run
2856 // within the alieve environment.
2857 // It will initialize AliEVE in a way
2858 // so that it can visualize event processed
2859 // by AliReconstruction.
2860 // The return flag shows whenever the
2861 // AliEVE initialization was successful or not.
2864 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2865 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2866 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2868 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
2869 gROOT->ProcessLine("alieve_online_init()");
2874 //_____________________________________________________________________________
2875 void AliReconstruction::RunAliEVE()
2877 // Runs AliEVE visualisation of
2878 // the current event.
2879 // Should be executed only after
2880 // successful initialization of AliEVE.
2882 AliInfo("Running AliEVE...");
2883 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
2887 //_____________________________________________________________________________
2888 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2890 // Allows to run QA for a selected set of detectors
2891 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2892 // all selected detectors run the same selected tasks
2894 if (!detAndAction.Contains(":")) {
2895 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2899 Int_t colon = detAndAction.Index(":") ;
2900 fQADetectors = detAndAction(0, colon) ;
2901 if (fQADetectors.Contains("ALL") )
2902 fQADetectors = fFillESD ;
2903 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2904 if (fQATasks.Contains("ALL") ) {
2905 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2907 fQATasks.ToUpper() ;
2909 if ( fQATasks.Contains("RAW") )
2910 tempo = Form("%d ", AliQA::kRAWS) ;
2911 if ( fQATasks.Contains("RECPOINT") )
2912 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2913 if ( fQATasks.Contains("ESD") )
2914 tempo += Form("%d ", AliQA::kESDS) ;
2916 if (fQATasks.IsNull()) {
2917 AliInfo("No QA requested\n") ;
2922 TString tempo(fQATasks) ;
2923 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2924 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2925 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2926 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2931 //_____________________________________________________________________________
2932 Bool_t AliReconstruction::InitRecoParams()
2934 // The method accesses OCDB and retrieves all
2935 // the available reco-param objects from there.
2937 Bool_t isOK = kTRUE;
2939 TString detStr = fLoadCDB;
2940 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2942 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2944 if (fRecoParam.GetDetRecoParamArray(iDet)) {
2945 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
2949 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
2951 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
2952 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
2954 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
2958 TObject *recoParamObj = entry->GetObject();
2959 if (dynamic_cast<TObjArray*>(recoParamObj)) {
2960 // The detector has a normal TobjArray of AliDetectorRecoParam objects
2961 // Registering them in AliRecoParam
2962 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
2964 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
2965 // The detector has only onse set of reco parameters
2966 // Registering it in AliRecoParam
2967 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
2968 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
2969 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
2972 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
2976 AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
2980 if (AliDebugLevel() > 0) fRecoParam.Print();
2985 //_____________________________________________________________________________
2986 Bool_t AliReconstruction::GetEventInfo()
2988 // Fill the event info object
2990 AliCodeTimerAuto("")
2992 AliCentralTrigger *aCTP = NULL;
2994 fEventInfo.SetEventType(fRawReader->GetType());
2996 ULong64_t mask = fRawReader->GetClassMask();
2997 fEventInfo.SetTriggerMask(mask);
2998 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
2999 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3001 aCTP = new AliCentralTrigger();
3002 TString configstr("");
3003 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3004 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3008 aCTP->SetClassMask(mask);
3009 aCTP->SetClusterMask(clmask);
3012 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3014 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3015 aCTP = fRunLoader->GetTrigger();
3016 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3017 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3020 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3025 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3027 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3028 if (fRawReader) delete aCTP;
3032 UChar_t clustmask = 0;
3034 ULong64_t trmask = fEventInfo.GetTriggerMask();
3035 const TObjArray& classesArray = config->GetClasses();
3036 Int_t nclasses = classesArray.GetEntriesFast();
3037 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3038 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3040 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3041 fesd->SetTriggerClass(trclass->GetName(),trindex);
3042 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3043 if (trmask & (1 << trindex)) {
3045 trclasses += trclass->GetName();
3047 clustmask |= trclass->GetCluster()->GetClusterMask();
3051 fEventInfo.SetTriggerClasses(trclasses);
3053 // Set the information in ESD
3054 fesd->SetTriggerMask(trmask);
3055 fesd->SetTriggerCluster(clustmask);
3057 if (!aCTP->CheckTriggeredDetectors()) {
3058 if (fRawReader) delete aCTP;
3062 if (fRawReader) delete aCTP;
3064 // We have to fill also the HLT decision here!!
3070 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3072 // Match the detector list found in the rec.C or the default 'ALL'
3073 // to the list found in the GRP (stored there by the shuttle PP which
3074 // gets the information from ECS)
3075 static TString resultList;
3076 TString detList = detectorList;
3080 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3081 if ((detectorMask >> iDet) & 0x1) {
3082 TString det = AliDAQ::OfflineModuleName(iDet);
3083 if ((detList.CompareTo("ALL") == 0) ||
3084 ((detList.BeginsWith("ALL ") ||
3085 detList.EndsWith(" ALL") ||
3086 detList.Contains(" ALL ")) &&
3087 !(detList.BeginsWith("-"+det+" ") ||
3088 detList.EndsWith(" -"+det) ||
3089 detList.Contains(" -"+det+" "))) ||
3090 (detList.CompareTo(det) == 0) ||
3091 detList.BeginsWith(det+" ") ||
3092 detList.EndsWith(" "+det) ||
3093 detList.Contains( " "+det+" " )) {
3094 if (!resultList.EndsWith(det + " ")) {
3103 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3104 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3105 if ((detList.CompareTo("ALL") == 0) ||
3106 ((detList.BeginsWith("ALL ") ||
3107 detList.EndsWith(" ALL") ||
3108 detList.Contains(" ALL ")) &&
3109 !(detList.BeginsWith("-"+hltDet+" ") ||
3110 detList.EndsWith(" -"+hltDet) ||
3111 detList.Contains(" -"+hltDet+" "))) ||
3112 (detList.CompareTo(hltDet) == 0) ||
3113 detList.BeginsWith(hltDet+" ") ||
3114 detList.EndsWith(" "+hltDet) ||
3115 detList.Contains( " "+hltDet+" " )) {
3116 resultList += hltDet;
3120 return resultList.Data();
3124 //______________________________________________________________________________
3125 void AliReconstruction::Abort(const char *method, EAbort what)
3127 // Abort processing. If what = kAbortProcess, the Process() loop will be
3128 // aborted. If what = kAbortFile, the current file in a chain will be
3129 // aborted and the processing will continue with the next file, if there
3130 // is no next file then Process() will be aborted. Abort() can also be
3131 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3132 // the SlaveTerminate() and Terminate() are always called. The abort flag
3133 // can be checked in these methods using GetAbort().
3135 // The method is overwritten in AliReconstruction for better handling of
3136 // reco specific errors
3138 if (!fStopOnError) return;
3142 TString whyMess = method;
3143 whyMess += " failed! Aborting...";
3145 AliError(whyMess.Data());
3148 TString mess = "Abort";
3149 if (fAbort == kAbortProcess)
3150 mess = "AbortProcess";
3151 else if (fAbort == kAbortFile)
3154 Info(mess, whyMess.Data());