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::fgkNDetectors] = {"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(1),
236 fLoadAlignFromCDB(kTRUE),
237 fLoadAlignData("ALL"),
245 fParentRawReader(NULL),
250 fDiamondProfile(NULL),
251 fDiamondProfileTPC(NULL),
252 fMeanVertexConstraint(kTRUE),
256 fAlignObjArray(NULL),
259 fInitCDBCalled(kFALSE),
260 fSetRunNumberFromDataCalled(kFALSE),
266 fSameQACycle(kFALSE),
268 fRunPlaneEff(kFALSE),
277 fIsNewRunLoader(kFALSE),
281 // create reconstruction object with default parameters
284 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
285 fReconstructor[iDet] = NULL;
286 fLoader[iDet] = NULL;
287 fTracker[iDet] = NULL;
289 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
290 fQACycles[iDet] = 999999 ;
291 fQAWriteExpert[iDet] = kFALSE ;
297 //_____________________________________________________________________________
298 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
300 fUniformField(rec.fUniformField),
301 fForcedFieldMap(NULL),
302 fRunVertexFinder(rec.fRunVertexFinder),
303 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
304 fRunHLTTracking(rec.fRunHLTTracking),
305 fRunMuonTracking(rec.fRunMuonTracking),
306 fRunV0Finder(rec.fRunV0Finder),
307 fRunCascadeFinder(rec.fRunCascadeFinder),
308 fStopOnError(rec.fStopOnError),
309 fWriteAlignmentData(rec.fWriteAlignmentData),
310 fWriteESDfriend(rec.fWriteESDfriend),
311 fFillTriggerESD(rec.fFillTriggerESD),
313 fCleanESD(rec.fCleanESD),
314 fV0DCAmax(rec.fV0DCAmax),
315 fV0CsPmin(rec.fV0CsPmin),
319 fRunLocalReconstruction(rec.fRunLocalReconstruction),
320 fRunTracking(rec.fRunTracking),
321 fFillESD(rec.fFillESD),
322 fLoadCDB(rec.fLoadCDB),
323 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
324 fGAliceFileName(rec.fGAliceFileName),
325 fRawInput(rec.fRawInput),
326 fEquipIdMap(rec.fEquipIdMap),
327 fFirstEvent(rec.fFirstEvent),
328 fLastEvent(rec.fLastEvent),
329 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
331 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
332 fLoadAlignData(rec.fLoadAlignData),
333 fESDPar(rec.fESDPar),
334 fUseHLTData(rec.fUseHLTData),
340 fParentRawReader(NULL),
342 fRecoParam(rec.fRecoParam),
345 fDiamondProfile(rec.fDiamondProfile),
346 fDiamondProfileTPC(rec.fDiamondProfileTPC),
347 fMeanVertexConstraint(rec.fMeanVertexConstraint),
351 fAlignObjArray(rec.fAlignObjArray),
352 fCDBUri(rec.fCDBUri),
354 fInitCDBCalled(rec.fInitCDBCalled),
355 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
356 fQADetectors(rec.fQADetectors),
358 fQATasks(rec.fQATasks),
360 fRunGlobalQA(rec.fRunGlobalQA),
361 fSameQACycle(rec.fSameQACycle),
362 fRunPlaneEff(rec.fRunPlaneEff),
371 fIsNewRunLoader(rec.fIsNewRunLoader),
377 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
378 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
380 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
381 fReconstructor[iDet] = NULL;
382 fLoader[iDet] = NULL;
383 fTracker[iDet] = NULL;
386 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
387 fQACycles[iDet] = rec.fQACycles[iDet];
388 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
391 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
392 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
396 //_____________________________________________________________________________
397 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
399 // assignment operator
400 // Used in PROOF mode
401 // Be very careful while modifing it!
402 // Simple rules to follow:
403 // for persistent data members - use their assignment operators
404 // for non-persistent ones - do nothing or take the default values from constructor
405 // TSelector members should not be touched
406 if(&rec == this) return *this;
408 fUniformField = rec.fUniformField;
409 fForcedFieldMap = NULL;
410 fRunVertexFinder = rec.fRunVertexFinder;
411 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
412 fRunHLTTracking = rec.fRunHLTTracking;
413 fRunMuonTracking = rec.fRunMuonTracking;
414 fRunV0Finder = rec.fRunV0Finder;
415 fRunCascadeFinder = rec.fRunCascadeFinder;
416 fStopOnError = rec.fStopOnError;
417 fWriteAlignmentData = rec.fWriteAlignmentData;
418 fWriteESDfriend = rec.fWriteESDfriend;
419 fFillTriggerESD = rec.fFillTriggerESD;
421 fCleanESD = rec.fCleanESD;
422 fV0DCAmax = rec.fV0DCAmax;
423 fV0CsPmin = rec.fV0CsPmin;
427 fRunLocalReconstruction = rec.fRunLocalReconstruction;
428 fRunTracking = rec.fRunTracking;
429 fFillESD = rec.fFillESD;
430 fLoadCDB = rec.fLoadCDB;
431 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
432 fGAliceFileName = rec.fGAliceFileName;
433 fRawInput = rec.fRawInput;
434 fEquipIdMap = rec.fEquipIdMap;
435 fFirstEvent = rec.fFirstEvent;
436 fLastEvent = rec.fLastEvent;
437 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
439 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
440 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
443 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
444 fLoadAlignData = rec.fLoadAlignData;
445 fESDPar = rec.fESDPar;
446 fUseHLTData = rec.fUseHLTData;
448 delete fRunInfo; fRunInfo = NULL;
449 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
451 fEventInfo = rec.fEventInfo;
455 fParentRawReader = NULL;
457 fRecoParam = rec.fRecoParam;
459 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
460 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
461 delete fLoader[iDet]; fLoader[iDet] = NULL;
462 delete fTracker[iDet]; fTracker[iDet] = NULL;
465 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
466 fQACycles[iDet] = rec.fQACycles[iDet];
467 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
471 delete fDiamondProfile; fDiamondProfile = NULL;
472 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
473 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
474 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
475 fMeanVertexConstraint = rec.fMeanVertexConstraint;
477 delete fGRPData; fGRPData = NULL;
478 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
479 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
481 delete fAlignObjArray; fAlignObjArray = NULL;
484 fSpecCDBUri.Delete();
485 fInitCDBCalled = rec.fInitCDBCalled;
486 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
487 fQADetectors = rec.fQADetectors;
489 fQATasks = rec.fQATasks;
491 fRunGlobalQA = rec.fRunGlobalQA;
492 fSameQACycle = rec.fSameQACycle;
493 fRunPlaneEff = rec.fRunPlaneEff;
502 fIsNewRunLoader = rec.fIsNewRunLoader;
509 //_____________________________________________________________________________
510 AliReconstruction::~AliReconstruction()
516 delete fForcedFieldMap;
518 if (fAlignObjArray) {
519 fAlignObjArray->Delete();
520 delete fAlignObjArray;
522 fSpecCDBUri.Delete();
524 AliCodeTimer::Instance()->Print();
527 //_____________________________________________________________________________
528 void AliReconstruction::InitCDB()
530 // activate a default CDB storage
531 // First check if we have any CDB storage set, because it is used
532 // to retrieve the calibration and alignment constants
533 AliCodeTimerAuto("");
535 if (fInitCDBCalled) return;
536 fInitCDBCalled = kTRUE;
538 AliCDBManager* man = AliCDBManager::Instance();
539 if (man->IsDefaultStorageSet())
541 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
542 AliWarning("Default CDB storage has been already set !");
543 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
544 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
545 fCDBUri = man->GetDefaultStorage()->GetURI();
548 if (fCDBUri.Length() > 0)
550 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
551 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
552 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
554 fCDBUri="local://$ALICE_ROOT";
555 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
556 AliWarning("Default CDB storage not yet set !!!!");
557 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
558 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
561 man->SetDefaultStorage(fCDBUri);
564 // Now activate the detector specific CDB storage locations
565 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
566 TObject* obj = fSpecCDBUri[i];
568 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
569 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
570 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
571 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
573 AliSysInfo::AddStamp("InitCDB");
576 //_____________________________________________________________________________
577 void AliReconstruction::SetDefaultStorage(const char* uri) {
578 // Store the desired default CDB storage location
579 // Activate it later within the Run() method
585 //_____________________________________________________________________________
586 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
587 // Store a detector-specific CDB storage location
588 // Activate it later within the Run() method
590 AliCDBPath aPath(calibType);
591 if(!aPath.IsValid()){
592 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
593 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
594 if(!strcmp(calibType, fgkDetectorName[iDet])) {
595 aPath.SetPath(Form("%s/*", calibType));
596 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
600 if(!aPath.IsValid()){
601 AliError(Form("Not a valid path or detector: %s", calibType));
606 // // check that calibType refers to a "valid" detector name
607 // Bool_t isDetector = kFALSE;
608 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
609 // TString detName = fgkDetectorName[iDet];
610 // if(aPath.GetLevel0() == detName) {
611 // isDetector = kTRUE;
617 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
621 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
622 if (obj) fSpecCDBUri.Remove(obj);
623 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
627 //_____________________________________________________________________________
628 Bool_t AliReconstruction::SetRunNumberFromData()
630 // The method is called in Run() in order
631 // to set a correct run number.
632 // In case of raw data reconstruction the
633 // run number is taken from the raw data header
635 if (fSetRunNumberFromDataCalled) return kTRUE;
636 fSetRunNumberFromDataCalled = kTRUE;
638 AliCDBManager* man = AliCDBManager::Instance();
641 if(fRawReader->NextEvent()) {
642 if(man->GetRun() > 0) {
643 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
645 man->SetRun(fRawReader->GetRunNumber());
646 fRawReader->RewindEvents();
649 if(man->GetRun() > 0) {
650 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
653 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
659 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
661 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
666 // read run number from gAlice
667 if(rl->GetHeader()) {
668 man->SetRun(rl->GetHeader()->GetRun());
673 AliError("Neither run-loader header nor RawReader objects are found !");
685 //_____________________________________________________________________________
686 void AliReconstruction::SetCDBLock() {
687 // Set CDB lock: from now on it is forbidden to reset the run number
688 // or the default storage or to activate any further storage!
690 AliCDBManager::Instance()->SetLock(1);
693 //_____________________________________________________________________________
694 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
696 // Read the alignment objects from CDB.
697 // Each detector is supposed to have the
698 // alignment objects in DET/Align/Data CDB path.
699 // All the detector objects are then collected,
700 // sorted by geometry level (starting from ALIC) and
701 // then applied to the TGeo geometry.
702 // Finally an overlaps check is performed.
704 // Load alignment data from CDB and fill fAlignObjArray
705 if(fLoadAlignFromCDB){
707 TString detStr = detectors;
708 TString loadAlObjsListOfDets = "";
710 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
711 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
712 if(fgkDetectorName[iDet]=="HLT") continue;
713 if(AliGeomManager::IsModuleInGeom(fgkDetectorName[iDet]))
715 loadAlObjsListOfDets += fgkDetectorName[iDet];
716 loadAlObjsListOfDets += " ";
718 } // end loop over detectors
719 if(AliGeomManager::IsModuleInGeom("FRAME"))
720 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
721 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
722 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
724 // Check if the array with alignment objects was
725 // provided by the user. If yes, apply the objects
726 // to the present TGeo geometry
727 if (fAlignObjArray) {
728 if (gGeoManager && gGeoManager->IsClosed()) {
729 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
730 AliError("The misalignment of one or more volumes failed!"
731 "Compare the list of simulated detectors and the list of detector alignment data!");
736 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
742 if (fAlignObjArray) {
743 fAlignObjArray->Delete();
744 delete fAlignObjArray; fAlignObjArray=NULL;
750 //_____________________________________________________________________________
751 void AliReconstruction::SetGAliceFile(const char* fileName)
753 // set the name of the galice file
755 fGAliceFileName = fileName;
758 //_____________________________________________________________________________
759 void AliReconstruction::SetInput(const char* input)
761 // In case the input string starts with 'mem://', we run in an online mode
762 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
763 // file is assumed. One can give as an input:
764 // mem://: - events taken from DAQ monitoring libs online
766 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
767 if (input) fRawInput = input;
770 //_____________________________________________________________________________
771 void AliReconstruction::SetOption(const char* detector, const char* option)
773 // set options for the reconstruction of a detector
775 TObject* obj = fOptions.FindObject(detector);
776 if (obj) fOptions.Remove(obj);
777 fOptions.Add(new TNamed(detector, option));
780 //_____________________________________________________________________________
781 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
783 // Set custom reconstruction parameters for a given detector
784 // Single set of parameters for all the events
786 // First check if the reco-params are global
787 if(!strcmp(detector, "GRP")) {
789 fRecoParam.AddDetRecoParam(fgkNDetectors,par);
793 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
794 if(!strcmp(detector, fgkDetectorName[iDet])) {
796 fRecoParam.AddDetRecoParam(iDet,par);
803 //_____________________________________________________________________________
804 Bool_t AliReconstruction::SetFieldMap(Float_t l3Current, Float_t diCurrent, Float_t factor, const char *path) {
805 //------------------------------------------------
806 // The magnetic field map, defined externally...
807 // L3 current 30000 A -> 0.5 T
808 // L3 current 12000 A -> 0.2 T
809 // dipole current 6000 A
810 // The polarities must be the same
811 //------------------------------------------------
812 const Float_t l3NominalCurrent1=30000.; // (A)
813 const Float_t l3NominalCurrent2=12000.; // (A)
814 const Float_t diNominalCurrent =6000. ; // (A)
816 const Float_t tolerance=0.03; // relative current tolerance
817 const Float_t zero=77.; // "zero" current (A)
820 Bool_t dipoleON=kFALSE;
822 TString s=(factor < 0) ? "L3: -" : "L3: +";
824 l3Current = TMath::Abs(l3Current);
825 if (TMath::Abs(l3Current-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
826 map=AliMagWrapCheb::k5kG;
829 if (TMath::Abs(l3Current-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
830 map=AliMagWrapCheb::k2kG;
833 if (l3Current < zero) {
834 map=AliMagWrapCheb::k2kG;
836 factor=0.; // in fact, this is a global factor...
837 fUniformField=kTRUE; // track with the uniform (zero) B field
839 AliError(Form("Wrong L3 current (%f A)!",l3Current));
843 diCurrent = TMath::Abs(diCurrent);
844 if (TMath::Abs(diCurrent-diNominalCurrent)/diNominalCurrent < tolerance) {
845 // 3% current tolerance...
849 if (diCurrent < zero) { // some small current..
853 AliError(Form("Wrong dipole current (%f A)!",diCurrent));
857 delete fForcedFieldMap;
859 new AliMagWrapCheb("B field map ",s,2,factor,10.,map,dipoleON,path);
861 fForcedFieldMap->Print();
863 AliTracker::SetFieldMap(fForcedFieldMap,fUniformField);
869 Bool_t AliReconstruction::InitGRP() {
870 //------------------------------------
871 // Initialization of the GRP entry
872 //------------------------------------
873 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
877 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
880 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
882 fGRPData = new AliGRPObject();
883 fGRPData->ReadValuesFromMap(m);
887 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
888 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
892 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
896 AliError("No GRP entry found in OCDB!");
900 TString lhcState = fGRPData->GetLHCState();
901 if (lhcState==AliGRPObject::GetInvalidString()) {
902 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
903 lhcState = "UNKNOWN";
906 TString beamType = fGRPData->GetBeamType();
907 if (beamType==AliGRPObject::GetInvalidString()) {
908 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
909 beamType = "UNKNOWN";
912 Float_t beamEnergy = fGRPData->GetBeamEnergy();
913 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
914 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
918 TString runType = fGRPData->GetRunType();
919 if (runType==AliGRPObject::GetInvalidString()) {
920 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
924 Int_t activeDetectors = fGRPData->GetDetectorMask();
925 if (activeDetectors==AliGRPObject::GetInvalidInt()) {
926 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
927 activeDetectors = 1074790399;
930 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
935 // Process the list of active detectors
936 if (activeDetectors) {
937 UInt_t detMask = activeDetectors;
938 fLoadCDB.Form("%s %s %s %s",
939 fRunLocalReconstruction.Data(),
942 fQADetectors.Data());
943 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
944 fRunTracking = MatchDetectorList(fRunTracking,detMask);
945 fFillESD = MatchDetectorList(fFillESD,detMask);
946 fQADetectors = MatchDetectorList(fQADetectors,detMask);
947 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
950 AliInfo("===================================================================================");
951 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
952 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
953 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
954 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
955 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
956 AliInfo("===================================================================================");
958 //*** Dealing with the magnetic field map
959 if (AliTracker::GetFieldMap()) {
960 AliInfo("Running with the externally set B field !");
962 // Construct the field map out of the information retrieved from GRP.
967 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
968 if (l3Current == AliGRPObject::GetInvalidFloat()) {
969 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
973 Char_t l3Polarity = fGRPData->GetL3Polarity();
974 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
975 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
980 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
981 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
982 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
986 Char_t diPolarity = fGRPData->GetDipolePolarity();
987 if (diPolarity == AliGRPObject::GetInvalidChar()) {
988 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
993 TObjString *l3Current=
994 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
996 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
999 TObjString *l3Polarity=
1000 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
1002 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1007 TObjString *diCurrent=
1008 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
1010 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1013 TObjString *diPolarity=
1014 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
1016 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1022 Float_t l3Cur=TMath::Abs(l3Current);
1023 Float_t diCur=TMath::Abs(diCurrent);
1024 Float_t l3Pol=l3Polarity;
1025 // Float_t l3Cur=TMath::Abs(atof(l3Current->GetName()));
1026 //Float_t diCur=TMath::Abs(atof(diCurrent->GetName()));
1027 //Float_t l3Pol=atof(l3Polarity->GetName());
1029 if (l3Pol != 0.) factor=-1.;
1032 if (!SetFieldMap(l3Cur, diCur, factor)) {
1033 AliFatal("Failed to creat a B field map ! Exiting...");
1035 AliInfo("Running with the B field constructed out of GRP !");
1038 AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1043 //*** Get the diamond profile from OCDB
1044 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1046 if (fMeanVertexConstraint)
1047 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1049 AliError("No diamond profile found in OCDB!");
1052 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1054 if (fMeanVertexConstraint)
1055 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1057 AliError("No diamond profile found in OCDB!");
1063 //_____________________________________________________________________________
1064 Bool_t AliReconstruction::LoadCDB()
1066 AliCodeTimerAuto("");
1068 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1070 TString detStr = fLoadCDB;
1071 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1072 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1073 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1078 //_____________________________________________________________________________
1079 Bool_t AliReconstruction::Run(const char* input)
1082 AliCodeTimerAuto("");
1085 if (GetAbort() != TSelector::kContinue) return kFALSE;
1087 TChain *chain = NULL;
1088 if (fRawReader && (chain = fRawReader->GetChain())) {
1091 gProof->AddInput(this);
1093 outputFile.SetProtocol("root",kTRUE);
1094 outputFile.SetHost(gSystem->HostName());
1095 outputFile.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1096 AliInfo(Form("Output file with ESDs is %s",outputFile.GetUrl()));
1097 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.GetUrl()));
1099 chain->Process("AliReconstruction");
1102 chain->Process(this);
1107 if (GetAbort() != TSelector::kContinue) return kFALSE;
1109 if (GetAbort() != TSelector::kContinue) return kFALSE;
1110 //******* The loop over events
1111 AliInfo("Starting looping over events");
1113 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1114 (fRawReader && fRawReader->NextEvent())) {
1115 if (!ProcessEvent(iEvent)) {
1116 Abort("ProcessEvent",TSelector::kAbortFile);
1122 if (GetAbort() != TSelector::kContinue) return kFALSE;
1124 if (GetAbort() != TSelector::kContinue) return kFALSE;
1130 //_____________________________________________________________________________
1131 void AliReconstruction::InitRawReader(const char* input)
1133 AliCodeTimerAuto("");
1135 // Init raw-reader and
1136 // set the input in case of raw data
1137 if (input) fRawInput = input;
1138 fRawReader = AliRawReader::Create(fRawInput.Data());
1140 AliInfo("Reconstruction will run over digits");
1142 if (!fEquipIdMap.IsNull() && fRawReader)
1143 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1145 if (!fUseHLTData.IsNull()) {
1146 // create the RawReaderHLT which performs redirection of HLT input data for
1147 // the specified detectors
1148 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1150 fParentRawReader=fRawReader;
1151 fRawReader=pRawReader;
1153 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1156 AliSysInfo::AddStamp("CreateRawReader");
1159 //_____________________________________________________________________________
1160 void AliReconstruction::InitRun(const char* input)
1162 // Initialization of raw-reader,
1163 // run number, CDB etc.
1164 AliCodeTimerAuto("");
1165 AliSysInfo::AddStamp("Start");
1167 // Initialize raw-reader if any
1168 InitRawReader(input);
1170 // Initialize the CDB storage
1173 // Set run number in CDBManager (if it is not already set by the user)
1174 if (!SetRunNumberFromData()) {
1175 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1179 // Set CDB lock: from now on it is forbidden to reset the run number
1180 // or the default storage or to activate any further storage!
1185 //_____________________________________________________________________________
1186 void AliReconstruction::Begin(TTree *)
1188 // Initialize AlReconstruction before
1189 // going into the event loop
1190 // Should follow the TSelector convention
1191 // i.e. initialize only the object on the client side
1192 AliCodeTimerAuto("");
1194 AliReconstruction *reco = NULL;
1196 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1199 AliSysInfo::AddStamp("ReadInputInBegin");
1202 // Import ideal TGeo geometry and apply misalignment
1204 TString geom(gSystem->DirName(fGAliceFileName));
1205 geom += "/geometry.root";
1206 AliGeomManager::LoadGeometry(geom.Data());
1208 Abort("LoadGeometry", TSelector::kAbortProcess);
1211 AliSysInfo::AddStamp("LoadGeom");
1212 TString detsToCheck=fRunLocalReconstruction;
1213 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1214 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1217 AliSysInfo::AddStamp("CheckGeom");
1220 if (!MisalignGeometry(fLoadAlignData)) {
1221 Abort("MisalignGeometry", TSelector::kAbortProcess);
1224 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1225 AliSysInfo::AddStamp("MisalignGeom");
1228 Abort("InitGRP", TSelector::kAbortProcess);
1231 AliSysInfo::AddStamp("InitGRP");
1234 Abort("LoadCDB", TSelector::kAbortProcess);
1237 AliSysInfo::AddStamp("LoadCDB");
1239 // Read the reconstruction parameters from OCDB
1240 if (!InitRecoParams()) {
1241 AliWarning("Not all detectors have correct RecoParam objects initialized");
1243 AliSysInfo::AddStamp("InitRecoParams");
1246 if (reco) *reco = *this;
1247 fInput->Add(gGeoManager);
1249 fInput->Add(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()));
1250 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1251 AliMagF *magFieldMap = (AliMagF*)AliTracker::GetFieldMap();
1252 magFieldMap->SetName("MagneticFieldMap");
1253 fInput->Add(magFieldMap);
1258 //_____________________________________________________________________________
1259 void AliReconstruction::SlaveBegin(TTree*)
1261 // Initialization related to run-loader,
1262 // vertexer, trackers, recontructors
1263 // In proof mode it is executed on the slave
1264 AliCodeTimerAuto("");
1266 TProofOutputFile *outProofFile = NULL;
1268 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1271 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1273 AliGeomManager::SetGeometry(tgeo);
1275 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1276 Int_t runNumber = -1;
1277 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1278 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1279 man->SetCacheFlag(kTRUE);
1280 man->SetLock(kTRUE);
1284 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1285 AliTracker::SetFieldMap(map,fUniformField);
1287 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1288 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1289 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1290 fOutput->Add(outProofFile);
1292 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1295 // get the run loader
1296 if (!InitRunLoader()) {
1297 Abort("InitRunLoader", TSelector::kAbortProcess);
1300 AliSysInfo::AddStamp("LoadLoader");
1302 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1303 if(fDiamondProfile && fMeanVertexConstraint) ftVertexer->SetVtxStart(fDiamondProfile);
1306 if (fRunVertexFinder && !CreateVertexer()) {
1307 Abort("CreateVertexer", TSelector::kAbortProcess);
1310 AliSysInfo::AddStamp("CreateVertexer");
1313 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1314 Abort("CreateTrackers", TSelector::kAbortProcess);
1317 AliSysInfo::AddStamp("CreateTrackers");
1319 // create the ESD output file and tree
1320 if (!outProofFile) {
1321 ffile = TFile::Open("AliESDs.root", "RECREATE");
1322 ffile->SetCompressionLevel(2);
1323 if (!ffile->IsOpen()) {
1324 Abort("OpenESDFile", TSelector::kAbortProcess);
1329 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1330 Abort(Form("Problems opening output PROOF file: %s/%s",
1331 outProofFile->GetDir(), outProofFile->GetFileName()),
1332 TSelector::kAbortProcess);
1337 ftree = new TTree("esdTree", "Tree with ESD objects");
1338 fesd = new AliESDEvent();
1339 fesd->CreateStdContent();
1340 if (fWriteESDfriend) {
1341 fesdf = new AliESDfriend();
1342 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1343 br->SetFile("AliESDfriends.root");
1344 fesd->AddObject(fesdf);
1346 fesd->WriteToTree(ftree);
1347 ftree->GetUserInfo()->Add(fesd);
1349 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1350 fhltesd = new AliESDEvent();
1351 fhltesd->CreateStdContent();
1352 fhltesd->WriteToTree(fhlttree);
1353 fhlttree->GetUserInfo()->Add(fhltesd);
1355 ProcInfo_t ProcInfo;
1356 gSystem->GetProcInfo(&ProcInfo);
1357 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1360 //Initialize the QA and start of cycle
1362 fQASteer = new AliQADataMakerSteer("rec") ;
1363 fQASteer->SetActiveDetectors(fQADetectors) ;
1364 for (Int_t det = 0 ; det < AliQA::kNDET ; det++) {
1365 fQASteer->SetCycleLength(AliQA::DETECTORINDEX_t(det), fQACycles[det]) ;
1366 if (fQAWriteExpert[det])
1367 fQASteer->SetWriteExpert(AliQA::DETECTORINDEX_t(det)) ;
1370 if (!fRawReader && fQATasks.Contains(AliQA::kRAWS))
1371 fQATasks.ReplaceAll(Form("%d",AliQA::kRAWS), "") ;
1372 fQASteer->SetTasks(fQATasks) ;
1373 fQASteer->InitQADataMaker(AliCDBManager::Instance()->GetRun(), fRecoParam) ;
1377 Bool_t sameCycle = kFALSE ;
1378 if (!fQASteer) fQASteer = new AliQADataMakerSteer("rec") ;
1379 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1380 AliInfo(Form("Initializing the global QA data maker"));
1381 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1382 qadm->StartOfCycle(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1383 TObjArray *arr=qadm->Init(AliQA::kRECPOINTS);
1384 AliTracker::SetResidualsArray(arr);
1387 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1388 qadm->StartOfCycle(AliQA::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1389 qadm->Init(AliQA::kESDS);
1393 //Initialize the Plane Efficiency framework
1394 if (fRunPlaneEff && !InitPlaneEff()) {
1395 Abort("InitPlaneEff", TSelector::kAbortProcess);
1399 if (strcmp(gProgName,"alieve") == 0)
1400 fRunAliEVE = InitAliEVE();
1405 //_____________________________________________________________________________
1406 Bool_t AliReconstruction::Process(Long64_t entry)
1408 // run the reconstruction over a single entry
1409 // from the chain with raw data
1410 AliCodeTimerAuto("");
1412 TTree *currTree = fChain->GetTree();
1413 AliRawEvent *event = new AliRawEvent;
1414 currTree->SetBranchAddress("rawevent",&event);
1415 currTree->GetEntry(entry);
1416 fRawReader = new AliRawReaderRoot(event);
1417 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1425 //_____________________________________________________________________________
1426 void AliReconstruction::Init(TTree *tree)
1429 AliError("The input tree is not found!");
1435 //_____________________________________________________________________________
1436 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1438 // run the reconstruction over a single event
1439 // The event loop is steered in Run method
1441 AliCodeTimerAuto("");
1443 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1444 fRunLoader->SetEventNumber(iEvent);
1445 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1447 fRunLoader->TreeE()->Fill();
1448 if (fRawReader && fRawReader->UseAutoSaveESD())
1449 fRunLoader->TreeE()->AutoSave("SaveSelf");
1452 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1456 AliInfo(Form("processing event %d", iEvent));
1458 // Fill Event-info object
1460 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo);
1462 // Set the reco-params
1464 TString detStr = fLoadCDB;
1465 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1466 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1467 AliReconstructor *reconstructor = GetReconstructor(iDet);
1468 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1469 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1470 reconstructor->SetRecoParam(par);
1475 fRunLoader->GetEvent(iEvent);
1479 fQASteer->RunOneEvent(fRawReader) ;
1481 // local single event reconstruction
1482 if (!fRunLocalReconstruction.IsNull()) {
1483 TString detectors=fRunLocalReconstruction;
1484 // run HLT event reconstruction first
1485 // ;-( IsSelected changes the string
1486 if (IsSelected("HLT", detectors) &&
1487 !RunLocalEventReconstruction("HLT")) {
1488 if (fStopOnError) {CleanUp(); return kFALSE;}
1490 detectors=fRunLocalReconstruction;
1491 detectors.ReplaceAll("HLT", "");
1492 if (!RunLocalEventReconstruction(detectors)) {
1493 if (fStopOnError) {CleanUp(); return kFALSE;}
1497 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1498 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1499 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1500 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1502 // Set magnetic field from the tracker
1503 fesd->SetMagneticField(AliTracker::GetBz());
1504 fhltesd->SetMagneticField(AliTracker::GetBz());
1506 // Set most probable pt, for B=0 tracking
1507 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(fgkNDetectors));
1508 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1510 // Fill raw-data error log into the ESD
1511 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1514 if (fRunVertexFinder) {
1515 if (!RunVertexFinder(fesd)) {
1516 if (fStopOnError) {CleanUp(); return kFALSE;}
1521 if (!fRunTracking.IsNull()) {
1522 if (fRunMuonTracking) {
1523 if (!RunMuonTracking(fesd)) {
1524 if (fStopOnError) {CleanUp(); return kFALSE;}
1530 if (!fRunTracking.IsNull()) {
1531 if (!RunTracking(fesd)) {
1532 if (fStopOnError) {CleanUp(); return kFALSE;}
1537 if (!fFillESD.IsNull()) {
1538 TString detectors=fFillESD;
1539 // run HLT first and on hltesd
1540 // ;-( IsSelected changes the string
1541 if (IsSelected("HLT", detectors) &&
1542 !FillESD(fhltesd, "HLT")) {
1543 if (fStopOnError) {CleanUp(); return kFALSE;}
1546 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1547 if (detectors.Contains("ALL")) {
1549 for (Int_t idet=0; idet<fgkNDetectors; ++idet){
1550 detectors += fgkDetectorName[idet];
1554 detectors.ReplaceAll("HLT", "");
1555 if (!FillESD(fesd, detectors)) {
1556 if (fStopOnError) {CleanUp(); return kFALSE;}
1560 // fill Event header information from the RawEventHeader
1561 if (fRawReader){FillRawEventHeaderESD(fesd);}
1564 AliESDpid::MakePID(fesd);
1566 if (fFillTriggerESD) {
1567 if (!FillTriggerESD(fesd)) {
1568 if (fStopOnError) {CleanUp(); return kFALSE;}
1575 // Propagate track to the beam pipe (if not already done by ITS)
1577 const Int_t ntracks = fesd->GetNumberOfTracks();
1578 const Double_t kBz = fesd->GetMagneticField();
1579 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1582 UShort_t *selectedIdx=new UShort_t[ntracks];
1584 for (Int_t itrack=0; itrack<ntracks; itrack++){
1585 const Double_t kMaxStep = 5; //max step over the material
1588 AliESDtrack *track = fesd->GetTrack(itrack);
1589 if (!track) continue;
1591 AliExternalTrackParam *tpcTrack =
1592 (AliExternalTrackParam *)track->GetTPCInnerParam();
1596 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
1599 Int_t n=trkArray.GetEntriesFast();
1600 selectedIdx[n]=track->GetID();
1601 trkArray.AddLast(tpcTrack);
1604 //Tracks refitted by ITS should already be at the SPD vertex
1605 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1608 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1609 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1614 // Improve the reconstructed primary vertex position using the tracks
1616 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1617 if(fesd->GetPrimaryVertexSPD()) {
1618 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1619 if(vtitle.Contains("cosmics")) {
1620 runVertexFinderTracks=kFALSE;
1624 if (runVertexFinderTracks) {
1625 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1626 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(fgkNDetectors));
1628 // TPC + ITS primary vertex
1629 ftVertexer->SetITSMode();
1630 // get cuts for vertexer from AliGRPRecoParam
1632 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1633 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1634 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1635 ftVertexer->SetCuts(cutsVertexer);
1636 delete [] cutsVertexer; cutsVertexer = NULL;
1638 if(fDiamondProfile && fMeanVertexConstraint) {
1639 ftVertexer->SetVtxStart(fDiamondProfile);
1641 ftVertexer->SetConstraintOff();
1643 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1645 if (pvtx->GetStatus()) {
1646 fesd->SetPrimaryVertex(pvtx);
1647 for (Int_t i=0; i<ntracks; i++) {
1648 AliESDtrack *t = fesd->GetTrack(i);
1649 t->RelateToVertex(pvtx, kBz, kVeryBig);
1654 // TPC-only primary vertex
1655 ftVertexer->SetTPCMode();
1656 // get cuts for vertexer from AliGRPRecoParam
1658 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1659 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1660 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1661 ftVertexer->SetCuts(cutsVertexer);
1662 delete [] cutsVertexer; cutsVertexer = NULL;
1664 if(fDiamondProfileTPC && fMeanVertexConstraint) {
1665 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1667 ftVertexer->SetConstraintOff();
1669 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1671 if (pvtx->GetStatus()) {
1672 fesd->SetPrimaryVertexTPC(pvtx);
1673 for (Int_t i=0; i<ntracks; i++) {
1674 AliESDtrack *t = fesd->GetTrack(i);
1675 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1681 delete[] selectedIdx;
1683 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1688 AliV0vertexer vtxer;
1689 vtxer.Tracks2V0vertices(fesd);
1691 if (fRunCascadeFinder) {
1693 AliCascadeVertexer cvtxer;
1694 cvtxer.V0sTracks2CascadeVertices(fesd);
1699 if (fCleanESD) CleanESD(fesd);
1702 fQASteer->RunOneEvent(fesd) ;
1705 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1706 if (qadm && fQATasks.Contains(Form("%d", AliQA::kESDS)))
1707 qadm->Exec(AliQA::kESDS, fesd);
1710 if (fWriteESDfriend) {
1711 // fesdf->~AliESDfriend();
1712 // new (fesdf) AliESDfriend(); // Reset...
1713 fesd->GetESDfriend(fesdf);
1717 // Auto-save the ESD tree in case of prompt reco @P2
1718 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1719 ftree->AutoSave("SaveSelf");
1720 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1721 if (friendfile) friendfile->Save();
1728 if (fRunAliEVE) RunAliEVE();
1732 if (fWriteESDfriend) {
1733 fesdf->~AliESDfriend();
1734 new (fesdf) AliESDfriend(); // Reset...
1737 ProcInfo_t ProcInfo;
1738 gSystem->GetProcInfo(&ProcInfo);
1739 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1742 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1743 if (fReconstructor[iDet])
1744 fReconstructor[iDet]->SetRecoParam(NULL);
1747 if (fRunQA || fRunGlobalQA)
1748 fQASteer->Increment() ;
1753 //_____________________________________________________________________________
1754 void AliReconstruction::SlaveTerminate()
1756 // Finalize the run on the slave side
1757 // Called after the exit
1758 // from the event loop
1759 AliCodeTimerAuto("");
1761 if (fIsNewRunLoader) { // galice.root didn't exist
1762 fRunLoader->WriteHeader("OVERWRITE");
1763 fRunLoader->CdGAFile();
1764 fRunLoader->Write(0, TObject::kOverwrite);
1767 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1768 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1770 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1771 cdbMapCopy->SetOwner(1);
1772 cdbMapCopy->SetName("cdbMap");
1773 TIter iter(cdbMap->GetTable());
1776 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1777 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1778 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1779 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1782 TList *cdbListCopy = new TList();
1783 cdbListCopy->SetOwner(1);
1784 cdbListCopy->SetName("cdbList");
1786 TIter iter2(cdbList);
1789 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1790 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1793 ftree->GetUserInfo()->Add(cdbMapCopy);
1794 ftree->GetUserInfo()->Add(cdbListCopy);
1797 if(fESDPar.Contains("ESD.par")){
1798 AliInfo("Attaching ESD.par to Tree");
1799 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1800 ftree->GetUserInfo()->Add(fn);
1806 if (fWriteESDfriend)
1807 ftree->SetBranchStatus("ESDfriend*",0);
1808 // we want to have only one tree version number
1809 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1812 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1813 if (fRunPlaneEff && !FinishPlaneEff()) {
1814 AliWarning("Finish PlaneEff evaluation failed");
1817 // End of cycle for the in-loop
1819 fQASteer->EndOfCycle() ;
1822 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1824 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1825 qadm->EndOfCycle(AliQA::kRECPOINTS);
1826 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1827 qadm->EndOfCycle(AliQA::kESDS);
1835 //_____________________________________________________________________________
1836 void AliReconstruction::Terminate()
1838 // Create tags for the events in the ESD tree (the ESD tree is always present)
1839 // In case of empty events the tags will contain dummy values
1840 AliCodeTimerAuto("");
1842 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1843 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1845 // Cleanup of CDB manager: cache and active storages!
1846 AliCDBManager::Instance()->ClearCache();
1849 //_____________________________________________________________________________
1850 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1852 // run the local reconstruction
1854 static Int_t eventNr=0;
1855 AliCodeTimerAuto("")
1857 TString detStr = detectors;
1858 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1859 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1860 AliReconstructor* reconstructor = GetReconstructor(iDet);
1861 if (!reconstructor) continue;
1862 AliLoader* loader = fLoader[iDet];
1863 // Matthias April 2008: temporary fix to run HLT reconstruction
1864 // although the HLT loader is missing
1865 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1867 reconstructor->Reconstruct(fRawReader, NULL);
1870 reconstructor->Reconstruct(dummy, NULL);
1875 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1878 // conversion of digits
1879 if (fRawReader && reconstructor->HasDigitConversion()) {
1880 AliInfo(Form("converting raw data digits into root objects for %s",
1881 fgkDetectorName[iDet]));
1882 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1883 // fgkDetectorName[iDet]));
1884 loader->LoadDigits("update");
1885 loader->CleanDigits();
1886 loader->MakeDigitsContainer();
1887 TTree* digitsTree = loader->TreeD();
1888 reconstructor->ConvertDigits(fRawReader, digitsTree);
1889 loader->WriteDigits("OVERWRITE");
1890 loader->UnloadDigits();
1892 // local reconstruction
1893 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1894 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1895 loader->LoadRecPoints("update");
1896 loader->CleanRecPoints();
1897 loader->MakeRecPointsContainer();
1898 TTree* clustersTree = loader->TreeR();
1899 if (fRawReader && !reconstructor->HasDigitConversion()) {
1900 reconstructor->Reconstruct(fRawReader, clustersTree);
1902 loader->LoadDigits("read");
1903 TTree* digitsTree = loader->TreeD();
1905 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1906 if (fStopOnError) return kFALSE;
1908 reconstructor->Reconstruct(digitsTree, clustersTree);
1910 loader->UnloadDigits();
1913 TString detQAStr(fQADetectors) ;
1915 fQASteer->RunOneEventInOneDetector(iDet, clustersTree) ;
1917 loader->WriteRecPoints("OVERWRITE");
1918 loader->UnloadRecPoints();
1919 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1921 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1922 AliError(Form("the following detectors were not found: %s",
1924 if (fStopOnError) return kFALSE;
1930 //_____________________________________________________________________________
1931 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1933 // run the barrel tracking
1935 AliCodeTimerAuto("")
1937 AliESDVertex* vertex = NULL;
1938 Double_t vtxPos[3] = {0, 0, 0};
1939 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1940 TArrayF mcVertex(3);
1941 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1942 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1943 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1947 AliInfo("running the ITS vertex finder");
1949 fLoader[0]->LoadRecPoints();
1950 TTree* cltree = fLoader[0]->TreeR();
1952 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1953 vertex = fVertexer->FindVertexForCurrentEvent(cltree);
1956 AliError("Can't get the ITS cluster tree");
1958 fLoader[0]->UnloadRecPoints();
1961 AliError("Can't get the ITS loader");
1964 AliWarning("Vertex not found");
1965 vertex = new AliESDVertex();
1966 vertex->SetName("default");
1969 vertex->SetName("reconstructed");
1973 AliInfo("getting the primary vertex from MC");
1974 vertex = new AliESDVertex(vtxPos, vtxErr);
1978 vertex->GetXYZ(vtxPos);
1979 vertex->GetSigmaXYZ(vtxErr);
1981 AliWarning("no vertex reconstructed");
1982 vertex = new AliESDVertex(vtxPos, vtxErr);
1984 esd->SetPrimaryVertexSPD(vertex);
1985 // if SPD multiplicity has been determined, it is stored in the ESD
1986 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1987 if(mult)esd->SetMultiplicity(mult);
1989 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1990 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1997 //_____________________________________________________________________________
1998 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2000 // run the HLT barrel tracking
2002 AliCodeTimerAuto("")
2005 AliError("Missing runLoader!");
2009 AliInfo("running HLT tracking");
2011 // Get a pointer to the HLT reconstructor
2012 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
2013 if (!reconstructor) return kFALSE;
2016 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2017 TString detName = fgkDetectorName[iDet];
2018 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2019 reconstructor->SetOption(detName.Data());
2020 AliTracker *tracker = reconstructor->CreateTracker();
2022 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2023 if (fStopOnError) return kFALSE;
2027 Double_t vtxErr[3]={0.005,0.005,0.010};
2028 const AliESDVertex *vertex = esd->GetVertex();
2029 vertex->GetXYZ(vtxPos);
2030 tracker->SetVertex(vtxPos,vtxErr);
2032 fLoader[iDet]->LoadRecPoints("read");
2033 TTree* tree = fLoader[iDet]->TreeR();
2035 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2038 tracker->LoadClusters(tree);
2040 if (tracker->Clusters2Tracks(esd) != 0) {
2041 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2045 tracker->UnloadClusters();
2053 //_____________________________________________________________________________
2054 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2056 // run the muon spectrometer tracking
2058 AliCodeTimerAuto("")
2061 AliError("Missing runLoader!");
2064 Int_t iDet = 7; // for MUON
2066 AliInfo("is running...");
2068 // Get a pointer to the MUON reconstructor
2069 AliReconstructor *reconstructor = GetReconstructor(iDet);
2070 if (!reconstructor) return kFALSE;
2073 TString detName = fgkDetectorName[iDet];
2074 AliDebug(1, Form("%s tracking", detName.Data()));
2075 AliTracker *tracker = reconstructor->CreateTracker();
2077 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2082 fLoader[iDet]->LoadRecPoints("read");
2084 tracker->LoadClusters(fLoader[iDet]->TreeR());
2086 Int_t rv = tracker->Clusters2Tracks(esd);
2090 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2094 fLoader[iDet]->UnloadRecPoints();
2096 tracker->UnloadClusters();
2104 //_____________________________________________________________________________
2105 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2107 // run the barrel tracking
2108 static Int_t eventNr=0;
2109 AliCodeTimerAuto("")
2111 AliInfo("running tracking");
2113 //Fill the ESD with the T0 info (will be used by the TOF)
2114 if (fReconstructor[11] && fLoader[11]) {
2115 fLoader[11]->LoadRecPoints("READ");
2116 TTree *treeR = fLoader[11]->TreeR();
2117 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2120 // pass 1: TPC + ITS inwards
2121 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2122 if (!fTracker[iDet]) continue;
2123 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2126 fLoader[iDet]->LoadRecPoints("read");
2127 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2128 TTree* tree = fLoader[iDet]->TreeR();
2130 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2133 fTracker[iDet]->LoadClusters(tree);
2134 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2136 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2137 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2140 // preliminary PID in TPC needed by the ITS tracker
2142 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2143 AliESDpid::MakePID(esd);
2145 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2148 // pass 2: ALL backwards
2150 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2151 if (!fTracker[iDet]) continue;
2152 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2155 if (iDet > 1) { // all except ITS, TPC
2157 fLoader[iDet]->LoadRecPoints("read");
2158 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2159 tree = fLoader[iDet]->TreeR();
2161 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2164 fTracker[iDet]->LoadClusters(tree);
2165 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2169 if (iDet>1) // start filling residuals for the "outer" detectors
2170 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2172 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2173 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2178 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2179 fTracker[iDet]->UnloadClusters();
2180 fLoader[iDet]->UnloadRecPoints();
2182 // updated PID in TPC needed by the ITS tracker -MI
2184 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2185 AliESDpid::MakePID(esd);
2187 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2189 //stop filling residuals for the "outer" detectors
2190 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2192 // pass 3: TRD + TPC + ITS refit inwards
2194 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2195 if (!fTracker[iDet]) continue;
2196 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2199 if (iDet<2) // start filling residuals for TPC and ITS
2200 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2202 if (fTracker[iDet]->RefitInward(esd) != 0) {
2203 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2206 // run postprocessing
2207 if (fTracker[iDet]->PostProcess(esd) != 0) {
2208 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2211 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2214 // write space-points to the ESD in case alignment data output
2216 if (fWriteAlignmentData)
2217 WriteAlignmentData(esd);
2219 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2220 if (!fTracker[iDet]) continue;
2222 fTracker[iDet]->UnloadClusters();
2223 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2224 fLoader[iDet]->UnloadRecPoints();
2225 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2227 // stop filling residuals for TPC and ITS
2228 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2234 //_____________________________________________________________________________
2235 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2237 // Remove the data which are not needed for the physics analysis.
2240 Int_t nTracks=esd->GetNumberOfTracks();
2241 Int_t nV0s=esd->GetNumberOfV0s();
2243 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2245 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2246 Bool_t rc=esd->Clean(cleanPars);
2248 nTracks=esd->GetNumberOfTracks();
2249 nV0s=esd->GetNumberOfV0s();
2251 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2256 //_____________________________________________________________________________
2257 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2259 // fill the event summary data
2261 AliCodeTimerAuto("")
2262 static Int_t eventNr=0;
2263 TString detStr = detectors;
2265 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2266 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2267 AliReconstructor* reconstructor = GetReconstructor(iDet);
2268 if (!reconstructor) continue;
2269 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2270 TTree* clustersTree = NULL;
2271 if (fLoader[iDet]) {
2272 fLoader[iDet]->LoadRecPoints("read");
2273 clustersTree = fLoader[iDet]->TreeR();
2274 if (!clustersTree) {
2275 AliError(Form("Can't get the %s clusters tree",
2276 fgkDetectorName[iDet]));
2277 if (fStopOnError) return kFALSE;
2280 if (fRawReader && !reconstructor->HasDigitConversion()) {
2281 reconstructor->FillESD(fRawReader, clustersTree, esd);
2283 TTree* digitsTree = NULL;
2284 if (fLoader[iDet]) {
2285 fLoader[iDet]->LoadDigits("read");
2286 digitsTree = fLoader[iDet]->TreeD();
2288 AliError(Form("Can't get the %s digits tree",
2289 fgkDetectorName[iDet]));
2290 if (fStopOnError) return kFALSE;
2293 reconstructor->FillESD(digitsTree, clustersTree, esd);
2294 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2296 if (fLoader[iDet]) {
2297 fLoader[iDet]->UnloadRecPoints();
2301 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2302 AliError(Form("the following detectors were not found: %s",
2304 if (fStopOnError) return kFALSE;
2306 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2311 //_____________________________________________________________________________
2312 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2314 // Reads the trigger decision which is
2315 // stored in Trigger.root file and fills
2316 // the corresponding esd entries
2318 AliCodeTimerAuto("")
2320 AliInfo("Filling trigger information into the ESD");
2323 AliCTPRawStream input(fRawReader);
2324 if (!input.Next()) {
2325 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2328 if (esd->GetTriggerMask() != input.GetClassMask())
2329 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2330 input.GetClassMask(),esd->GetTriggerMask()));
2331 if (esd->GetOrbitNumber() != input.GetOrbitID())
2332 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2333 input.GetOrbitID(),esd->GetOrbitNumber()));
2334 if (esd->GetBunchCrossNumber() != input.GetBCID())
2335 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2336 input.GetBCID(),esd->GetBunchCrossNumber()));
2339 // Here one has to add the filling of trigger inputs and
2340 // interaction records
2350 //_____________________________________________________________________________
2351 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2354 // Filling information from RawReader Header
2357 if (!fRawReader) return kFALSE;
2359 AliInfo("Filling information from RawReader Header");
2361 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2362 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2363 esd->SetPeriodNumber(fRawReader->GetPeriod());
2365 esd->SetTimeStamp(fRawReader->GetTimestamp());
2366 esd->SetEventType(fRawReader->GetType());
2372 //_____________________________________________________________________________
2373 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2375 // check whether detName is contained in detectors
2376 // if yes, it is removed from detectors
2378 // check if all detectors are selected
2379 if ((detectors.CompareTo("ALL") == 0) ||
2380 detectors.BeginsWith("ALL ") ||
2381 detectors.EndsWith(" ALL") ||
2382 detectors.Contains(" ALL ")) {
2387 // search for the given detector
2388 Bool_t result = kFALSE;
2389 if ((detectors.CompareTo(detName) == 0) ||
2390 detectors.BeginsWith(detName+" ") ||
2391 detectors.EndsWith(" "+detName) ||
2392 detectors.Contains(" "+detName+" ")) {
2393 detectors.ReplaceAll(detName, "");
2397 // clean up the detectors string
2398 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2399 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2400 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2405 //_____________________________________________________________________________
2406 Bool_t AliReconstruction::InitRunLoader()
2408 // get or create the run loader
2410 if (gAlice) delete gAlice;
2413 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2414 // load all base libraries to get the loader classes
2415 TString libs = gSystem->GetLibraries();
2416 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2417 TString detName = fgkDetectorName[iDet];
2418 if (detName == "HLT") continue;
2419 if (libs.Contains("lib" + detName + "base.so")) continue;
2420 gSystem->Load("lib" + detName + "base.so");
2422 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2424 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2429 fRunLoader->CdGAFile();
2430 fRunLoader->LoadgAlice();
2432 //PH This is a temporary fix to give access to the kinematics
2433 //PH that is needed for the labels of ITS clusters
2434 fRunLoader->LoadHeader();
2435 fRunLoader->LoadKinematics();
2437 } else { // galice.root does not exist
2439 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2441 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2442 AliConfig::GetDefaultEventFolderName(),
2445 AliError(Form("could not create run loader in file %s",
2446 fGAliceFileName.Data()));
2450 fIsNewRunLoader = kTRUE;
2451 fRunLoader->MakeTree("E");
2453 if (fNumberOfEventsPerFile > 0)
2454 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2456 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2462 //_____________________________________________________________________________
2463 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2465 // get the reconstructor object and the loader for a detector
2467 if (fReconstructor[iDet]) {
2468 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2469 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2470 fReconstructor[iDet]->SetRecoParam(par);
2472 return fReconstructor[iDet];
2475 // load the reconstructor object
2476 TPluginManager* pluginManager = gROOT->GetPluginManager();
2477 TString detName = fgkDetectorName[iDet];
2478 TString recName = "Ali" + detName + "Reconstructor";
2480 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2482 AliReconstructor* reconstructor = NULL;
2483 // first check if a plugin is defined for the reconstructor
2484 TPluginHandler* pluginHandler =
2485 pluginManager->FindHandler("AliReconstructor", detName);
2486 // if not, add a plugin for it
2487 if (!pluginHandler) {
2488 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2489 TString libs = gSystem->GetLibraries();
2490 if (libs.Contains("lib" + detName + "base.so") ||
2491 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2492 pluginManager->AddHandler("AliReconstructor", detName,
2493 recName, detName + "rec", recName + "()");
2495 pluginManager->AddHandler("AliReconstructor", detName,
2496 recName, detName, recName + "()");
2498 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2500 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2501 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2503 if (reconstructor) {
2504 TObject* obj = fOptions.FindObject(detName.Data());
2505 if (obj) reconstructor->SetOption(obj->GetTitle());
2506 reconstructor->Init();
2507 fReconstructor[iDet] = reconstructor;
2510 // get or create the loader
2511 if (detName != "HLT") {
2512 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2513 if (!fLoader[iDet]) {
2514 AliConfig::Instance()
2515 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2517 // first check if a plugin is defined for the loader
2519 pluginManager->FindHandler("AliLoader", detName);
2520 // if not, add a plugin for it
2521 if (!pluginHandler) {
2522 TString loaderName = "Ali" + detName + "Loader";
2523 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2524 pluginManager->AddHandler("AliLoader", detName,
2525 loaderName, detName + "base",
2526 loaderName + "(const char*, TFolder*)");
2527 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2529 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2531 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2532 fRunLoader->GetEventFolder());
2534 if (!fLoader[iDet]) { // use default loader
2535 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2537 if (!fLoader[iDet]) {
2538 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2539 if (fStopOnError) return NULL;
2541 fRunLoader->AddLoader(fLoader[iDet]);
2542 fRunLoader->CdGAFile();
2543 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2544 fRunLoader->Write(0, TObject::kOverwrite);
2549 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2550 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2551 reconstructor->SetRecoParam(par);
2553 return reconstructor;
2556 //_____________________________________________________________________________
2557 Bool_t AliReconstruction::CreateVertexer()
2559 // create the vertexer
2562 AliReconstructor* itsReconstructor = GetReconstructor(0);
2563 if (itsReconstructor) {
2564 fVertexer = itsReconstructor->CreateVertexer();
2567 AliWarning("couldn't create a vertexer for ITS");
2568 if (fStopOnError) return kFALSE;
2574 //_____________________________________________________________________________
2575 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2577 // create the trackers
2578 AliInfo("Creating trackers");
2580 TString detStr = detectors;
2581 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2582 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2583 AliReconstructor* reconstructor = GetReconstructor(iDet);
2584 if (!reconstructor) continue;
2585 TString detName = fgkDetectorName[iDet];
2586 if (detName == "HLT") {
2587 fRunHLTTracking = kTRUE;
2590 if (detName == "MUON") {
2591 fRunMuonTracking = kTRUE;
2596 fTracker[iDet] = reconstructor->CreateTracker();
2597 if (!fTracker[iDet] && (iDet < 7)) {
2598 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2599 if (fStopOnError) return kFALSE;
2601 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2607 //_____________________________________________________________________________
2608 void AliReconstruction::CleanUp()
2610 // delete trackers and the run loader and close and delete the file
2612 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2613 delete fReconstructor[iDet];
2614 fReconstructor[iDet] = NULL;
2615 fLoader[iDet] = NULL;
2616 delete fTracker[iDet];
2617 fTracker[iDet] = NULL;
2628 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2629 delete fDiamondProfile;
2630 fDiamondProfile = NULL;
2631 delete fDiamondProfileTPC;
2632 fDiamondProfileTPC = NULL;
2638 delete fParentRawReader;
2639 fParentRawReader=NULL;
2648 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2650 // Write space-points which are then used in the alignment procedures
2651 // For the moment only ITS, TPC, TRD and TOF
2653 Int_t ntracks = esd->GetNumberOfTracks();
2654 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2656 AliESDtrack *track = esd->GetTrack(itrack);
2659 for (Int_t iDet = 3; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2660 nsp += track->GetNcls(iDet);
2662 if (iDet==0) { // ITS "extra" clusters
2663 track->GetClusters(iDet,idx);
2664 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2669 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2670 track->SetTrackPointArray(sp);
2672 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2673 AliTracker *tracker = fTracker[iDet];
2674 if (!tracker) continue;
2675 Int_t nspdet = track->GetClusters(iDet,idx);
2677 if (iDet==0) // ITS "extra" clusters
2678 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2680 if (nspdet <= 0) continue;
2684 while (isp2 < nspdet) {
2685 Bool_t isvalid=kTRUE;
2687 Int_t index=idx[isp++];
2688 if (index < 0) continue;
2690 TString dets = fgkDetectorName[iDet];
2691 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2692 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2693 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2694 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2695 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2697 isvalid = tracker->GetTrackPoint(index,p);
2700 if (!isvalid) continue;
2701 sp->AddPoint(isptrack,&p); isptrack++;
2708 //_____________________________________________________________________________
2709 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2711 // The method reads the raw-data error log
2712 // accumulated within the rawReader.
2713 // It extracts the raw-data errors related to
2714 // the current event and stores them into
2715 // a TClonesArray inside the esd object.
2717 if (!fRawReader) return;
2719 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2721 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2723 if (iEvent != log->GetEventNumber()) continue;
2725 esd->AddRawDataErrorLog(log);
2730 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2731 // Dump a file content into a char in TNamed
2733 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2734 Int_t kBytes = (Int_t)in.tellg();
2735 printf("Size: %d \n",kBytes);
2738 char* memblock = new char [kBytes];
2739 in.seekg (0, ios::beg);
2740 in.read (memblock, kBytes);
2742 TString fData(memblock,kBytes);
2743 fn = new TNamed(pName,fData);
2744 printf("fData Size: %d \n",fData.Sizeof());
2745 printf("pName Size: %d \n",pName.Sizeof());
2746 printf("fn Size: %d \n",fn->Sizeof());
2750 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2756 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2757 // This is not really needed in AliReconstruction at the moment
2758 // but can serve as a template
2760 TList *fList = fTree->GetUserInfo();
2761 TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2762 printf("fn Size: %d \n",fn->Sizeof());
2764 TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2765 const char* cdata = fn->GetTitle();
2766 printf("fTmp Size %d\n",fTmp.Sizeof());
2768 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2769 printf("calculated size %d\n",size);
2770 ofstream out(pName.Data(),ios::out | ios::binary);
2771 out.write(cdata,size);
2776 //_____________________________________________________________________________
2777 void AliReconstruction::CheckQA()
2779 // check the QA of SIM for this run and remove the detectors
2780 // with status Fatal
2782 TString newRunLocalReconstruction ;
2783 TString newRunTracking ;
2784 TString newFillESD ;
2786 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2787 TString detName(AliQA::GetDetName(iDet)) ;
2788 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2789 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2790 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2792 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2793 fRunLocalReconstruction.Contains("ALL") ) {
2794 newRunLocalReconstruction += detName ;
2795 newRunLocalReconstruction += " " ;
2797 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2798 fRunTracking.Contains("ALL") ) {
2799 newRunTracking += detName ;
2800 newRunTracking += " " ;
2802 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2803 fFillESD.Contains("ALL") ) {
2804 newFillESD += detName ;
2809 fRunLocalReconstruction = newRunLocalReconstruction ;
2810 fRunTracking = newRunTracking ;
2811 fFillESD = newFillESD ;
2814 //_____________________________________________________________________________
2815 Int_t AliReconstruction::GetDetIndex(const char* detector)
2817 // return the detector index corresponding to detector
2819 for (index = 0; index < fgkNDetectors ; index++) {
2820 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2825 //_____________________________________________________________________________
2826 Bool_t AliReconstruction::FinishPlaneEff() {
2828 // Here execute all the necessary operationis, at the end of the tracking phase,
2829 // in case that evaluation of PlaneEfficiencies was required for some detector.
2830 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2832 // This Preliminary version works only FOR ITS !!!!!
2833 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2836 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2839 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2840 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2841 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2842 if(fTracker[iDet]) {
2843 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2844 TString name=planeeff->GetName();
2846 TFile* pefile = TFile::Open(name, "RECREATE");
2847 ret=(Bool_t)planeeff->Write();
2849 if(planeeff->GetCreateHistos()) {
2850 TString hname=planeeff->GetName();
2851 hname+="Histo.root";
2852 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
2858 //_____________________________________________________________________________
2859 Bool_t AliReconstruction::InitPlaneEff() {
2861 // Here execute all the necessary operations, before of the tracking phase,
2862 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2863 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2864 // which should be updated/recalculated.
2866 // This Preliminary version will work only FOR ITS !!!!!
2867 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2870 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2872 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2876 //_____________________________________________________________________________
2877 Bool_t AliReconstruction::InitAliEVE()
2879 // This method should be called only in case
2880 // AliReconstruction is run
2881 // within the alieve environment.
2882 // It will initialize AliEVE in a way
2883 // so that it can visualize event processed
2884 // by AliReconstruction.
2885 // The return flag shows whenever the
2886 // AliEVE initialization was successful or not.
2889 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2890 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2891 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2893 gROOT->ProcessLine("if (!gAliEveEvent) {gAliEveEvent = new AliEveEventManager();gAliEveEvent->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(gAliEveEvent);};");
2894 gROOT->ProcessLine("alieve_online_init()");
2899 //_____________________________________________________________________________
2900 void AliReconstruction::RunAliEVE()
2902 // Runs AliEVE visualisation of
2903 // the current event.
2904 // Should be executed only after
2905 // successful initialization of AliEVE.
2907 AliInfo("Running AliEVE...");
2908 gROOT->ProcessLine(Form("gAliEveEvent->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2912 //_____________________________________________________________________________
2913 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2915 // Allows to run QA for a selected set of detectors
2916 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2917 // all selected detectors run the same selected tasks
2919 if (!detAndAction.Contains(":")) {
2920 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2924 Int_t colon = detAndAction.Index(":") ;
2925 fQADetectors = detAndAction(0, colon) ;
2926 if (fQADetectors.Contains("ALL") )
2927 fQADetectors = fFillESD ;
2928 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2929 if (fQATasks.Contains("ALL") ) {
2930 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2932 fQATasks.ToUpper() ;
2934 if ( fQATasks.Contains("RAW") )
2935 tempo = Form("%d ", AliQA::kRAWS) ;
2936 if ( fQATasks.Contains("RECPOINT") )
2937 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2938 if ( fQATasks.Contains("ESD") )
2939 tempo += Form("%d ", AliQA::kESDS) ;
2941 if (fQATasks.IsNull()) {
2942 AliInfo("No QA requested\n") ;
2947 TString tempo(fQATasks) ;
2948 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2949 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2950 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2951 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2956 //_____________________________________________________________________________
2957 Bool_t AliReconstruction::InitRecoParams()
2959 // The method accesses OCDB and retrieves all
2960 // the available reco-param objects from there.
2962 Bool_t isOK = kTRUE;
2964 TString detStr = fLoadCDB;
2965 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2967 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2969 if (fRecoParam.GetDetRecoParamArray(iDet)) {
2970 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
2974 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
2976 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
2977 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
2979 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
2983 TObject *recoParamObj = entry->GetObject();
2984 if (dynamic_cast<TObjArray*>(recoParamObj)) {
2985 // The detector has a normal TobjArray of AliDetectorRecoParam objects
2986 // Registering them in AliRecoParam
2987 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
2989 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
2990 // The detector has only onse set of reco parameters
2991 // Registering it in AliRecoParam
2992 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
2993 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
2994 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
2997 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3001 AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3005 if (AliDebugLevel() > 0) fRecoParam.Print();
3010 //_____________________________________________________________________________
3011 Bool_t AliReconstruction::GetEventInfo()
3013 // Fill the event info object
3015 AliCodeTimerAuto("")
3017 AliCentralTrigger *aCTP = NULL;
3019 fEventInfo.SetEventType(fRawReader->GetType());
3021 ULong64_t mask = fRawReader->GetClassMask();
3022 fEventInfo.SetTriggerMask(mask);
3023 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3024 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3026 aCTP = new AliCentralTrigger();
3027 TString configstr("");
3028 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3029 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3033 aCTP->SetClassMask(mask);
3034 aCTP->SetClusterMask(clmask);
3037 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3039 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3040 aCTP = fRunLoader->GetTrigger();
3041 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3042 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3045 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3050 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3052 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3053 if (fRawReader) delete aCTP;
3057 UChar_t clustmask = 0;
3059 ULong64_t trmask = fEventInfo.GetTriggerMask();
3060 const TObjArray& classesArray = config->GetClasses();
3061 Int_t nclasses = classesArray.GetEntriesFast();
3062 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3063 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3065 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
3066 fesd->SetTriggerClass(trclass->GetName(),trindex);
3067 if (trmask & (1 << trindex)) {
3069 trclasses += trclass->GetName();
3071 clustmask |= trclass->GetCluster()->GetClusterMask();
3075 fEventInfo.SetTriggerClasses(trclasses);
3077 // Set the information in ESD
3078 fesd->SetTriggerMask(trmask);
3079 fesd->SetTriggerCluster(clustmask);
3081 if (!aCTP->CheckTriggeredDetectors()) {
3082 if (fRawReader) delete aCTP;
3086 if (fRawReader) delete aCTP;
3088 // We have to fill also the HLT decision here!!
3094 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3096 // Match the detector list found in the rec.C or the default 'ALL'
3097 // to the list found in the GRP (stored there by the shuttle PP which
3098 // gets the information from ECS)
3099 static TString resultList;
3100 TString detList = detectorList;
3104 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3105 if ((detectorMask >> iDet) & 0x1) {
3106 TString det = AliDAQ::OfflineModuleName(iDet);
3107 if ((detList.CompareTo("ALL") == 0) ||
3108 detList.BeginsWith("ALL ") ||
3109 detList.EndsWith(" ALL") ||
3110 detList.Contains(" ALL ") ||
3111 (detList.CompareTo(det) == 0) ||
3112 detList.BeginsWith(det) ||
3113 detList.EndsWith(det) ||
3114 detList.Contains( " "+det+" " )) {
3115 if (!resultList.EndsWith(det + " ")) {
3124 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3125 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3126 if ((detList.CompareTo("ALL") == 0) ||
3127 detList.BeginsWith("ALL ") ||
3128 detList.EndsWith(" ALL") ||
3129 detList.Contains(" ALL ") ||
3130 (detList.CompareTo(hltDet) == 0) ||
3131 detList.BeginsWith(hltDet) ||
3132 detList.EndsWith(hltDet) ||
3133 detList.Contains( " "+hltDet+" " )) {
3134 resultList += hltDet;
3138 return resultList.Data();
3142 //______________________________________________________________________________
3143 void AliReconstruction::Abort(const char *method, EAbort what)
3145 // Abort processing. If what = kAbortProcess, the Process() loop will be
3146 // aborted. If what = kAbortFile, the current file in a chain will be
3147 // aborted and the processing will continue with the next file, if there
3148 // is no next file then Process() will be aborted. Abort() can also be
3149 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3150 // the SlaveTerminate() and Terminate() are always called. The abort flag
3151 // can be checked in these methods using GetAbort().
3153 // The method is overwritten in AliReconstruction for better handling of
3154 // reco specific errors
3156 if (!fStopOnError) return;
3160 TString whyMess = method;
3161 whyMess += " failed! Aborting...";
3163 AliError(whyMess.Data());
3166 TString mess = "Abort";
3167 if (fAbort == kAbortProcess)
3168 mess = "AbortProcess";
3169 else if (fAbort == kAbortFile)
3172 Info(mess, whyMess.Data());