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 <TPluginManager.h>
117 #include <TGeoManager.h>
118 #include <TLorentzVector.h>
121 #include <TObjArray.h>
125 #include <TProofOutputFile.h>
126 #include <TParameter.h>
128 #include "AliReconstruction.h"
129 #include "AliCodeTimer.h"
130 #include "AliReconstructor.h"
132 #include "AliRunLoader.h"
134 #include "AliRawReaderFile.h"
135 #include "AliRawReaderDate.h"
136 #include "AliRawReaderRoot.h"
137 #include "AliRawEventHeaderBase.h"
138 #include "AliRawEvent.h"
139 #include "AliESDEvent.h"
140 #include "AliESDMuonTrack.h"
141 #include "AliESDfriend.h"
142 #include "AliESDVertex.h"
143 #include "AliESDcascade.h"
144 #include "AliESDkink.h"
145 #include "AliESDtrack.h"
146 #include "AliESDCaloCluster.h"
147 #include "AliESDCaloCells.h"
148 #include "AliMultiplicity.h"
149 #include "AliTracker.h"
150 #include "AliVertexer.h"
151 #include "AliVertexerTracks.h"
152 #include "AliV0vertexer.h"
153 #include "AliCascadeVertexer.h"
154 #include "AliHeader.h"
155 #include "AliGenEventHeader.h"
157 #include "AliESDpid.h"
158 #include "AliESDtrack.h"
159 #include "AliESDPmdTrack.h"
161 #include "AliESDTagCreator.h"
163 #include "AliGeomManager.h"
164 #include "AliTrackPointArray.h"
165 #include "AliCDBManager.h"
166 #include "AliCDBStorage.h"
167 #include "AliCDBEntry.h"
168 #include "AliAlignObj.h"
170 #include "AliCentralTrigger.h"
171 #include "AliTriggerConfiguration.h"
172 #include "AliTriggerClass.h"
173 #include "AliTriggerCluster.h"
174 #include "AliCTPRawStream.h"
176 #include "AliQADataMakerRec.h"
177 #include "AliGlobalQADataMaker.h"
179 #include "AliQADataMakerSteer.h"
181 #include "AliPlaneEff.h"
183 #include "AliSysInfo.h" // memory snapshots
184 #include "AliRawHLTManager.h"
186 #include "AliMagWrapCheb.h"
188 #include "AliDetectorRecoParam.h"
189 #include "AliGRPRecoParam.h"
190 #include "AliRunInfo.h"
191 #include "AliEventInfo.h"
195 #include "AliGRPObject.h"
197 ClassImp(AliReconstruction)
199 //_____________________________________________________________________________
200 const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
202 //_____________________________________________________________________________
203 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
205 fUniformField(kFALSE),
206 fForcedFieldMap(NULL),
207 fRunVertexFinder(kTRUE),
208 fRunVertexFinderTracks(kTRUE),
209 fRunHLTTracking(kFALSE),
210 fRunMuonTracking(kFALSE),
212 fRunCascadeFinder(kTRUE),
213 fStopOnError(kFALSE),
214 fWriteAlignmentData(kFALSE),
215 fWriteESDfriend(kFALSE),
216 fFillTriggerESD(kTRUE),
224 fRunLocalReconstruction("ALL"),
228 fUseTrackingErrorsForAlignment(""),
229 fGAliceFileName(gAliceFilename),
234 fNumberOfEventsPerFile((UInt_t)-1),
236 fLoadAlignFromCDB(kTRUE),
237 fLoadAlignData("ALL"),
244 fParentRawReader(NULL),
248 fDiamondProfileSPD(NULL),
249 fDiamondProfile(NULL),
250 fDiamondProfileTPC(NULL),
254 fAlignObjArray(NULL),
257 fInitCDBCalled(kFALSE),
258 fSetRunNumberFromDataCalled(kFALSE),
264 fSameQACycle(kFALSE),
266 fRunPlaneEff(kFALSE),
275 fIsNewRunLoader(kFALSE),
279 // create reconstruction object with default parameters
282 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
283 fReconstructor[iDet] = NULL;
284 fLoader[iDet] = NULL;
285 fTracker[iDet] = NULL;
287 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
288 fQACycles[iDet] = 999999 ;
289 fQAWriteExpert[iDet] = kFALSE ;
295 //_____________________________________________________________________________
296 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
298 fUniformField(rec.fUniformField),
299 fForcedFieldMap(NULL),
300 fRunVertexFinder(rec.fRunVertexFinder),
301 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
302 fRunHLTTracking(rec.fRunHLTTracking),
303 fRunMuonTracking(rec.fRunMuonTracking),
304 fRunV0Finder(rec.fRunV0Finder),
305 fRunCascadeFinder(rec.fRunCascadeFinder),
306 fStopOnError(rec.fStopOnError),
307 fWriteAlignmentData(rec.fWriteAlignmentData),
308 fWriteESDfriend(rec.fWriteESDfriend),
309 fFillTriggerESD(rec.fFillTriggerESD),
311 fCleanESD(rec.fCleanESD),
312 fV0DCAmax(rec.fV0DCAmax),
313 fV0CsPmin(rec.fV0CsPmin),
317 fRunLocalReconstruction(rec.fRunLocalReconstruction),
318 fRunTracking(rec.fRunTracking),
319 fFillESD(rec.fFillESD),
320 fLoadCDB(rec.fLoadCDB),
321 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
322 fGAliceFileName(rec.fGAliceFileName),
323 fRawInput(rec.fRawInput),
324 fEquipIdMap(rec.fEquipIdMap),
325 fFirstEvent(rec.fFirstEvent),
326 fLastEvent(rec.fLastEvent),
327 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
329 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
330 fLoadAlignData(rec.fLoadAlignData),
331 fUseHLTData(rec.fUseHLTData),
337 fParentRawReader(NULL),
339 fRecoParam(rec.fRecoParam),
341 fDiamondProfileSPD(rec.fDiamondProfileSPD),
342 fDiamondProfile(rec.fDiamondProfile),
343 fDiamondProfileTPC(rec.fDiamondProfileTPC),
347 fAlignObjArray(rec.fAlignObjArray),
348 fCDBUri(rec.fCDBUri),
350 fInitCDBCalled(rec.fInitCDBCalled),
351 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
352 fQADetectors(rec.fQADetectors),
354 fQATasks(rec.fQATasks),
356 fRunGlobalQA(rec.fRunGlobalQA),
357 fSameQACycle(rec.fSameQACycle),
358 fRunPlaneEff(rec.fRunPlaneEff),
367 fIsNewRunLoader(rec.fIsNewRunLoader),
373 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
374 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
376 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
377 fReconstructor[iDet] = NULL;
378 fLoader[iDet] = NULL;
379 fTracker[iDet] = NULL;
382 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
383 fQACycles[iDet] = rec.fQACycles[iDet];
384 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
387 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
388 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
392 //_____________________________________________________________________________
393 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
395 // assignment operator
396 // Used in PROOF mode
397 // Be very careful while modifing it!
398 // Simple rules to follow:
399 // for persistent data members - use their assignment operators
400 // for non-persistent ones - do nothing or take the default values from constructor
401 // TSelector members should not be touched
402 if(&rec == this) return *this;
404 fUniformField = rec.fUniformField;
405 fForcedFieldMap = NULL;
406 fRunVertexFinder = rec.fRunVertexFinder;
407 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
408 fRunHLTTracking = rec.fRunHLTTracking;
409 fRunMuonTracking = rec.fRunMuonTracking;
410 fRunV0Finder = rec.fRunV0Finder;
411 fRunCascadeFinder = rec.fRunCascadeFinder;
412 fStopOnError = rec.fStopOnError;
413 fWriteAlignmentData = rec.fWriteAlignmentData;
414 fWriteESDfriend = rec.fWriteESDfriend;
415 fFillTriggerESD = rec.fFillTriggerESD;
417 fCleanESD = rec.fCleanESD;
418 fV0DCAmax = rec.fV0DCAmax;
419 fV0CsPmin = rec.fV0CsPmin;
423 fRunLocalReconstruction = rec.fRunLocalReconstruction;
424 fRunTracking = rec.fRunTracking;
425 fFillESD = rec.fFillESD;
426 fLoadCDB = rec.fLoadCDB;
427 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
428 fGAliceFileName = rec.fGAliceFileName;
429 fRawInput = rec.fRawInput;
430 fEquipIdMap = rec.fEquipIdMap;
431 fFirstEvent = rec.fFirstEvent;
432 fLastEvent = rec.fLastEvent;
433 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
435 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
436 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
439 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
440 fLoadAlignData = rec.fLoadAlignData;
441 fUseHLTData = rec.fUseHLTData;
443 delete fRunInfo; fRunInfo = NULL;
444 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
446 fEventInfo = rec.fEventInfo;
450 fParentRawReader = NULL;
452 fRecoParam = rec.fRecoParam;
454 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
455 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
456 delete fLoader[iDet]; fLoader[iDet] = NULL;
457 delete fTracker[iDet]; fTracker[iDet] = NULL;
460 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
461 fQACycles[iDet] = rec.fQACycles[iDet];
462 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
465 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
466 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
467 delete fDiamondProfile; fDiamondProfile = NULL;
468 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
469 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
470 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
472 delete fGRPData; fGRPData = NULL;
473 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
474 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
476 delete fAlignObjArray; fAlignObjArray = NULL;
479 fSpecCDBUri.Delete();
480 fInitCDBCalled = rec.fInitCDBCalled;
481 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
482 fQADetectors = rec.fQADetectors;
484 fQATasks = rec.fQATasks;
486 fRunGlobalQA = rec.fRunGlobalQA;
487 fSameQACycle = rec.fSameQACycle;
488 fRunPlaneEff = rec.fRunPlaneEff;
497 fIsNewRunLoader = rec.fIsNewRunLoader;
504 //_____________________________________________________________________________
505 AliReconstruction::~AliReconstruction()
511 delete fForcedFieldMap;
513 if (fAlignObjArray) {
514 fAlignObjArray->Delete();
515 delete fAlignObjArray;
517 fSpecCDBUri.Delete();
519 AliCodeTimer::Instance()->Print();
522 //_____________________________________________________________________________
523 void AliReconstruction::InitCDB()
525 // activate a default CDB storage
526 // First check if we have any CDB storage set, because it is used
527 // to retrieve the calibration and alignment constants
528 AliCodeTimerAuto("");
530 if (fInitCDBCalled) return;
531 fInitCDBCalled = kTRUE;
533 AliCDBManager* man = AliCDBManager::Instance();
534 if (man->IsDefaultStorageSet())
536 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
537 AliWarning("Default CDB storage has been already set !");
538 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
539 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
540 fCDBUri = man->GetDefaultStorage()->GetURI();
543 if (fCDBUri.Length() > 0)
545 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
546 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
547 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
549 fCDBUri="local://$ALICE_ROOT";
550 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
551 AliWarning("Default CDB storage not yet set !!!!");
552 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
553 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
556 man->SetDefaultStorage(fCDBUri);
559 // Now activate the detector specific CDB storage locations
560 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
561 TObject* obj = fSpecCDBUri[i];
563 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
564 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
565 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
566 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
568 AliSysInfo::AddStamp("InitCDB");
571 //_____________________________________________________________________________
572 void AliReconstruction::SetDefaultStorage(const char* uri) {
573 // Store the desired default CDB storage location
574 // Activate it later within the Run() method
580 //_____________________________________________________________________________
581 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
582 // Store a detector-specific CDB storage location
583 // Activate it later within the Run() method
585 AliCDBPath aPath(calibType);
586 if(!aPath.IsValid()){
587 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
588 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
589 if(!strcmp(calibType, fgkDetectorName[iDet])) {
590 aPath.SetPath(Form("%s/*", calibType));
591 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
595 if(!aPath.IsValid()){
596 AliError(Form("Not a valid path or detector: %s", calibType));
601 // // check that calibType refers to a "valid" detector name
602 // Bool_t isDetector = kFALSE;
603 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
604 // TString detName = fgkDetectorName[iDet];
605 // if(aPath.GetLevel0() == detName) {
606 // isDetector = kTRUE;
612 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
616 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
617 if (obj) fSpecCDBUri.Remove(obj);
618 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
622 //_____________________________________________________________________________
623 Bool_t AliReconstruction::SetRunNumberFromData()
625 // The method is called in Run() in order
626 // to set a correct run number.
627 // In case of raw data reconstruction the
628 // run number is taken from the raw data header
630 if (fSetRunNumberFromDataCalled) return kTRUE;
631 fSetRunNumberFromDataCalled = kTRUE;
633 AliCDBManager* man = AliCDBManager::Instance();
636 if(fRawReader->NextEvent()) {
637 if(man->GetRun() > 0) {
638 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
640 man->SetRun(fRawReader->GetRunNumber());
641 fRawReader->RewindEvents();
644 if(man->GetRun() > 0) {
645 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
648 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
654 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
656 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
661 // read run number from gAlice
662 if(rl->GetHeader()) {
663 man->SetRun(rl->GetHeader()->GetRun());
668 AliError("Neither run-loader header nor RawReader objects are found !");
680 //_____________________________________________________________________________
681 void AliReconstruction::SetCDBLock() {
682 // Set CDB lock: from now on it is forbidden to reset the run number
683 // or the default storage or to activate any further storage!
685 AliCDBManager::Instance()->SetLock(1);
688 //_____________________________________________________________________________
689 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
691 // Read the alignment objects from CDB.
692 // Each detector is supposed to have the
693 // alignment objects in DET/Align/Data CDB path.
694 // All the detector objects are then collected,
695 // sorted by geometry level (starting from ALIC) and
696 // then applied to the TGeo geometry.
697 // Finally an overlaps check is performed.
699 // Load alignment data from CDB and fill fAlignObjArray
700 if(fLoadAlignFromCDB){
702 TString detStr = detectors;
703 TString loadAlObjsListOfDets = "";
705 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
706 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
707 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
709 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
711 loadAlObjsListOfDets += fgkDetectorName[iDet];
712 loadAlObjsListOfDets += " ";
714 } // end loop over detectors
716 if(AliGeomManager::GetNalignable("GRP") != 0)
717 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
718 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
719 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
721 // Check if the array with alignment objects was
722 // provided by the user. If yes, apply the objects
723 // to the present TGeo geometry
724 if (fAlignObjArray) {
725 if (gGeoManager && gGeoManager->IsClosed()) {
726 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
727 AliError("The misalignment of one or more volumes failed!"
728 "Compare the list of simulated detectors and the list of detector alignment data!");
733 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
739 if (fAlignObjArray) {
740 fAlignObjArray->Delete();
741 delete fAlignObjArray; fAlignObjArray=NULL;
747 //_____________________________________________________________________________
748 void AliReconstruction::SetGAliceFile(const char* fileName)
750 // set the name of the galice file
752 fGAliceFileName = fileName;
755 //_____________________________________________________________________________
756 void AliReconstruction::SetInput(const char* input)
758 // In case the input string starts with 'mem://', we run in an online mode
759 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
760 // file is assumed. One can give as an input:
761 // mem://: - events taken from DAQ monitoring libs online
763 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
764 if (input) fRawInput = input;
767 //_____________________________________________________________________________
768 void AliReconstruction::SetOption(const char* detector, const char* option)
770 // set options for the reconstruction of a detector
772 TObject* obj = fOptions.FindObject(detector);
773 if (obj) fOptions.Remove(obj);
774 fOptions.Add(new TNamed(detector, option));
777 //_____________________________________________________________________________
778 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
780 // Set custom reconstruction parameters for a given detector
781 // Single set of parameters for all the events
783 // First check if the reco-params are global
784 if(!strcmp(detector, "GRP")) {
786 fRecoParam.AddDetRecoParam(kNDetectors,par);
790 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
791 if(!strcmp(detector, fgkDetectorName[iDet])) {
793 fRecoParam.AddDetRecoParam(iDet,par);
800 //_____________________________________________________________________________
801 Bool_t AliReconstruction::SetFieldMap(Float_t l3Current, Float_t diCurrent, Float_t factor, const char *path) {
802 //------------------------------------------------
803 // The magnetic field map, defined externally...
804 // L3 current 30000 A -> 0.5 T
805 // L3 current 12000 A -> 0.2 T
806 // dipole current 6000 A
807 // The polarities must be the same
808 //------------------------------------------------
809 const Float_t l3NominalCurrent1=30000.; // (A)
810 const Float_t l3NominalCurrent2=12000.; // (A)
811 const Float_t diNominalCurrent =6000. ; // (A)
813 const Float_t tolerance=0.03; // relative current tolerance
814 const Float_t zero=77.; // "zero" current (A)
817 Bool_t dipoleON=kFALSE;
819 TString s=(factor < 0) ? "L3: -" : "L3: +";
821 l3Current = TMath::Abs(l3Current);
822 if (TMath::Abs(l3Current-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
823 map=AliMagWrapCheb::k5kG;
826 if (TMath::Abs(l3Current-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
827 map=AliMagWrapCheb::k2kG;
830 if (l3Current < zero) {
831 map=AliMagWrapCheb::k2kG;
833 factor=0.; // in fact, this is a global factor...
834 fUniformField=kTRUE; // track with the uniform (zero) B field
836 AliError(Form("Wrong L3 current (%f A)!",l3Current));
840 diCurrent = TMath::Abs(diCurrent);
841 if (TMath::Abs(diCurrent-diNominalCurrent)/diNominalCurrent < tolerance) {
842 // 3% current tolerance...
846 if (diCurrent < zero) { // some small current..
850 AliError(Form("Wrong dipole current (%f A)!",diCurrent));
854 delete fForcedFieldMap;
856 new AliMagWrapCheb("B field map ",s,2,factor,10.,map,dipoleON,path);
858 fForcedFieldMap->Print();
860 AliTracker::SetFieldMap(fForcedFieldMap,fUniformField);
866 Bool_t AliReconstruction::InitGRP() {
867 //------------------------------------
868 // Initialization of the GRP entry
869 //------------------------------------
870 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
874 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
877 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
879 fGRPData = new AliGRPObject();
880 fGRPData->ReadValuesFromMap(m);
884 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
885 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
889 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
893 AliError("No GRP entry found in OCDB!");
897 TString lhcState = fGRPData->GetLHCState();
898 if (lhcState==AliGRPObject::GetInvalidString()) {
899 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
900 lhcState = "UNKNOWN";
903 TString beamType = fGRPData->GetBeamType();
904 if (beamType==AliGRPObject::GetInvalidString()) {
905 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
906 beamType = "UNKNOWN";
909 Float_t beamEnergy = fGRPData->GetBeamEnergy();
910 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
911 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
915 TString runType = fGRPData->GetRunType();
916 if (runType==AliGRPObject::GetInvalidString()) {
917 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
921 Int_t activeDetectors = fGRPData->GetDetectorMask();
922 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
923 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
924 activeDetectors = 1074790399;
927 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
932 // Process the list of active detectors
933 if (activeDetectors) {
934 UInt_t detMask = activeDetectors;
935 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
936 fRunTracking = MatchDetectorList(fRunTracking,detMask);
937 fFillESD = MatchDetectorList(fFillESD,detMask);
938 fQADetectors = MatchDetectorList(fQADetectors,detMask);
939 fLoadCDB.Form("%s %s %s %s",
940 fRunLocalReconstruction.Data(),
943 fQADetectors.Data());
944 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
945 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1)) {
946 // switch off the vertexer
947 AliInfo("SPD is not in the list of active detectors. Vertexer switched off.");
948 fRunVertexFinder = kFALSE;
950 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
951 // switch off the reading of CTP raw-data payload
952 if (fFillTriggerESD) {
953 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
954 fFillTriggerESD = kFALSE;
959 AliInfo("===================================================================================");
960 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
961 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
962 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
963 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
964 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
965 AliInfo("===================================================================================");
967 //*** Dealing with the magnetic field map
968 if (AliTracker::GetFieldMap()) {
969 AliInfo("Running with the externally set B field !");
971 // Construct the field map out of the information retrieved from GRP.
976 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
977 if (l3Current == AliGRPObject::GetInvalidFloat()) {
978 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
982 Char_t l3Polarity = fGRPData->GetL3Polarity();
983 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
984 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
989 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
990 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
991 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
995 Char_t diPolarity = fGRPData->GetDipolePolarity();
996 if (diPolarity == AliGRPObject::GetInvalidChar()) {
997 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1002 TObjString *l3Current=
1003 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
1005 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1008 TObjString *l3Polarity=
1009 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
1011 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1016 TObjString *diCurrent=
1017 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
1019 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1022 TObjString *diPolarity=
1023 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
1025 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1031 Float_t l3Cur=TMath::Abs(l3Current);
1032 Float_t diCur=TMath::Abs(diCurrent);
1033 Float_t l3Pol=l3Polarity;
1034 // Float_t l3Cur=TMath::Abs(atof(l3Current->GetName()));
1035 //Float_t diCur=TMath::Abs(atof(diCurrent->GetName()));
1036 //Float_t l3Pol=atof(l3Polarity->GetName());
1038 if (l3Pol != 0.) factor=-1.;
1041 if (!SetFieldMap(l3Cur, diCur, factor)) {
1042 AliFatal("Failed to creat a B field map ! Exiting...");
1044 AliInfo("Running with the B field constructed out of GRP !");
1047 AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1052 //*** Get the diamond profiles from OCDB
1053 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1055 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1057 AliError("No SPD diamond profile found in OCDB!");
1060 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1062 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1064 AliError("No diamond profile found in OCDB!");
1067 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1069 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1071 AliError("No TPC diamond profile found in OCDB!");
1077 //_____________________________________________________________________________
1078 Bool_t AliReconstruction::LoadCDB()
1080 AliCodeTimerAuto("");
1082 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1084 TString detStr = fLoadCDB;
1085 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1086 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1087 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1092 //_____________________________________________________________________________
1093 Bool_t AliReconstruction::Run(const char* input)
1096 AliCodeTimerAuto("");
1099 if (GetAbort() != TSelector::kContinue) return kFALSE;
1101 TChain *chain = NULL;
1102 if (fRawReader && (chain = fRawReader->GetChain())) {
1105 gProof->AddInput(this);
1107 outputFile.SetProtocol("root",kTRUE);
1108 outputFile.SetHost(gSystem->HostName());
1109 outputFile.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1110 AliInfo(Form("Output file with ESDs is %s",outputFile.GetUrl()));
1111 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.GetUrl()));
1113 chain->Process("AliReconstruction");
1116 chain->Process(this);
1121 if (GetAbort() != TSelector::kContinue) return kFALSE;
1123 if (GetAbort() != TSelector::kContinue) return kFALSE;
1124 //******* The loop over events
1125 AliInfo("Starting looping over events");
1127 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1128 (fRawReader && fRawReader->NextEvent())) {
1129 if (!ProcessEvent(iEvent)) {
1130 Abort("ProcessEvent",TSelector::kAbortFile);
1136 if (GetAbort() != TSelector::kContinue) return kFALSE;
1138 if (GetAbort() != TSelector::kContinue) return kFALSE;
1144 //_____________________________________________________________________________
1145 void AliReconstruction::InitRawReader(const char* input)
1147 AliCodeTimerAuto("");
1149 // Init raw-reader and
1150 // set the input in case of raw data
1151 if (input) fRawInput = input;
1152 fRawReader = AliRawReader::Create(fRawInput.Data());
1154 AliInfo("Reconstruction will run over digits");
1156 if (!fEquipIdMap.IsNull() && fRawReader)
1157 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1159 if (!fUseHLTData.IsNull()) {
1160 // create the RawReaderHLT which performs redirection of HLT input data for
1161 // the specified detectors
1162 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1164 fParentRawReader=fRawReader;
1165 fRawReader=pRawReader;
1167 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1170 AliSysInfo::AddStamp("CreateRawReader");
1173 //_____________________________________________________________________________
1174 void AliReconstruction::InitRun(const char* input)
1176 // Initialization of raw-reader,
1177 // run number, CDB etc.
1178 AliCodeTimerAuto("");
1179 AliSysInfo::AddStamp("Start");
1181 // Initialize raw-reader if any
1182 InitRawReader(input);
1184 // Initialize the CDB storage
1187 // Set run number in CDBManager (if it is not already set by the user)
1188 if (!SetRunNumberFromData()) {
1189 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1193 // Set CDB lock: from now on it is forbidden to reset the run number
1194 // or the default storage or to activate any further storage!
1199 //_____________________________________________________________________________
1200 void AliReconstruction::Begin(TTree *)
1202 // Initialize AlReconstruction before
1203 // going into the event loop
1204 // Should follow the TSelector convention
1205 // i.e. initialize only the object on the client side
1206 AliCodeTimerAuto("");
1208 AliReconstruction *reco = NULL;
1210 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1213 AliSysInfo::AddStamp("ReadInputInBegin");
1216 // Import ideal TGeo geometry and apply misalignment
1218 TString geom(gSystem->DirName(fGAliceFileName));
1219 geom += "/geometry.root";
1220 AliGeomManager::LoadGeometry(geom.Data());
1222 Abort("LoadGeometry", TSelector::kAbortProcess);
1225 AliSysInfo::AddStamp("LoadGeom");
1226 TString detsToCheck=fRunLocalReconstruction;
1227 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1228 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1231 AliSysInfo::AddStamp("CheckGeom");
1234 if (!MisalignGeometry(fLoadAlignData)) {
1235 Abort("MisalignGeometry", TSelector::kAbortProcess);
1238 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1239 AliSysInfo::AddStamp("MisalignGeom");
1242 Abort("InitGRP", TSelector::kAbortProcess);
1245 AliSysInfo::AddStamp("InitGRP");
1248 Abort("LoadCDB", TSelector::kAbortProcess);
1251 AliSysInfo::AddStamp("LoadCDB");
1253 // Read the reconstruction parameters from OCDB
1254 if (!InitRecoParams()) {
1255 AliWarning("Not all detectors have correct RecoParam objects initialized");
1257 AliSysInfo::AddStamp("InitRecoParams");
1260 if (reco) *reco = *this;
1261 fInput->Add(gGeoManager);
1263 fInput->Add(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()));
1264 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1265 AliMagF *magFieldMap = (AliMagF*)AliTracker::GetFieldMap();
1266 magFieldMap->SetName("MagneticFieldMap");
1267 fInput->Add(magFieldMap);
1272 //_____________________________________________________________________________
1273 void AliReconstruction::SlaveBegin(TTree*)
1275 // Initialization related to run-loader,
1276 // vertexer, trackers, recontructors
1277 // In proof mode it is executed on the slave
1278 AliCodeTimerAuto("");
1280 TProofOutputFile *outProofFile = NULL;
1282 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1285 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1287 AliGeomManager::SetGeometry(tgeo);
1289 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1290 Int_t runNumber = -1;
1291 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1292 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1293 man->SetCacheFlag(kTRUE);
1294 man->SetLock(kTRUE);
1298 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1299 AliTracker::SetFieldMap(map,fUniformField);
1301 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1302 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1303 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1304 fOutput->Add(outProofFile);
1306 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1309 // get the run loader
1310 if (!InitRunLoader()) {
1311 Abort("InitRunLoader", TSelector::kAbortProcess);
1314 AliSysInfo::AddStamp("LoadLoader");
1316 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1319 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1320 Abort("CreateTrackers", TSelector::kAbortProcess);
1323 AliSysInfo::AddStamp("CreateTrackers");
1325 // create the ESD output file and tree
1326 if (!outProofFile) {
1327 ffile = TFile::Open("AliESDs.root", "RECREATE");
1328 ffile->SetCompressionLevel(2);
1329 if (!ffile->IsOpen()) {
1330 Abort("OpenESDFile", TSelector::kAbortProcess);
1335 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1336 Abort(Form("Problems opening output PROOF file: %s/%s",
1337 outProofFile->GetDir(), outProofFile->GetFileName()),
1338 TSelector::kAbortProcess);
1343 ftree = new TTree("esdTree", "Tree with ESD objects");
1344 fesd = new AliESDEvent();
1345 fesd->CreateStdContent();
1347 fesd->WriteToTree(ftree);
1348 if (fWriteESDfriend) {
1350 // Since we add the branch manually we must
1351 // book and add it after WriteToTree
1352 // otherwise it is created twice,
1353 // once via writetotree and once here.
1354 // The case for AliESDfriend is now
1355 // caught also in AlIESDEvent::WriteToTree but
1356 // be careful when changing the name (AliESDfriend is not
1357 // a TNamed so we had to hardwire it)
1358 fesdf = new AliESDfriend();
1359 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1360 br->SetFile("AliESDfriends.root");
1361 fesd->AddObject(fesdf);
1363 ftree->GetUserInfo()->Add(fesd);
1365 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1366 fhltesd = new AliESDEvent();
1367 fhltesd->CreateStdContent();
1369 // read the ESD template from CDB
1370 // HLT is allowed to put non-std content to its ESD, the non-std
1371 // objects need to be created before invocation of WriteToTree in
1372 // order to create all branches. Initialization is done from an
1373 // ESD layout template in CDB
1374 AliCDBManager* man = AliCDBManager::Instance();
1375 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1376 AliCDBEntry* hltESDConfig=NULL;
1377 if (man->GetId(hltESDConfigPath)!=NULL &&
1378 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1379 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1381 // init all internal variables from the list of objects
1382 pESDLayout->GetStdContent();
1384 // copy content and create non-std objects
1385 *fhltesd=*pESDLayout;
1388 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1389 hltESDConfigPath.GetPath().Data()));
1393 fhltesd->WriteToTree(fhlttree);
1394 fhlttree->GetUserInfo()->Add(fhltesd);
1396 ProcInfo_t procInfo;
1397 gSystem->GetProcInfo(&procInfo);
1398 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1401 //Initialize the QA and start of cycle
1403 fQASteer = new AliQADataMakerSteer("rec") ;
1404 fQASteer->SetActiveDetectors(fQADetectors) ;
1405 for (Int_t det = 0 ; det < AliQA::kNDET ; det++) {
1406 fQASteer->SetCycleLength(AliQA::DETECTORINDEX_t(det), fQACycles[det]) ;
1407 if (fQAWriteExpert[det])
1408 fQASteer->SetWriteExpert(AliQA::DETECTORINDEX_t(det)) ;
1411 if (!fRawReader && fQATasks.Contains(AliQA::kRAWS))
1412 fQATasks.ReplaceAll(Form("%d",AliQA::kRAWS), "") ;
1413 fQASteer->SetTasks(fQATasks) ;
1414 fQASteer->InitQADataMaker(AliCDBManager::Instance()->GetRun(), fRecoParam) ;
1418 Bool_t sameCycle = kFALSE ;
1419 if (!fQASteer) fQASteer = new AliQADataMakerSteer("rec") ;
1420 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1421 AliInfo(Form("Initializing the global QA data maker"));
1422 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1423 qadm->StartOfCycle(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1424 TObjArray *arr=qadm->Init(AliQA::kRECPOINTS);
1425 AliTracker::SetResidualsArray(arr);
1428 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1429 qadm->StartOfCycle(AliQA::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1430 qadm->Init(AliQA::kESDS);
1434 //Initialize the Plane Efficiency framework
1435 if (fRunPlaneEff && !InitPlaneEff()) {
1436 Abort("InitPlaneEff", TSelector::kAbortProcess);
1440 if (strcmp(gProgName,"alieve") == 0)
1441 fRunAliEVE = InitAliEVE();
1446 //_____________________________________________________________________________
1447 Bool_t AliReconstruction::Process(Long64_t entry)
1449 // run the reconstruction over a single entry
1450 // from the chain with raw data
1451 AliCodeTimerAuto("");
1453 TTree *currTree = fChain->GetTree();
1454 AliRawEvent *event = new AliRawEvent;
1455 currTree->SetBranchAddress("rawevent",&event);
1456 currTree->GetEntry(entry);
1457 fRawReader = new AliRawReaderRoot(event);
1458 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1466 //_____________________________________________________________________________
1467 void AliReconstruction::Init(TTree *tree)
1470 AliError("The input tree is not found!");
1476 //_____________________________________________________________________________
1477 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1479 // run the reconstruction over a single event
1480 // The event loop is steered in Run method
1482 AliCodeTimerAuto("");
1484 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1485 fRunLoader->SetEventNumber(iEvent);
1486 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1488 fRunLoader->TreeE()->Fill();
1489 if (fRawReader && fRawReader->UseAutoSaveESD())
1490 fRunLoader->TreeE()->AutoSave("SaveSelf");
1493 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1497 AliInfo(Form("processing event %d", iEvent));
1499 // Fill Event-info object
1501 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo);
1502 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1504 // Set the reco-params
1506 TString detStr = fLoadCDB;
1507 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1508 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1509 AliReconstructor *reconstructor = GetReconstructor(iDet);
1510 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1511 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1512 reconstructor->SetRecoParam(par);
1517 fRunLoader->GetEvent(iEvent);
1521 fQASteer->RunOneEvent(fRawReader) ;
1523 // local single event reconstruction
1524 if (!fRunLocalReconstruction.IsNull()) {
1525 TString detectors=fRunLocalReconstruction;
1526 // run HLT event reconstruction first
1527 // ;-( IsSelected changes the string
1528 if (IsSelected("HLT", detectors) &&
1529 !RunLocalEventReconstruction("HLT")) {
1530 if (fStopOnError) {CleanUp(); return kFALSE;}
1532 detectors=fRunLocalReconstruction;
1533 detectors.ReplaceAll("HLT", "");
1534 if (!RunLocalEventReconstruction(detectors)) {
1535 if (fStopOnError) {CleanUp(); return kFALSE;}
1539 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1540 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1541 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1542 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1544 // Set magnetic field from the tracker
1545 fesd->SetMagneticField(AliTracker::GetBz());
1546 fhltesd->SetMagneticField(AliTracker::GetBz());
1548 // Set most probable pt, for B=0 tracking
1549 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1550 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1551 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1553 // Fill raw-data error log into the ESD
1554 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1557 if (fRunVertexFinder) {
1558 if (!RunVertexFinder(fesd)) {
1559 if (fStopOnError) {CleanUp(); return kFALSE;}
1564 if (!fRunTracking.IsNull()) {
1565 if (fRunMuonTracking) {
1566 if (!RunMuonTracking(fesd)) {
1567 if (fStopOnError) {CleanUp(); return kFALSE;}
1573 if (!fRunTracking.IsNull()) {
1574 if (!RunTracking(fesd)) {
1575 if (fStopOnError) {CleanUp(); return kFALSE;}
1580 if (!fFillESD.IsNull()) {
1581 TString detectors=fFillESD;
1582 // run HLT first and on hltesd
1583 // ;-( IsSelected changes the string
1584 if (IsSelected("HLT", detectors) &&
1585 !FillESD(fhltesd, "HLT")) {
1586 if (fStopOnError) {CleanUp(); return kFALSE;}
1589 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1590 if (detectors.Contains("ALL")) {
1592 for (Int_t idet=0; idet<kNDetectors; ++idet){
1593 detectors += fgkDetectorName[idet];
1597 detectors.ReplaceAll("HLT", "");
1598 if (!FillESD(fesd, detectors)) {
1599 if (fStopOnError) {CleanUp(); return kFALSE;}
1603 // fill Event header information from the RawEventHeader
1604 if (fRawReader){FillRawEventHeaderESD(fesd);}
1607 AliESDpid::MakePID(fesd);
1609 if (fFillTriggerESD) {
1610 if (!FillTriggerESD(fesd)) {
1611 if (fStopOnError) {CleanUp(); return kFALSE;}
1618 // Propagate track to the beam pipe (if not already done by ITS)
1620 const Int_t ntracks = fesd->GetNumberOfTracks();
1621 const Double_t kBz = fesd->GetMagneticField();
1622 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1625 UShort_t *selectedIdx=new UShort_t[ntracks];
1627 for (Int_t itrack=0; itrack<ntracks; itrack++){
1628 const Double_t kMaxStep = 1; //max step over the material
1631 AliESDtrack *track = fesd->GetTrack(itrack);
1632 if (!track) continue;
1634 AliExternalTrackParam *tpcTrack =
1635 (AliExternalTrackParam *)track->GetTPCInnerParam();
1639 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1642 Int_t n=trkArray.GetEntriesFast();
1643 selectedIdx[n]=track->GetID();
1644 trkArray.AddLast(tpcTrack);
1647 //Tracks refitted by ITS should already be at the SPD vertex
1648 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1651 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1652 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1657 // Improve the reconstructed primary vertex position using the tracks
1659 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1660 if(fesd->GetPrimaryVertexSPD()) {
1661 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1662 if(vtitle.Contains("cosmics")) {
1663 runVertexFinderTracks=kFALSE;
1667 if (runVertexFinderTracks) {
1668 // TPC + ITS primary vertex
1669 ftVertexer->SetITSMode();
1670 ftVertexer->SetConstraintOff();
1671 // get cuts for vertexer from AliGRPRecoParam
1673 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1674 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1675 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1676 ftVertexer->SetCuts(cutsVertexer);
1677 delete [] cutsVertexer; cutsVertexer = NULL;
1678 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1679 ftVertexer->SetVtxStart(fDiamondProfile);
1681 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1683 if (pvtx->GetStatus()) {
1684 fesd->SetPrimaryVertexTracks(pvtx);
1685 for (Int_t i=0; i<ntracks; i++) {
1686 AliESDtrack *t = fesd->GetTrack(i);
1687 t->RelateToVertex(pvtx, kBz, kVeryBig);
1692 // TPC-only primary vertex
1693 ftVertexer->SetTPCMode();
1694 ftVertexer->SetConstraintOff();
1695 // get cuts for vertexer from AliGRPRecoParam
1697 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1698 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1699 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1700 ftVertexer->SetCuts(cutsVertexer);
1701 delete [] cutsVertexer; cutsVertexer = NULL;
1702 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1703 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1705 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1707 if (pvtx->GetStatus()) {
1708 fesd->SetPrimaryVertexTPC(pvtx);
1709 for (Int_t i=0; i<ntracks; i++) {
1710 AliESDtrack *t = fesd->GetTrack(i);
1711 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1717 delete[] selectedIdx;
1719 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1724 AliV0vertexer vtxer;
1725 vtxer.Tracks2V0vertices(fesd);
1727 if (fRunCascadeFinder) {
1729 AliCascadeVertexer cvtxer;
1730 cvtxer.V0sTracks2CascadeVertices(fesd);
1735 if (fCleanESD) CleanESD(fesd);
1738 fQASteer->RunOneEvent(fesd) ;
1741 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1742 if (qadm && fQATasks.Contains(Form("%d", AliQA::kESDS)))
1743 qadm->Exec(AliQA::kESDS, fesd);
1746 if (fWriteESDfriend) {
1747 // fesdf->~AliESDfriend();
1748 // new (fesdf) AliESDfriend(); // Reset...
1749 fesd->GetESDfriend(fesdf);
1753 // Auto-save the ESD tree in case of prompt reco @P2
1754 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1755 ftree->AutoSave("SaveSelf");
1756 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1757 if (friendfile) friendfile->Save();
1764 if (fRunAliEVE) RunAliEVE();
1768 if (fWriteESDfriend) {
1769 fesdf->~AliESDfriend();
1770 new (fesdf) AliESDfriend(); // Reset...
1773 ProcInfo_t procInfo;
1774 gSystem->GetProcInfo(&procInfo);
1775 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1778 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1779 if (fReconstructor[iDet])
1780 fReconstructor[iDet]->SetRecoParam(NULL);
1783 if (fRunQA || fRunGlobalQA)
1784 fQASteer->Increment() ;
1789 //_____________________________________________________________________________
1790 void AliReconstruction::SlaveTerminate()
1792 // Finalize the run on the slave side
1793 // Called after the exit
1794 // from the event loop
1795 AliCodeTimerAuto("");
1797 if (fIsNewRunLoader) { // galice.root didn't exist
1798 fRunLoader->WriteHeader("OVERWRITE");
1799 fRunLoader->CdGAFile();
1800 fRunLoader->Write(0, TObject::kOverwrite);
1803 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1804 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1806 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1807 cdbMapCopy->SetOwner(1);
1808 cdbMapCopy->SetName("cdbMap");
1809 TIter iter(cdbMap->GetTable());
1812 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1813 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1814 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1815 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1818 TList *cdbListCopy = new TList();
1819 cdbListCopy->SetOwner(1);
1820 cdbListCopy->SetName("cdbList");
1822 TIter iter2(cdbList);
1825 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1826 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1829 ftree->GetUserInfo()->Add(cdbMapCopy);
1830 ftree->GetUserInfo()->Add(cdbListCopy);
1835 if (fWriteESDfriend)
1836 ftree->SetBranchStatus("ESDfriend*",0);
1837 // we want to have only one tree version number
1838 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1841 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1842 if (fRunPlaneEff && !FinishPlaneEff()) {
1843 AliWarning("Finish PlaneEff evaluation failed");
1846 // End of cycle for the in-loop
1848 fQASteer->EndOfCycle() ;
1851 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1853 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1854 qadm->EndOfCycle(AliQA::kRECPOINTS);
1855 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1856 qadm->EndOfCycle(AliQA::kESDS);
1864 //_____________________________________________________________________________
1865 void AliReconstruction::Terminate()
1867 // Create tags for the events in the ESD tree (the ESD tree is always present)
1868 // In case of empty events the tags will contain dummy values
1869 AliCodeTimerAuto("");
1871 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1872 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1874 // Cleanup of CDB manager: cache and active storages!
1875 AliCDBManager::Instance()->ClearCache();
1878 //_____________________________________________________________________________
1879 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1881 // run the local reconstruction
1883 static Int_t eventNr=0;
1884 AliCodeTimerAuto("")
1886 TString detStr = detectors;
1887 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1888 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1889 AliReconstructor* reconstructor = GetReconstructor(iDet);
1890 if (!reconstructor) continue;
1891 AliLoader* loader = fLoader[iDet];
1892 // Matthias April 2008: temporary fix to run HLT reconstruction
1893 // although the HLT loader is missing
1894 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1896 reconstructor->Reconstruct(fRawReader, NULL);
1899 reconstructor->Reconstruct(dummy, NULL);
1904 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1907 // conversion of digits
1908 if (fRawReader && reconstructor->HasDigitConversion()) {
1909 AliInfo(Form("converting raw data digits into root objects for %s",
1910 fgkDetectorName[iDet]));
1911 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1912 // fgkDetectorName[iDet]));
1913 loader->LoadDigits("update");
1914 loader->CleanDigits();
1915 loader->MakeDigitsContainer();
1916 TTree* digitsTree = loader->TreeD();
1917 reconstructor->ConvertDigits(fRawReader, digitsTree);
1918 loader->WriteDigits("OVERWRITE");
1919 loader->UnloadDigits();
1921 // local reconstruction
1922 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1923 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1924 loader->LoadRecPoints("update");
1925 loader->CleanRecPoints();
1926 loader->MakeRecPointsContainer();
1927 TTree* clustersTree = loader->TreeR();
1928 if (fRawReader && !reconstructor->HasDigitConversion()) {
1929 reconstructor->Reconstruct(fRawReader, clustersTree);
1931 loader->LoadDigits("read");
1932 TTree* digitsTree = loader->TreeD();
1934 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1935 if (fStopOnError) return kFALSE;
1937 reconstructor->Reconstruct(digitsTree, clustersTree);
1939 loader->UnloadDigits();
1942 TString detQAStr(fQADetectors) ;
1944 fQASteer->RunOneEventInOneDetector(iDet, clustersTree) ;
1946 loader->WriteRecPoints("OVERWRITE");
1947 loader->UnloadRecPoints();
1948 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1950 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1951 AliError(Form("the following detectors were not found: %s",
1953 if (fStopOnError) return kFALSE;
1959 //_____________________________________________________________________________
1960 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1962 // run the barrel tracking
1964 AliCodeTimerAuto("")
1966 AliVertexer *vertexer = CreateVertexer();
1967 if (!vertexer) return kFALSE;
1969 AliInfo("running the ITS vertex finder");
1970 AliESDVertex* vertex = NULL;
1972 fLoader[0]->LoadRecPoints();
1973 TTree* cltree = fLoader[0]->TreeR();
1975 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
1976 vertex = vertexer->FindVertexForCurrentEvent(cltree);
1979 AliError("Can't get the ITS cluster tree");
1981 fLoader[0]->UnloadRecPoints();
1984 AliError("Can't get the ITS loader");
1987 AliWarning("Vertex not found");
1988 vertex = new AliESDVertex();
1989 vertex->SetName("default");
1992 vertex->SetName("reconstructed");
1997 vertex->GetXYZ(vtxPos);
1998 vertex->GetSigmaXYZ(vtxErr);
2000 esd->SetPrimaryVertexSPD(vertex);
2001 // if SPD multiplicity has been determined, it is stored in the ESD
2002 AliMultiplicity *mult = vertexer->GetMultiplicity();
2003 if(mult)esd->SetMultiplicity(mult);
2005 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2006 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2015 //_____________________________________________________________________________
2016 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2018 // run the HLT barrel tracking
2020 AliCodeTimerAuto("")
2023 AliError("Missing runLoader!");
2027 AliInfo("running HLT tracking");
2029 // Get a pointer to the HLT reconstructor
2030 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2031 if (!reconstructor) return kFALSE;
2034 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2035 TString detName = fgkDetectorName[iDet];
2036 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2037 reconstructor->SetOption(detName.Data());
2038 AliTracker *tracker = reconstructor->CreateTracker();
2040 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2041 if (fStopOnError) return kFALSE;
2045 Double_t vtxErr[3]={0.005,0.005,0.010};
2046 const AliESDVertex *vertex = esd->GetVertex();
2047 vertex->GetXYZ(vtxPos);
2048 tracker->SetVertex(vtxPos,vtxErr);
2050 fLoader[iDet]->LoadRecPoints("read");
2051 TTree* tree = fLoader[iDet]->TreeR();
2053 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2056 tracker->LoadClusters(tree);
2058 if (tracker->Clusters2Tracks(esd) != 0) {
2059 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2063 tracker->UnloadClusters();
2071 //_____________________________________________________________________________
2072 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2074 // run the muon spectrometer tracking
2076 AliCodeTimerAuto("")
2079 AliError("Missing runLoader!");
2082 Int_t iDet = 7; // for MUON
2084 AliInfo("is running...");
2086 // Get a pointer to the MUON reconstructor
2087 AliReconstructor *reconstructor = GetReconstructor(iDet);
2088 if (!reconstructor) return kFALSE;
2091 TString detName = fgkDetectorName[iDet];
2092 AliDebug(1, Form("%s tracking", detName.Data()));
2093 AliTracker *tracker = reconstructor->CreateTracker();
2095 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2100 fLoader[iDet]->LoadRecPoints("read");
2102 tracker->LoadClusters(fLoader[iDet]->TreeR());
2104 Int_t rv = tracker->Clusters2Tracks(esd);
2108 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2112 fLoader[iDet]->UnloadRecPoints();
2114 tracker->UnloadClusters();
2122 //_____________________________________________________________________________
2123 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2125 // run the barrel tracking
2126 static Int_t eventNr=0;
2127 AliCodeTimerAuto("")
2129 AliInfo("running tracking");
2131 //Fill the ESD with the T0 info (will be used by the TOF)
2132 if (fReconstructor[11] && fLoader[11]) {
2133 fLoader[11]->LoadRecPoints("READ");
2134 TTree *treeR = fLoader[11]->TreeR();
2136 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2140 // pass 1: TPC + ITS inwards
2141 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2142 if (!fTracker[iDet]) continue;
2143 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2146 fLoader[iDet]->LoadRecPoints("read");
2147 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2148 TTree* tree = fLoader[iDet]->TreeR();
2150 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2153 fTracker[iDet]->LoadClusters(tree);
2154 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2156 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2157 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2160 // preliminary PID in TPC needed by the ITS tracker
2162 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2163 AliESDpid::MakePID(esd);
2165 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2168 // pass 2: ALL backwards
2170 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2171 if (!fTracker[iDet]) continue;
2172 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2175 if (iDet > 1) { // all except ITS, TPC
2177 fLoader[iDet]->LoadRecPoints("read");
2178 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2179 tree = fLoader[iDet]->TreeR();
2181 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2184 fTracker[iDet]->LoadClusters(tree);
2185 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2189 if (iDet>1) // start filling residuals for the "outer" detectors
2190 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2192 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2193 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2198 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2199 fTracker[iDet]->UnloadClusters();
2200 fLoader[iDet]->UnloadRecPoints();
2202 // updated PID in TPC needed by the ITS tracker -MI
2204 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2205 AliESDpid::MakePID(esd);
2207 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2209 //stop filling residuals for the "outer" detectors
2210 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2212 // pass 3: TRD + TPC + ITS refit inwards
2214 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2215 if (!fTracker[iDet]) continue;
2216 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2219 if (iDet<2) // start filling residuals for TPC and ITS
2220 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2222 if (fTracker[iDet]->RefitInward(esd) != 0) {
2223 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2226 // run postprocessing
2227 if (fTracker[iDet]->PostProcess(esd) != 0) {
2228 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2231 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2234 // write space-points to the ESD in case alignment data output
2236 if (fWriteAlignmentData)
2237 WriteAlignmentData(esd);
2239 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2240 if (!fTracker[iDet]) continue;
2242 fTracker[iDet]->UnloadClusters();
2243 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2244 fLoader[iDet]->UnloadRecPoints();
2245 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2247 // stop filling residuals for TPC and ITS
2248 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2254 //_____________________________________________________________________________
2255 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2257 // Remove the data which are not needed for the physics analysis.
2260 Int_t nTracks=esd->GetNumberOfTracks();
2261 Int_t nV0s=esd->GetNumberOfV0s();
2263 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2265 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2266 Bool_t rc=esd->Clean(cleanPars);
2268 nTracks=esd->GetNumberOfTracks();
2269 nV0s=esd->GetNumberOfV0s();
2271 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2276 //_____________________________________________________________________________
2277 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2279 // fill the event summary data
2281 AliCodeTimerAuto("")
2282 static Int_t eventNr=0;
2283 TString detStr = detectors;
2285 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2286 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2287 AliReconstructor* reconstructor = GetReconstructor(iDet);
2288 if (!reconstructor) continue;
2289 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2290 TTree* clustersTree = NULL;
2291 if (fLoader[iDet]) {
2292 fLoader[iDet]->LoadRecPoints("read");
2293 clustersTree = fLoader[iDet]->TreeR();
2294 if (!clustersTree) {
2295 AliError(Form("Can't get the %s clusters tree",
2296 fgkDetectorName[iDet]));
2297 if (fStopOnError) return kFALSE;
2300 if (fRawReader && !reconstructor->HasDigitConversion()) {
2301 reconstructor->FillESD(fRawReader, clustersTree, esd);
2303 TTree* digitsTree = NULL;
2304 if (fLoader[iDet]) {
2305 fLoader[iDet]->LoadDigits("read");
2306 digitsTree = fLoader[iDet]->TreeD();
2308 AliError(Form("Can't get the %s digits tree",
2309 fgkDetectorName[iDet]));
2310 if (fStopOnError) return kFALSE;
2313 reconstructor->FillESD(digitsTree, clustersTree, esd);
2314 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2316 if (fLoader[iDet]) {
2317 fLoader[iDet]->UnloadRecPoints();
2321 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2322 AliError(Form("the following detectors were not found: %s",
2324 if (fStopOnError) return kFALSE;
2326 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2331 //_____________________________________________________________________________
2332 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2334 // Reads the trigger decision which is
2335 // stored in Trigger.root file and fills
2336 // the corresponding esd entries
2338 AliCodeTimerAuto("")
2340 AliInfo("Filling trigger information into the ESD");
2343 AliCTPRawStream input(fRawReader);
2344 if (!input.Next()) {
2345 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2348 if (esd->GetTriggerMask() != input.GetClassMask())
2349 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2350 input.GetClassMask(),esd->GetTriggerMask()));
2351 if (esd->GetOrbitNumber() != input.GetOrbitID())
2352 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2353 input.GetOrbitID(),esd->GetOrbitNumber()));
2354 if (esd->GetBunchCrossNumber() != input.GetBCID())
2355 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2356 input.GetBCID(),esd->GetBunchCrossNumber()));
2359 // Here one has to add the filling of trigger inputs and
2360 // interaction records
2370 //_____________________________________________________________________________
2371 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2374 // Filling information from RawReader Header
2377 if (!fRawReader) return kFALSE;
2379 AliInfo("Filling information from RawReader Header");
2381 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2382 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2383 esd->SetPeriodNumber(fRawReader->GetPeriod());
2385 esd->SetTimeStamp(fRawReader->GetTimestamp());
2386 esd->SetEventType(fRawReader->GetType());
2392 //_____________________________________________________________________________
2393 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2395 // check whether detName is contained in detectors
2396 // if yes, it is removed from detectors
2398 // check if all detectors are selected
2399 if ((detectors.CompareTo("ALL") == 0) ||
2400 detectors.BeginsWith("ALL ") ||
2401 detectors.EndsWith(" ALL") ||
2402 detectors.Contains(" ALL ")) {
2407 // search for the given detector
2408 Bool_t result = kFALSE;
2409 if ((detectors.CompareTo(detName) == 0) ||
2410 detectors.BeginsWith(detName+" ") ||
2411 detectors.EndsWith(" "+detName) ||
2412 detectors.Contains(" "+detName+" ")) {
2413 detectors.ReplaceAll(detName, "");
2417 // clean up the detectors string
2418 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2419 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2420 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2425 //_____________________________________________________________________________
2426 Bool_t AliReconstruction::InitRunLoader()
2428 // get or create the run loader
2430 if (gAlice) delete gAlice;
2433 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2434 // load all base libraries to get the loader classes
2435 TString libs = gSystem->GetLibraries();
2436 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2437 TString detName = fgkDetectorName[iDet];
2438 if (detName == "HLT") continue;
2439 if (libs.Contains("lib" + detName + "base.so")) continue;
2440 gSystem->Load("lib" + detName + "base.so");
2442 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2444 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2449 fRunLoader->CdGAFile();
2450 fRunLoader->LoadgAlice();
2452 //PH This is a temporary fix to give access to the kinematics
2453 //PH that is needed for the labels of ITS clusters
2454 fRunLoader->LoadHeader();
2455 fRunLoader->LoadKinematics();
2457 } else { // galice.root does not exist
2459 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2461 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2462 AliConfig::GetDefaultEventFolderName(),
2465 AliError(Form("could not create run loader in file %s",
2466 fGAliceFileName.Data()));
2470 fIsNewRunLoader = kTRUE;
2471 fRunLoader->MakeTree("E");
2473 if (fNumberOfEventsPerFile > 0)
2474 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2476 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2482 //_____________________________________________________________________________
2483 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2485 // get the reconstructor object and the loader for a detector
2487 if (fReconstructor[iDet]) {
2488 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2489 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2490 fReconstructor[iDet]->SetRecoParam(par);
2492 return fReconstructor[iDet];
2495 // load the reconstructor object
2496 TPluginManager* pluginManager = gROOT->GetPluginManager();
2497 TString detName = fgkDetectorName[iDet];
2498 TString recName = "Ali" + detName + "Reconstructor";
2500 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2502 AliReconstructor* reconstructor = NULL;
2503 // first check if a plugin is defined for the reconstructor
2504 TPluginHandler* pluginHandler =
2505 pluginManager->FindHandler("AliReconstructor", detName);
2506 // if not, add a plugin for it
2507 if (!pluginHandler) {
2508 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2509 TString libs = gSystem->GetLibraries();
2510 if (libs.Contains("lib" + detName + "base.so") ||
2511 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2512 pluginManager->AddHandler("AliReconstructor", detName,
2513 recName, detName + "rec", recName + "()");
2515 pluginManager->AddHandler("AliReconstructor", detName,
2516 recName, detName, recName + "()");
2518 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2520 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2521 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2523 if (reconstructor) {
2524 TObject* obj = fOptions.FindObject(detName.Data());
2525 if (obj) reconstructor->SetOption(obj->GetTitle());
2526 reconstructor->Init();
2527 fReconstructor[iDet] = reconstructor;
2530 // get or create the loader
2531 if (detName != "HLT") {
2532 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2533 if (!fLoader[iDet]) {
2534 AliConfig::Instance()
2535 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2537 // first check if a plugin is defined for the loader
2539 pluginManager->FindHandler("AliLoader", detName);
2540 // if not, add a plugin for it
2541 if (!pluginHandler) {
2542 TString loaderName = "Ali" + detName + "Loader";
2543 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2544 pluginManager->AddHandler("AliLoader", detName,
2545 loaderName, detName + "base",
2546 loaderName + "(const char*, TFolder*)");
2547 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2549 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2551 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2552 fRunLoader->GetEventFolder());
2554 if (!fLoader[iDet]) { // use default loader
2555 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2557 if (!fLoader[iDet]) {
2558 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2559 if (fStopOnError) return NULL;
2561 fRunLoader->AddLoader(fLoader[iDet]);
2562 fRunLoader->CdGAFile();
2563 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2564 fRunLoader->Write(0, TObject::kOverwrite);
2569 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2570 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2571 reconstructor->SetRecoParam(par);
2573 return reconstructor;
2576 //_____________________________________________________________________________
2577 AliVertexer* AliReconstruction::CreateVertexer()
2579 // create the vertexer
2580 // Please note that the caller is the owner of the
2583 AliVertexer* vertexer = NULL;
2584 AliReconstructor* itsReconstructor = GetReconstructor(0);
2585 if (itsReconstructor) {
2586 vertexer = itsReconstructor->CreateVertexer();
2589 AliWarning("couldn't create a vertexer for ITS");
2595 //_____________________________________________________________________________
2596 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2598 // create the trackers
2599 AliInfo("Creating trackers");
2601 TString detStr = detectors;
2602 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2603 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2604 AliReconstructor* reconstructor = GetReconstructor(iDet);
2605 if (!reconstructor) continue;
2606 TString detName = fgkDetectorName[iDet];
2607 if (detName == "HLT") {
2608 fRunHLTTracking = kTRUE;
2611 if (detName == "MUON") {
2612 fRunMuonTracking = kTRUE;
2617 fTracker[iDet] = reconstructor->CreateTracker();
2618 if (!fTracker[iDet] && (iDet < 7)) {
2619 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2620 if (fStopOnError) return kFALSE;
2622 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2628 //_____________________________________________________________________________
2629 void AliReconstruction::CleanUp()
2631 // delete trackers and the run loader and close and delete the file
2633 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2634 delete fReconstructor[iDet];
2635 fReconstructor[iDet] = NULL;
2636 fLoader[iDet] = NULL;
2637 delete fTracker[iDet];
2638 fTracker[iDet] = NULL;
2646 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2647 delete fDiamondProfileSPD;
2648 fDiamondProfileSPD = NULL;
2649 delete fDiamondProfile;
2650 fDiamondProfile = NULL;
2651 delete fDiamondProfileTPC;
2652 fDiamondProfileTPC = NULL;
2658 delete fParentRawReader;
2659 fParentRawReader=NULL;
2668 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2670 // Write space-points which are then used in the alignment procedures
2671 // For the moment only ITS, TPC, TRD and TOF
2673 Int_t ntracks = esd->GetNumberOfTracks();
2674 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2676 AliESDtrack *track = esd->GetTrack(itrack);
2679 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2680 nsp += track->GetNcls(iDet);
2682 if (iDet==0) { // ITS "extra" clusters
2683 track->GetClusters(iDet,idx);
2684 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2689 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2690 track->SetTrackPointArray(sp);
2692 for (Int_t iDet = 5; iDet >= 0; iDet--) {
2693 AliTracker *tracker = fTracker[iDet];
2694 if (!tracker) continue;
2695 Int_t nspdet = track->GetClusters(iDet,idx);
2697 if (iDet==0) // ITS "extra" clusters
2698 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2700 if (nspdet <= 0) continue;
2704 while (isp2 < nspdet) {
2705 Bool_t isvalid=kTRUE;
2707 Int_t index=idx[isp++];
2708 if (index < 0) continue;
2710 TString dets = fgkDetectorName[iDet];
2711 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2712 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2713 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2714 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2715 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2717 isvalid = tracker->GetTrackPoint(index,p);
2720 if (!isvalid) continue;
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, AliQA::kFATAL)) {
2764 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2766 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2767 fRunLocalReconstruction.Contains("ALL") ) {
2768 newRunLocalReconstruction += detName ;
2769 newRunLocalReconstruction += " " ;
2771 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2772 fRunTracking.Contains("ALL") ) {
2773 newRunTracking += detName ;
2774 newRunTracking += " " ;
2776 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2777 fFillESD.Contains("ALL") ) {
2778 newFillESD += detName ;
2783 fRunLocalReconstruction = newRunLocalReconstruction ;
2784 fRunTracking = newRunTracking ;
2785 fFillESD = newFillESD ;
2788 //_____________________________________________________________________________
2789 Int_t AliReconstruction::GetDetIndex(const char* detector)
2791 // return the detector index corresponding to detector
2793 for (index = 0; index < kNDetectors ; index++) {
2794 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2799 //_____________________________________________________________________________
2800 Bool_t AliReconstruction::FinishPlaneEff() {
2802 // Here execute all the necessary operationis, at the end of the tracking phase,
2803 // in case that evaluation of PlaneEfficiencies was required for some detector.
2804 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2806 // This Preliminary version works only FOR ITS !!!!!
2807 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2810 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2813 //for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2814 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2815 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2816 if(fTracker[iDet]) {
2817 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2818 TString name=planeeff->GetName();
2820 TFile* pefile = TFile::Open(name, "RECREATE");
2821 ret=(Bool_t)planeeff->Write();
2823 if(planeeff->GetCreateHistos()) {
2824 TString hname=planeeff->GetName();
2825 hname+="Histo.root";
2826 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
2832 //_____________________________________________________________________________
2833 Bool_t AliReconstruction::InitPlaneEff() {
2835 // Here execute all the necessary operations, before of the tracking phase,
2836 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2837 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2838 // which should be updated/recalculated.
2840 // This Preliminary version will work only FOR ITS !!!!!
2841 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2844 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2846 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2850 //_____________________________________________________________________________
2851 Bool_t AliReconstruction::InitAliEVE()
2853 // This method should be called only in case
2854 // AliReconstruction is run
2855 // within the alieve environment.
2856 // It will initialize AliEVE in a way
2857 // so that it can visualize event processed
2858 // by AliReconstruction.
2859 // The return flag shows whenever the
2860 // AliEVE initialization was successful or not.
2863 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2864 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2865 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2867 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
2868 gROOT->ProcessLine("alieve_online_init()");
2873 //_____________________________________________________________________________
2874 void AliReconstruction::RunAliEVE()
2876 // Runs AliEVE visualisation of
2877 // the current event.
2878 // Should be executed only after
2879 // successful initialization of AliEVE.
2881 AliInfo("Running AliEVE...");
2882 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2886 //_____________________________________________________________________________
2887 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2889 // Allows to run QA for a selected set of detectors
2890 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2891 // all selected detectors run the same selected tasks
2893 if (!detAndAction.Contains(":")) {
2894 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2898 Int_t colon = detAndAction.Index(":") ;
2899 fQADetectors = detAndAction(0, colon) ;
2900 if (fQADetectors.Contains("ALL") )
2901 fQADetectors = fFillESD ;
2902 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2903 if (fQATasks.Contains("ALL") ) {
2904 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2906 fQATasks.ToUpper() ;
2908 if ( fQATasks.Contains("RAW") )
2909 tempo = Form("%d ", AliQA::kRAWS) ;
2910 if ( fQATasks.Contains("RECPOINT") )
2911 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2912 if ( fQATasks.Contains("ESD") )
2913 tempo += Form("%d ", AliQA::kESDS) ;
2915 if (fQATasks.IsNull()) {
2916 AliInfo("No QA requested\n") ;
2921 TString tempo(fQATasks) ;
2922 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2923 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2924 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2925 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2930 //_____________________________________________________________________________
2931 Bool_t AliReconstruction::InitRecoParams()
2933 // The method accesses OCDB and retrieves all
2934 // the available reco-param objects from there.
2936 Bool_t isOK = kTRUE;
2938 TString detStr = fLoadCDB;
2939 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2941 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2943 if (fRecoParam.GetDetRecoParamArray(iDet)) {
2944 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
2948 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
2950 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
2951 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
2953 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
2957 TObject *recoParamObj = entry->GetObject();
2958 if (dynamic_cast<TObjArray*>(recoParamObj)) {
2959 // The detector has a normal TobjArray of AliDetectorRecoParam objects
2960 // Registering them in AliRecoParam
2961 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
2963 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
2964 // The detector has only onse set of reco parameters
2965 // Registering it in AliRecoParam
2966 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
2967 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
2968 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
2971 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
2975 AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
2979 if (AliDebugLevel() > 0) fRecoParam.Print();
2984 //_____________________________________________________________________________
2985 Bool_t AliReconstruction::GetEventInfo()
2987 // Fill the event info object
2989 AliCodeTimerAuto("")
2991 AliCentralTrigger *aCTP = NULL;
2993 fEventInfo.SetEventType(fRawReader->GetType());
2995 ULong64_t mask = fRawReader->GetClassMask();
2996 fEventInfo.SetTriggerMask(mask);
2997 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
2998 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3000 aCTP = new AliCentralTrigger();
3001 TString configstr("");
3002 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3003 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3007 aCTP->SetClassMask(mask);
3008 aCTP->SetClusterMask(clmask);
3011 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3013 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3014 aCTP = fRunLoader->GetTrigger();
3015 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3016 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3019 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3024 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3026 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3027 if (fRawReader) delete aCTP;
3031 UChar_t clustmask = 0;
3033 ULong64_t trmask = fEventInfo.GetTriggerMask();
3034 const TObjArray& classesArray = config->GetClasses();
3035 Int_t nclasses = classesArray.GetEntriesFast();
3036 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3037 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3039 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3040 fesd->SetTriggerClass(trclass->GetName(),trindex);
3041 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3042 if (trmask & (1 << trindex)) {
3044 trclasses += trclass->GetName();
3046 clustmask |= trclass->GetCluster()->GetClusterMask();
3050 fEventInfo.SetTriggerClasses(trclasses);
3052 // Set the information in ESD
3053 fesd->SetTriggerMask(trmask);
3054 fesd->SetTriggerCluster(clustmask);
3056 if (!aCTP->CheckTriggeredDetectors()) {
3057 if (fRawReader) delete aCTP;
3061 if (fRawReader) delete aCTP;
3063 // We have to fill also the HLT decision here!!
3069 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3071 // Match the detector list found in the rec.C or the default 'ALL'
3072 // to the list found in the GRP (stored there by the shuttle PP which
3073 // gets the information from ECS)
3074 static TString resultList;
3075 TString detList = detectorList;
3079 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3080 if ((detectorMask >> iDet) & 0x1) {
3081 TString det = AliDAQ::OfflineModuleName(iDet);
3082 if ((detList.CompareTo("ALL") == 0) ||
3083 ((detList.BeginsWith("ALL ") ||
3084 detList.EndsWith(" ALL") ||
3085 detList.Contains(" ALL ")) &&
3086 !(detList.BeginsWith("-"+det+" ") ||
3087 detList.EndsWith(" -"+det) ||
3088 detList.Contains(" -"+det+" "))) ||
3089 (detList.CompareTo(det) == 0) ||
3090 detList.BeginsWith(det+" ") ||
3091 detList.EndsWith(" "+det) ||
3092 detList.Contains( " "+det+" " )) {
3093 if (!resultList.EndsWith(det + " ")) {
3102 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3103 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3104 if ((detList.CompareTo("ALL") == 0) ||
3105 ((detList.BeginsWith("ALL ") ||
3106 detList.EndsWith(" ALL") ||
3107 detList.Contains(" ALL ")) &&
3108 !(detList.BeginsWith("-"+hltDet+" ") ||
3109 detList.EndsWith(" -"+hltDet) ||
3110 detList.Contains(" -"+hltDet+" "))) ||
3111 (detList.CompareTo(hltDet) == 0) ||
3112 detList.BeginsWith(hltDet+" ") ||
3113 detList.EndsWith(" "+hltDet) ||
3114 detList.Contains( " "+hltDet+" " )) {
3115 resultList += hltDet;
3119 return resultList.Data();
3123 //______________________________________________________________________________
3124 void AliReconstruction::Abort(const char *method, EAbort what)
3126 // Abort processing. If what = kAbortProcess, the Process() loop will be
3127 // aborted. If what = kAbortFile, the current file in a chain will be
3128 // aborted and the processing will continue with the next file, if there
3129 // is no next file then Process() will be aborted. Abort() can also be
3130 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3131 // the SlaveTerminate() and Terminate() are always called. The abort flag
3132 // can be checked in these methods using GetAbort().
3134 // The method is overwritten in AliReconstruction for better handling of
3135 // reco specific errors
3137 if (!fStopOnError) return;
3141 TString whyMess = method;
3142 whyMess += " failed! Aborting...";
3144 AliError(whyMess.Data());
3147 TString mess = "Abort";
3148 if (fAbort == kAbortProcess)
3149 mess = "AbortProcess";
3150 else if (fAbort == kAbortFile)
3153 Info(mess, whyMess.Data());