1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for running the reconstruction //
22 // Clusters and tracks are created for all detectors and all events by //
25 // AliReconstruction rec; //
28 // The Run method returns kTRUE in case of successful execution. //
30 // If the input to the reconstruction are not simulated digits but raw data, //
31 // this can be specified by an argument of the Run method or by the method //
33 // rec.SetInput("..."); //
35 // The input formats and the corresponding argument are: //
36 // - DDL raw data files: directory name, ends with "/" //
37 // - raw data root file: root file name, extension ".root" //
38 // - raw data DATE file: DATE file name, any other non-empty string //
39 // - MC root files : empty string, default //
41 // By default all events are reconstructed. The reconstruction can be //
42 // limited to a range of events by giving the index of the first and the //
43 // last event as an argument to the Run method or by calling //
45 // rec.SetEventRange(..., ...); //
47 // The index -1 (default) can be used for the last event to indicate no //
48 // upper limit of the event range. //
50 // In case of raw-data reconstruction the user can modify the default //
51 // number of events per digits/clusters/tracks file. In case the option //
52 // is not used the number is set 1. In case the user provides 0, than //
53 // the number of events is equal to the number of events inside the //
54 // raw-data file (i.e. one digits/clusters/tracks file): //
56 // rec.SetNumberOfEventsPerFile(...); //
59 // The name of the galice file can be changed from the default //
60 // "galice.root" by passing it as argument to the AliReconstruction //
61 // constructor or by //
63 // rec.SetGAliceFile("..."); //
65 // The local reconstruction can be switched on or off for individual //
68 // rec.SetRunLocalReconstruction("..."); //
70 // The argument is a (case sensitive) string with the names of the //
71 // detectors separated by a space. The special string "ALL" selects all //
72 // available detectors. This is the default. //
74 // The reconstruction of the primary vertex position can be switched off by //
76 // rec.SetRunVertexFinder(kFALSE); //
78 // The tracking and the creation of ESD tracks can be switched on for //
79 // selected detectors by //
81 // rec.SetRunTracking("..."); //
83 // Uniform/nonuniform field tracking switches (default: uniform field) //
85 // rec.SetUniformFieldTracking(); ( rec.SetUniformFieldTracking(kFALSE); ) //
87 // The filling of additional ESD information can be steered by //
89 // rec.SetFillESD("..."); //
91 // Again, for both methods the string specifies the list of detectors. //
92 // The default is "ALL". //
94 // The call of the shortcut method //
96 // rec.SetRunReconstruction("..."); //
98 // is equivalent to calling SetRunLocalReconstruction, SetRunTracking and //
99 // SetFillESD with the same detector selecting string as argument. //
101 // The reconstruction requires digits or raw data as input. For the creation //
102 // of digits and raw data have a look at the class AliSimulation. //
104 // The input data of a detector can be replaced by the corresponding HLT //
105 // data by calling (usual detector string) //
106 // SetUseHLTData("..."); //
109 ///////////////////////////////////////////////////////////////////////////////
116 #include <TGeoGlobalMagField.h>
117 #include <TGeoManager.h>
119 #include <TLorentzVector.h>
121 #include <TObjArray.h>
122 #include <TPRegexp.h>
123 #include <TParameter.h>
124 #include <TPluginManager.h>
126 #include <TProofOutputFile.h>
129 #include <THashTable.h>
131 #include "AliAlignObj.h"
132 #include "AliCDBEntry.h"
133 #include "AliCDBManager.h"
134 #include "AliCDBStorage.h"
135 #include "AliCTPRawStream.h"
136 #include "AliCascadeVertexer.h"
137 #include "AliCentralTrigger.h"
138 #include "AliCodeTimer.h"
140 #include "AliDetectorRecoParam.h"
141 #include "AliESDCaloCells.h"
142 #include "AliESDCaloCluster.h"
143 #include "AliESDEvent.h"
144 #include "AliESDMuonTrack.h"
145 #include "AliESDPmdTrack.h"
146 #include "AliESDTagCreator.h"
147 #include "AliESDVertex.h"
148 #include "AliESDcascade.h"
149 #include "AliESDfriend.h"
150 #include "AliESDkink.h"
151 #include "AliESDpid.h"
152 #include "AliESDtrack.h"
153 #include "AliESDtrack.h"
154 #include "AliEventInfo.h"
155 #include "AliGRPObject.h"
156 #include "AliGRPRecoParam.h"
157 #include "AliGenEventHeader.h"
158 #include "AliGeomManager.h"
159 #include "AliGlobalQADataMaker.h"
160 #include "AliHeader.h"
163 #include "AliMultiplicity.h"
165 #include "AliPlaneEff.h"
167 #include "AliQADataMakerRec.h"
168 #include "AliQAManager.h"
169 #include "AliRawVEvent.h"
170 #include "AliRawEventHeaderBase.h"
171 #include "AliRawHLTManager.h"
172 #include "AliRawReaderDate.h"
173 #include "AliRawReaderFile.h"
174 #include "AliRawReaderRoot.h"
175 #include "AliReconstruction.h"
176 #include "AliReconstructor.h"
178 #include "AliRunInfo.h"
179 #include "AliRunLoader.h"
180 #include "AliSysInfo.h" // memory snapshots
181 #include "AliTrackPointArray.h"
182 #include "AliTracker.h"
183 #include "AliTriggerClass.h"
184 #include "AliTriggerCluster.h"
185 #include "AliTriggerIR.h"
186 #include "AliTriggerConfiguration.h"
187 #include "AliV0vertexer.h"
188 #include "AliVertexer.h"
189 #include "AliVertexerTracks.h"
191 ClassImp(AliReconstruction)
193 //_____________________________________________________________________________
194 const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
196 //_____________________________________________________________________________
197 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
199 fRunVertexFinder(kTRUE),
200 fRunVertexFinderTracks(kTRUE),
201 fRunHLTTracking(kFALSE),
202 fRunMuonTracking(kFALSE),
204 fRunCascadeFinder(kTRUE),
205 fStopOnError(kFALSE),
206 fWriteAlignmentData(kFALSE),
207 fWriteESDfriend(kFALSE),
208 fFillTriggerESD(kTRUE),
216 fRunLocalReconstruction("ALL"),
220 fUseTrackingErrorsForAlignment(""),
221 fGAliceFileName(gAliceFilename),
227 fNumberOfEventsPerFile((UInt_t)-1),
229 fLoadAlignFromCDB(kTRUE),
230 fLoadAlignData("ALL"),
238 fParentRawReader(NULL),
242 fSPDTrackleter(NULL),
244 fDiamondProfileSPD(NULL),
245 fDiamondProfile(NULL),
246 fDiamondProfileTPC(NULL),
247 fListOfCosmicTriggers(NULL),
251 fAlignObjArray(NULL),
255 fInitCDBCalled(kFALSE),
256 fSetRunNumberFromDataCalled(kFALSE),
261 fSameQACycle(kFALSE),
262 fInitQACalled(kFALSE),
263 fWriteQAExpertData(kTRUE),
264 fRunPlaneEff(kFALSE),
273 fIsNewRunLoader(kFALSE),
277 // create reconstruction object with default parameters
280 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
281 fReconstructor[iDet] = NULL;
282 fLoader[iDet] = NULL;
283 fTracker[iDet] = NULL;
285 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
286 fQACycles[iDet] = 999999 ;
287 fQAWriteExpert[iDet] = kFALSE ;
293 //_____________________________________________________________________________
294 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
296 fRunVertexFinder(rec.fRunVertexFinder),
297 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
298 fRunHLTTracking(rec.fRunHLTTracking),
299 fRunMuonTracking(rec.fRunMuonTracking),
300 fRunV0Finder(rec.fRunV0Finder),
301 fRunCascadeFinder(rec.fRunCascadeFinder),
302 fStopOnError(rec.fStopOnError),
303 fWriteAlignmentData(rec.fWriteAlignmentData),
304 fWriteESDfriend(rec.fWriteESDfriend),
305 fFillTriggerESD(rec.fFillTriggerESD),
307 fCleanESD(rec.fCleanESD),
308 fV0DCAmax(rec.fV0DCAmax),
309 fV0CsPmin(rec.fV0CsPmin),
313 fRunLocalReconstruction(rec.fRunLocalReconstruction),
314 fRunTracking(rec.fRunTracking),
315 fFillESD(rec.fFillESD),
316 fLoadCDB(rec.fLoadCDB),
317 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
318 fGAliceFileName(rec.fGAliceFileName),
319 fRawInput(rec.fRawInput),
320 fESDOutput(rec.fESDOutput),
321 fEquipIdMap(rec.fEquipIdMap),
322 fFirstEvent(rec.fFirstEvent),
323 fLastEvent(rec.fLastEvent),
324 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
326 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
327 fLoadAlignData(rec.fLoadAlignData),
328 fUseHLTData(rec.fUseHLTData),
335 fParentRawReader(NULL),
337 fRecoParam(rec.fRecoParam),
339 fSPDTrackleter(NULL),
341 fDiamondProfileSPD(rec.fDiamondProfileSPD),
342 fDiamondProfile(rec.fDiamondProfile),
343 fDiamondProfileTPC(rec.fDiamondProfileTPC),
344 fListOfCosmicTriggers(NULL),
348 fAlignObjArray(rec.fAlignObjArray),
349 fCDBUri(rec.fCDBUri),
350 fQARefUri(rec.fQARefUri),
352 fInitCDBCalled(rec.fInitCDBCalled),
353 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
354 fQADetectors(rec.fQADetectors),
355 fQATasks(rec.fQATasks),
357 fRunGlobalQA(rec.fRunGlobalQA),
358 fSameQACycle(rec.fSameQACycle),
359 fInitQACalled(rec.fInitQACalled),
360 fWriteQAExpertData(rec.fWriteQAExpertData),
361 fRunPlaneEff(rec.fRunPlaneEff),
370 fIsNewRunLoader(rec.fIsNewRunLoader),
376 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
377 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
379 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
380 fReconstructor[iDet] = NULL;
381 fLoader[iDet] = NULL;
382 fTracker[iDet] = NULL;
385 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
386 fQACycles[iDet] = rec.fQACycles[iDet];
387 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
390 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
391 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
395 //_____________________________________________________________________________
396 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
398 // assignment operator
399 // Used in PROOF mode
400 // Be very careful while modifing it!
401 // Simple rules to follow:
402 // for persistent data members - use their assignment operators
403 // for non-persistent ones - do nothing or take the default values from constructor
404 // TSelector members should not be touched
405 if(&rec == this) return *this;
407 fRunVertexFinder = rec.fRunVertexFinder;
408 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
409 fRunHLTTracking = rec.fRunHLTTracking;
410 fRunMuonTracking = rec.fRunMuonTracking;
411 fRunV0Finder = rec.fRunV0Finder;
412 fRunCascadeFinder = rec.fRunCascadeFinder;
413 fStopOnError = rec.fStopOnError;
414 fWriteAlignmentData = rec.fWriteAlignmentData;
415 fWriteESDfriend = rec.fWriteESDfriend;
416 fFillTriggerESD = rec.fFillTriggerESD;
418 fCleanESD = rec.fCleanESD;
419 fV0DCAmax = rec.fV0DCAmax;
420 fV0CsPmin = rec.fV0CsPmin;
424 fRunLocalReconstruction = rec.fRunLocalReconstruction;
425 fRunTracking = rec.fRunTracking;
426 fFillESD = rec.fFillESD;
427 fLoadCDB = rec.fLoadCDB;
428 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
429 fGAliceFileName = rec.fGAliceFileName;
430 fRawInput = rec.fRawInput;
431 fESDOutput = rec.fESDOutput;
432 fEquipIdMap = rec.fEquipIdMap;
433 fFirstEvent = rec.fFirstEvent;
434 fLastEvent = rec.fLastEvent;
435 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
437 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
438 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
441 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
442 fLoadAlignData = rec.fLoadAlignData;
443 fUseHLTData = rec.fUseHLTData;
445 delete fRunInfo; fRunInfo = NULL;
446 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
448 fEventInfo = rec.fEventInfo;
450 delete fRunScalers; fRunScalers = NULL;
451 if (rec.fRunScalers) fRunScalers = new AliTriggerRunScalers(*rec.fRunScalers);
456 fParentRawReader = NULL;
458 fRecoParam = rec.fRecoParam;
460 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
461 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
462 delete fLoader[iDet]; fLoader[iDet] = NULL;
463 delete fTracker[iDet]; fTracker[iDet] = NULL;
466 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
467 fQACycles[iDet] = rec.fQACycles[iDet];
468 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
471 delete fSPDTrackleter; fSPDTrackleter = NULL;
473 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
474 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
475 delete fDiamondProfile; fDiamondProfile = NULL;
476 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
477 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
478 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
480 delete fListOfCosmicTriggers; fListOfCosmicTriggers = NULL;
481 if (rec.fListOfCosmicTriggers) fListOfCosmicTriggers = (THashTable*)((rec.fListOfCosmicTriggers)->Clone());
483 delete fGRPData; fGRPData = NULL;
484 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
485 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
487 delete fAlignObjArray; fAlignObjArray = NULL;
490 fQARefUri = rec.fQARefUri;
491 fSpecCDBUri.Delete();
492 fInitCDBCalled = rec.fInitCDBCalled;
493 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
494 fQADetectors = rec.fQADetectors;
495 fQATasks = rec.fQATasks;
497 fRunGlobalQA = rec.fRunGlobalQA;
498 fSameQACycle = rec.fSameQACycle;
499 fInitQACalled = rec.fInitQACalled;
500 fWriteQAExpertData = rec.fWriteQAExpertData;
501 fRunPlaneEff = rec.fRunPlaneEff;
510 fIsNewRunLoader = rec.fIsNewRunLoader;
517 //_____________________________________________________________________________
518 AliReconstruction::~AliReconstruction()
523 if (fListOfCosmicTriggers) {
524 fListOfCosmicTriggers->Delete();
525 delete fListOfCosmicTriggers;
530 if (fAlignObjArray) {
531 fAlignObjArray->Delete();
532 delete fAlignObjArray;
534 fSpecCDBUri.Delete();
536 AliCodeTimer::Instance()->Print();
539 //_____________________________________________________________________________
540 void AliReconstruction::InitQA()
542 //Initialize the QA and start of cycle
543 AliCodeTimerAuto("");
545 if (fInitQACalled) return;
546 fInitQACalled = kTRUE;
548 AliQAManager * qam = AliQAManager::QAManager("rec") ;
549 if (fWriteQAExpertData)
550 qam->SetWriteExpert() ;
552 if (qam->IsDefaultStorageSet()) {
553 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
554 AliWarning("Default QA reference storage has been already set !");
555 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fQARefUri.Data()));
556 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
557 fQARefUri = qam->GetDefaultStorage()->GetURI();
559 if (fQARefUri.Length() > 0) {
560 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
561 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
562 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
564 fQARefUri="local://$ALICE_ROOT/QAref";
565 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
566 AliWarning("Default QA refeference storage not yet set !!!!");
567 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
568 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
571 qam->SetDefaultStorage(fQARefUri);
575 qam->SetActiveDetectors(fQADetectors) ;
576 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) {
577 qam->SetCycleLength(AliQAv1::DETECTORINDEX_t(det), fQACycles[det]) ;
578 qam->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
580 if (!fRawReader && !fInput && fQATasks.Contains(AliQAv1::kRAWS))
581 fQATasks.ReplaceAll(Form("%d",AliQAv1::kRAWS), "") ;
582 qam->SetTasks(fQATasks) ;
583 qam->InitQADataMaker(AliCDBManager::Instance()->GetRun()) ;
586 Bool_t sameCycle = kFALSE ;
587 AliQADataMaker *qadm = qam->GetQADataMaker(AliQAv1::kGLOBAL);
588 AliInfo(Form("Initializing the global QA data maker"));
589 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS))) {
590 qadm->StartOfCycle(AliQAv1::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
591 TObjArray **arr=qadm->Init(AliQAv1::kRECPOINTS);
592 AliTracker::SetResidualsArray(arr);
595 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS))) {
596 qadm->StartOfCycle(AliQAv1::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
597 qadm->Init(AliQAv1::kESDS);
600 AliSysInfo::AddStamp("InitQA") ;
603 //_____________________________________________________________________________
604 void AliReconstruction::MergeQA(const char *fileName)
606 //Initialize the QA and start of cycle
607 AliCodeTimerAuto("") ;
608 AliQAManager::QAManager()->Merge(AliCDBManager::Instance()->GetRun(),fileName) ;
609 AliSysInfo::AddStamp("MergeQA") ;
612 //_____________________________________________________________________________
613 void AliReconstruction::InitCDB()
615 // activate a default CDB storage
616 // First check if we have any CDB storage set, because it is used
617 // to retrieve the calibration and alignment constants
618 AliCodeTimerAuto("");
620 if (fInitCDBCalled) return;
621 fInitCDBCalled = kTRUE;
623 AliCDBManager* man = AliCDBManager::Instance();
624 if (man->IsDefaultStorageSet())
626 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
627 AliWarning("Default CDB storage has been already set !");
628 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
629 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
630 fCDBUri = man->GetDefaultStorage()->GetURI();
633 if (fCDBUri.Length() > 0)
635 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
636 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
637 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
639 fCDBUri="local://$ALICE_ROOT/OCDB";
640 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
641 AliWarning("Default CDB storage not yet set !!!!");
642 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
643 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
646 man->SetDefaultStorage(fCDBUri);
649 // Now activate the detector specific CDB storage locations
650 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
651 TObject* obj = fSpecCDBUri[i];
653 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
654 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
655 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
656 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
658 AliSysInfo::AddStamp("InitCDB");
661 //_____________________________________________________________________________
662 void AliReconstruction::SetDefaultStorage(const char* uri) {
663 // Store the desired default CDB storage location
664 // Activate it later within the Run() method
670 //_____________________________________________________________________________
671 void AliReconstruction::SetQARefDefaultStorage(const char* uri) {
672 // Store the desired default CDB storage location
673 // Activate it later within the Run() method
676 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
679 //_____________________________________________________________________________
680 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
681 // Store a detector-specific CDB storage location
682 // Activate it later within the Run() method
684 AliCDBPath aPath(calibType);
685 if(!aPath.IsValid()){
686 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
687 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
688 if(!strcmp(calibType, fgkDetectorName[iDet])) {
689 aPath.SetPath(Form("%s/*", calibType));
690 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
694 if(!aPath.IsValid()){
695 AliError(Form("Not a valid path or detector: %s", calibType));
700 // // check that calibType refers to a "valid" detector name
701 // Bool_t isDetector = kFALSE;
702 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
703 // TString detName = fgkDetectorName[iDet];
704 // if(aPath.GetLevel0() == detName) {
705 // isDetector = kTRUE;
711 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
715 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
716 if (obj) fSpecCDBUri.Remove(obj);
717 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
721 //_____________________________________________________________________________
722 Bool_t AliReconstruction::SetRunNumberFromData()
724 // The method is called in Run() in order
725 // to set a correct run number.
726 // In case of raw data reconstruction the
727 // run number is taken from the raw data header
729 if (fSetRunNumberFromDataCalled) return kTRUE;
730 fSetRunNumberFromDataCalled = kTRUE;
732 AliCDBManager* man = AliCDBManager::Instance();
735 if(fRawReader->NextEvent()) {
736 if(man->GetRun() > 0) {
737 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
739 man->SetRun(fRawReader->GetRunNumber());
740 fRawReader->RewindEvents();
743 if(man->GetRun() > 0) {
744 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
747 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
753 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
755 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
760 // read run number from gAlice
761 if(rl->GetHeader()) {
762 man->SetRun(rl->GetHeader()->GetRun());
767 AliError("Neither run-loader header nor RawReader objects are found !");
779 //_____________________________________________________________________________
780 void AliReconstruction::SetCDBLock() {
781 // Set CDB lock: from now on it is forbidden to reset the run number
782 // or the default storage or to activate any further storage!
784 AliCDBManager::Instance()->SetLock(1);
787 //_____________________________________________________________________________
788 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
790 // Read the alignment objects from CDB.
791 // Each detector is supposed to have the
792 // alignment objects in DET/Align/Data CDB path.
793 // All the detector objects are then collected,
794 // sorted by geometry level (starting from ALIC) and
795 // then applied to the TGeo geometry.
796 // Finally an overlaps check is performed.
798 // Load alignment data from CDB and fill fAlignObjArray
799 if(fLoadAlignFromCDB){
801 TString detStr = detectors;
802 TString loadAlObjsListOfDets = "";
804 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
805 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
806 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
808 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
810 loadAlObjsListOfDets += fgkDetectorName[iDet];
811 loadAlObjsListOfDets += " ";
813 } // end loop over detectors
815 if(AliGeomManager::GetNalignable("GRP") != 0)
816 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
817 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
818 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
820 // Check if the array with alignment objects was
821 // provided by the user. If yes, apply the objects
822 // to the present TGeo geometry
823 if (fAlignObjArray) {
824 if (gGeoManager && gGeoManager->IsClosed()) {
825 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
826 AliError("The misalignment of one or more volumes failed!"
827 "Compare the list of simulated detectors and the list of detector alignment data!");
832 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
838 if (fAlignObjArray) {
839 fAlignObjArray->Delete();
840 delete fAlignObjArray; fAlignObjArray=NULL;
846 //_____________________________________________________________________________
847 void AliReconstruction::SetGAliceFile(const char* fileName)
849 // set the name of the galice file
851 fGAliceFileName = fileName;
854 //_____________________________________________________________________________
855 void AliReconstruction::SetInput(const char* input)
857 // In case the input string starts with 'mem://', we run in an online mode
858 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
859 // file is assumed. One can give as an input:
860 // mem://: - events taken from DAQ monitoring libs online
862 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
863 if (input) fRawInput = input;
866 //_____________________________________________________________________________
867 void AliReconstruction::SetOutput(const char* output)
869 // Set the output ESD filename
870 // 'output' is a normalt ROOT url
871 // The method is used in case of raw-data reco with PROOF
872 if (output) fESDOutput.SetUrl(output);
875 //_____________________________________________________________________________
876 void AliReconstruction::SetOption(const char* detector, const char* option)
878 // set options for the reconstruction of a detector
880 TObject* obj = fOptions.FindObject(detector);
881 if (obj) fOptions.Remove(obj);
882 fOptions.Add(new TNamed(detector, option));
885 //_____________________________________________________________________________
886 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
888 // Set custom reconstruction parameters for a given detector
889 // Single set of parameters for all the events
891 // First check if the reco-params are global
892 if(!strcmp(detector, "GRP")) {
894 fRecoParam.AddDetRecoParam(kNDetectors,par);
898 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
899 if(!strcmp(detector, fgkDetectorName[iDet])) {
901 fRecoParam.AddDetRecoParam(iDet,par);
908 //_____________________________________________________________________________
909 Bool_t AliReconstruction::SetFieldMap(Float_t l3Cur, Float_t diCur, Float_t l3Pol,
910 Float_t diPol, Float_t beamenergy,
911 const Char_t *beamtype, const Char_t *path)
913 //------------------------------------------------
914 // The magnetic field map, defined externally...
915 // L3 current 30000 A -> 0.5 T
916 // L3 current 12000 A -> 0.2 T
917 // dipole current 6000 A
918 // The polarities must be the same
919 //------------------------------------------------
920 const Float_t l3NominalCurrent1=30000.; // (A)
921 const Float_t l3NominalCurrent2=12000.; // (A)
922 const Float_t diNominalCurrent =6000. ; // (A)
924 const Float_t tolerance=0.03; // relative current tolerance
925 const Float_t zero=77.; // "zero" current (A)
927 TString s=(l3Pol < 0) ? "L3: -" : "L3: +";
929 AliMagF::BMap_t map = AliMagF::k5kG;
933 l3Cur = TMath::Abs(l3Cur);
934 if (TMath::Abs(l3Cur-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
935 fcL3 = l3Cur/l3NominalCurrent1;
938 } else if (TMath::Abs(l3Cur-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
939 fcL3 = l3Cur/l3NominalCurrent2;
942 } else if (l3Cur <= zero) {
944 map = AliMagF::k5kGUniform;
947 AliError(Form("Wrong L3 current (%f A)!",l3Cur));
951 diCur = TMath::Abs(diCur);
952 if (TMath::Abs(diCur-diNominalCurrent)/diNominalCurrent < tolerance) {
953 // 3% current tolerance...
954 fcDip = diCur/diNominalCurrent;
956 } else if (diCur <= zero) { // some small current..
960 AliError(Form("Wrong dipole current (%f A)!",diCur));
964 if (l3Pol!=diPol && (map==AliMagF::k5kG || map==AliMagF::k2kG) && fcDip!=0) {
965 AliError("L3 and Dipole polarities must be the same");
969 if (l3Pol<0) fcL3 = -fcL3;
970 if (diPol<0) fcDip = -fcDip;
972 AliMagF::BeamType_t btype = AliMagF::kNoBeamField;
973 TString btypestr = beamtype;
975 TPRegexp protonBeam("(proton|p)\\s*-?\\s*\\1");
976 TPRegexp ionBeam("(lead|pb|ion|a)\\s*-?\\s*\\1");
977 if (btypestr.Contains(ionBeam)) btype = AliMagF::kBeamTypeAA;
978 else if (btypestr.Contains(protonBeam)) btype = AliMagF::kBeamTypepp;
980 AliInfo(Form("Cannot determine the beam type from %s, assume no LHC magnet field",beamtype));
983 AliMagF* fld = new AliMagF("MagneticFieldMap", s.Data(), 2, fcL3, fcDip, 10., map, path,
985 TGeoGlobalMagField::Instance()->SetField( fld );
986 TGeoGlobalMagField::Instance()->Lock();
992 Bool_t AliReconstruction::InitGRP() {
993 //------------------------------------
994 // Initialization of the GRP entry
995 //------------------------------------
996 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
1000 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
1003 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
1005 fGRPData = new AliGRPObject();
1006 fGRPData->ReadValuesFromMap(m);
1010 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
1011 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
1015 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
1019 AliError("No GRP entry found in OCDB!");
1023 TString lhcState = fGRPData->GetLHCState();
1024 if (lhcState==AliGRPObject::GetInvalidString()) {
1025 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
1026 lhcState = "UNKNOWN";
1029 TString beamType = fGRPData->GetBeamType();
1030 if (beamType==AliGRPObject::GetInvalidString()) {
1031 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
1032 beamType = "UNKNOWN";
1035 Float_t beamEnergy = fGRPData->GetBeamEnergy();
1036 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
1037 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
1040 // energy is provided in MeV*120
1041 beamEnergy /= 120E3;
1043 TString runType = fGRPData->GetRunType();
1044 if (runType==AliGRPObject::GetInvalidString()) {
1045 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
1046 runType = "UNKNOWN";
1049 Int_t activeDetectors = fGRPData->GetDetectorMask();
1050 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
1051 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
1052 activeDetectors = 1074790399;
1055 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1059 // Process the list of active detectors
1060 if (activeDetectors) {
1061 UInt_t detMask = activeDetectors;
1062 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1063 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1064 fFillESD = MatchDetectorList(fFillESD,detMask);
1065 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1066 fLoadCDB.Form("%s %s %s %s",
1067 fRunLocalReconstruction.Data(),
1068 fRunTracking.Data(),
1070 fQADetectors.Data());
1071 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1072 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1)) {
1073 // switch off the vertexer
1074 AliInfo("SPD is not in the list of active detectors. Vertexer switched off.");
1075 fRunVertexFinder = kFALSE;
1077 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1078 // switch off the reading of CTP raw-data payload
1079 if (fFillTriggerESD) {
1080 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1081 fFillTriggerESD = kFALSE;
1086 AliInfo("===================================================================================");
1087 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1088 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1089 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1090 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1091 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1092 AliInfo("===================================================================================");
1094 //*** Dealing with the magnetic field map
1095 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {AliInfo("Running with the externally locked B field !");}
1097 // Construct the field map out of the information retrieved from GRP.
1100 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1101 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1102 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1106 Char_t l3Polarity = fGRPData->GetL3Polarity();
1107 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1108 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1113 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1114 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1115 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1119 Char_t diPolarity = fGRPData->GetDipolePolarity();
1120 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1121 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1126 TObjString *l3Current=
1127 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
1129 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1132 TObjString *l3Polarity=
1133 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
1135 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1140 TObjString *diCurrent=
1141 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
1143 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1146 TObjString *diPolarity=
1147 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
1149 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1155 if ( !SetFieldMap(l3Current, diCurrent, l3Polarity ? -1:1, diPolarity ? -1:1) )
1156 AliFatal("Failed to creat a B field map ! Exiting...");
1157 AliInfo("Running with the B field constructed out of GRP !");
1159 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1163 //*** Get the diamond profiles from OCDB
1164 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1166 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1168 AliError("No SPD diamond profile found in OCDB!");
1171 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1173 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1175 AliError("No diamond profile found in OCDB!");
1178 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1180 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1182 AliError("No TPC diamond profile found in OCDB!");
1185 entry = AliCDBManager::Instance()->Get("GRP/Calib/CosmicTriggers");
1187 fListOfCosmicTriggers = dynamic_cast<THashTable*>(entry->GetObject());
1189 AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
1192 if (!fListOfCosmicTriggers) {
1193 AliWarning("Can not get list of cosmic triggers from OCDB! Cosmic event specie will be effectively disabled!");
1199 //_____________________________________________________________________________
1200 Bool_t AliReconstruction::LoadCDB()
1202 AliCodeTimerAuto("");
1204 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1206 TString detStr = fLoadCDB;
1207 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1208 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1209 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1213 //_____________________________________________________________________________
1214 Bool_t AliReconstruction::LoadTriggerScalersCDB()
1216 AliCodeTimerAuto("");
1218 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
1222 AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
1223 fRunScalers = dynamic_cast<AliTriggerRunScalers*> (entry->GetObject());
1228 //_____________________________________________________________________________
1229 Bool_t AliReconstruction::Run(const char* input)
1232 AliCodeTimerAuto("");
1235 if (GetAbort() != TSelector::kContinue) return kFALSE;
1237 TChain *chain = NULL;
1238 if (fRawReader && (chain = fRawReader->GetChain())) {
1241 gProof->AddInput(this);
1242 if (!fESDOutput.IsValid()) {
1243 fESDOutput.SetProtocol("root",kTRUE);
1244 fESDOutput.SetHost(gSystem->HostName());
1245 fESDOutput.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1247 AliInfo(Form("Output file with ESDs is %s",fESDOutput.GetUrl()));
1248 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",fESDOutput.GetUrl()));
1249 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1251 chain->Process("AliReconstruction");
1254 chain->Process(this);
1259 if (GetAbort() != TSelector::kContinue) return kFALSE;
1261 if (GetAbort() != TSelector::kContinue) return kFALSE;
1262 //******* The loop over events
1263 AliInfo("Starting looping over events");
1265 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1266 (fRawReader && fRawReader->NextEvent())) {
1267 if (!ProcessEvent(iEvent)) {
1268 Abort("ProcessEvent",TSelector::kAbortFile);
1274 if (GetAbort() != TSelector::kContinue) return kFALSE;
1276 if (GetAbort() != TSelector::kContinue) return kFALSE;
1282 //_____________________________________________________________________________
1283 void AliReconstruction::InitRawReader(const char* input)
1285 AliCodeTimerAuto("");
1287 // Init raw-reader and
1288 // set the input in case of raw data
1289 if (input) fRawInput = input;
1290 fRawReader = AliRawReader::Create(fRawInput.Data());
1292 AliInfo("Reconstruction will run over digits");
1294 if (!fEquipIdMap.IsNull() && fRawReader)
1295 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1297 if (!fUseHLTData.IsNull()) {
1298 // create the RawReaderHLT which performs redirection of HLT input data for
1299 // the specified detectors
1300 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1302 fParentRawReader=fRawReader;
1303 fRawReader=pRawReader;
1305 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1308 AliSysInfo::AddStamp("CreateRawReader");
1311 //_____________________________________________________________________________
1312 void AliReconstruction::InitRun(const char* input)
1314 // Initialization of raw-reader,
1315 // run number, CDB etc.
1316 AliCodeTimerAuto("");
1317 AliSysInfo::AddStamp("Start");
1319 // Initialize raw-reader if any
1320 InitRawReader(input);
1322 // Initialize the CDB storage
1325 // Set run number in CDBManager (if it is not already set by the user)
1326 if (!SetRunNumberFromData()) {
1327 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1331 // Set CDB lock: from now on it is forbidden to reset the run number
1332 // or the default storage or to activate any further storage!
1337 //_____________________________________________________________________________
1338 void AliReconstruction::Begin(TTree *)
1340 // Initialize AlReconstruction before
1341 // going into the event loop
1342 // Should follow the TSelector convention
1343 // i.e. initialize only the object on the client side
1344 AliCodeTimerAuto("");
1346 AliReconstruction *reco = NULL;
1348 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1351 AliSysInfo::AddStamp("ReadInputInBegin");
1354 // Import ideal TGeo geometry and apply misalignment
1356 TString geom(gSystem->DirName(fGAliceFileName));
1357 geom += "/geometry.root";
1358 AliGeomManager::LoadGeometry(geom.Data());
1360 Abort("LoadGeometry", TSelector::kAbortProcess);
1363 AliSysInfo::AddStamp("LoadGeom");
1364 TString detsToCheck=fRunLocalReconstruction;
1365 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1366 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1369 AliSysInfo::AddStamp("CheckGeom");
1372 if (!MisalignGeometry(fLoadAlignData)) {
1373 Abort("MisalignGeometry", TSelector::kAbortProcess);
1376 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1377 AliSysInfo::AddStamp("MisalignGeom");
1380 Abort("InitGRP", TSelector::kAbortProcess);
1383 AliSysInfo::AddStamp("InitGRP");
1386 Abort("LoadCDB", TSelector::kAbortProcess);
1389 AliSysInfo::AddStamp("LoadCDB");
1391 if (!LoadTriggerScalersCDB()) {
1392 Abort("LoadTriggerScalersCDB", TSelector::kAbortProcess);
1395 AliSysInfo::AddStamp("LoadTriggerScalersCDB");
1398 // Read the reconstruction parameters from OCDB
1399 if (!InitRecoParams()) {
1400 AliWarning("Not all detectors have correct RecoParam objects initialized");
1402 AliSysInfo::AddStamp("InitRecoParams");
1404 if (fInput && gProof) {
1405 if (reco) *reco = *this;
1407 gProof->AddInputData(gGeoManager,kTRUE);
1409 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1410 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1411 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1412 magFieldMap->SetName("MagneticFieldMap");
1413 gProof->AddInputData(magFieldMap,kTRUE);
1418 //_____________________________________________________________________________
1419 void AliReconstruction::SlaveBegin(TTree*)
1421 // Initialization related to run-loader,
1422 // vertexer, trackers, recontructors
1423 // In proof mode it is executed on the slave
1424 AliCodeTimerAuto("");
1426 TProofOutputFile *outProofFile = NULL;
1428 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1431 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1433 AliGeomManager::SetGeometry(tgeo);
1435 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1436 Int_t runNumber = -1;
1437 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1438 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1439 man->SetCacheFlag(kTRUE);
1440 man->SetLock(kTRUE);
1444 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1445 TGeoGlobalMagField::Instance()->SetField(map);
1447 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1448 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1449 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1450 fOutput->Add(outProofFile);
1452 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1455 // get the run loader
1456 if (!InitRunLoader()) {
1457 Abort("InitRunLoader", TSelector::kAbortProcess);
1460 AliSysInfo::AddStamp("LoadLoader");
1462 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1465 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1466 Abort("CreateTrackers", TSelector::kAbortProcess);
1469 AliSysInfo::AddStamp("CreateTrackers");
1471 // create the ESD output file and tree
1472 if (!outProofFile) {
1473 ffile = TFile::Open("AliESDs.root", "RECREATE");
1474 ffile->SetCompressionLevel(2);
1475 if (!ffile->IsOpen()) {
1476 Abort("OpenESDFile", TSelector::kAbortProcess);
1481 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1482 Abort(Form("Problems opening output PROOF file: %s/%s",
1483 outProofFile->GetDir(), outProofFile->GetFileName()),
1484 TSelector::kAbortProcess);
1489 ftree = new TTree("esdTree", "Tree with ESD objects");
1490 fesd = new AliESDEvent();
1491 fesd->CreateStdContent();
1493 fesd->WriteToTree(ftree);
1494 if (fWriteESDfriend) {
1496 // Since we add the branch manually we must
1497 // book and add it after WriteToTree
1498 // otherwise it is created twice,
1499 // once via writetotree and once here.
1500 // The case for AliESDfriend is now
1501 // caught also in AlIESDEvent::WriteToTree but
1502 // be careful when changing the name (AliESDfriend is not
1503 // a TNamed so we had to hardwire it)
1504 fesdf = new AliESDfriend();
1505 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1506 br->SetFile("AliESDfriends.root");
1507 fesd->AddObject(fesdf);
1509 ftree->GetUserInfo()->Add(fesd);
1511 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1512 fhltesd = new AliESDEvent();
1513 fhltesd->CreateStdContent();
1515 // read the ESD template from CDB
1516 // HLT is allowed to put non-std content to its ESD, the non-std
1517 // objects need to be created before invocation of WriteToTree in
1518 // order to create all branches. Initialization is done from an
1519 // ESD layout template in CDB
1520 AliCDBManager* man = AliCDBManager::Instance();
1521 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1522 AliCDBEntry* hltESDConfig=NULL;
1523 if (man->GetId(hltESDConfigPath)!=NULL &&
1524 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1525 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1527 // init all internal variables from the list of objects
1528 pESDLayout->GetStdContent();
1530 // copy content and create non-std objects
1531 *fhltesd=*pESDLayout;
1534 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1535 hltESDConfigPath.GetPath().Data()));
1539 fhltesd->WriteToTree(fhlttree);
1540 fhlttree->GetUserInfo()->Add(fhltesd);
1542 ProcInfo_t procInfo;
1543 gSystem->GetProcInfo(&procInfo);
1544 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1547 //Initialize the QA and start of cycle
1548 if (fRunQA || fRunGlobalQA)
1551 //Initialize the Plane Efficiency framework
1552 if (fRunPlaneEff && !InitPlaneEff()) {
1553 Abort("InitPlaneEff", TSelector::kAbortProcess);
1557 if (strcmp(gProgName,"alieve") == 0)
1558 fRunAliEVE = InitAliEVE();
1563 //_____________________________________________________________________________
1564 Bool_t AliReconstruction::Process(Long64_t entry)
1566 // run the reconstruction over a single entry
1567 // from the chain with raw data
1568 AliCodeTimerAuto("");
1570 TTree *currTree = fChain->GetTree();
1571 AliRawVEvent *event = NULL;
1572 currTree->SetBranchAddress("rawevent",&event);
1573 currTree->GetEntry(entry);
1574 fRawReader = new AliRawReaderRoot(event);
1575 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1583 //_____________________________________________________________________________
1584 void AliReconstruction::Init(TTree *tree)
1587 AliError("The input tree is not found!");
1593 //_____________________________________________________________________________
1594 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1596 // run the reconstruction over a single event
1597 // The event loop is steered in Run method
1599 AliCodeTimerAuto("");
1601 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1602 fRunLoader->SetEventNumber(iEvent);
1603 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1605 fRunLoader->TreeE()->Fill();
1606 if (fRawReader && fRawReader->UseAutoSaveESD())
1607 fRunLoader->TreeE()->AutoSave("SaveSelf");
1610 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1614 AliInfo(Form("processing event %d", iEvent));
1616 fRunLoader->GetEvent(iEvent);
1618 // Fill Event-info object
1620 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo,fListOfCosmicTriggers);
1621 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1623 // Set the reco-params
1625 TString detStr = fLoadCDB;
1626 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1627 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1628 AliReconstructor *reconstructor = GetReconstructor(iDet);
1629 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1630 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1631 reconstructor->SetRecoParam(par);
1633 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1634 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(par->GetEventSpecie())) ;
1642 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1643 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1645 // local single event reconstruction
1646 if (!fRunLocalReconstruction.IsNull()) {
1647 TString detectors=fRunLocalReconstruction;
1648 // run HLT event reconstruction first
1649 // ;-( IsSelected changes the string
1650 if (IsSelected("HLT", detectors) &&
1651 !RunLocalEventReconstruction("HLT")) {
1652 if (fStopOnError) {CleanUp(); return kFALSE;}
1654 detectors=fRunLocalReconstruction;
1655 detectors.ReplaceAll("HLT", "");
1656 if (!RunLocalEventReconstruction(detectors)) {
1657 if (fStopOnError) {CleanUp(); return kFALSE;}
1661 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1662 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1663 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1664 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1666 // Set magnetic field from the tracker
1667 fesd->SetMagneticField(AliTracker::GetBz());
1668 fhltesd->SetMagneticField(AliTracker::GetBz());
1670 // Set most probable pt, for B=0 tracking
1671 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1672 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1673 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1675 // Fill raw-data error log into the ESD
1676 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1679 if (fRunVertexFinder) {
1680 if (!RunVertexFinder(fesd)) {
1681 if (fStopOnError) {CleanUp(); return kFALSE;}
1685 // For Plane Efficiency: run the SPD trackleter
1686 if (fRunPlaneEff && fSPDTrackleter) {
1687 if (!RunSPDTrackleting(fesd)) {
1688 if (fStopOnError) {CleanUp(); return kFALSE;}
1693 if (!fRunTracking.IsNull()) {
1694 if (fRunMuonTracking) {
1695 if (!RunMuonTracking(fesd)) {
1696 if (fStopOnError) {CleanUp(); return kFALSE;}
1702 if (!fRunTracking.IsNull()) {
1703 if (!RunTracking(fesd)) {
1704 if (fStopOnError) {CleanUp(); return kFALSE;}
1709 if (!fFillESD.IsNull()) {
1710 TString detectors=fFillESD;
1711 // run HLT first and on hltesd
1712 // ;-( IsSelected changes the string
1713 if (IsSelected("HLT", detectors) &&
1714 !FillESD(fhltesd, "HLT")) {
1715 if (fStopOnError) {CleanUp(); return kFALSE;}
1718 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1719 if (detectors.Contains("ALL")) {
1721 for (Int_t idet=0; idet<kNDetectors; ++idet){
1722 detectors += fgkDetectorName[idet];
1726 detectors.ReplaceAll("HLT", "");
1727 if (!FillESD(fesd, detectors)) {
1728 if (fStopOnError) {CleanUp(); return kFALSE;}
1732 // fill Event header information from the RawEventHeader
1733 if (fRawReader){FillRawEventHeaderESD(fesd);}
1736 AliESDpid::MakePID(fesd);
1738 if (fFillTriggerESD) {
1739 if (!FillTriggerESD(fesd)) {
1740 if (fStopOnError) {CleanUp(); return kFALSE;}
1747 // Propagate track to the beam pipe (if not already done by ITS)
1749 const Int_t ntracks = fesd->GetNumberOfTracks();
1750 const Double_t kBz = fesd->GetMagneticField();
1751 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1754 UShort_t *selectedIdx=new UShort_t[ntracks];
1756 for (Int_t itrack=0; itrack<ntracks; itrack++){
1757 const Double_t kMaxStep = 1; //max step over the material
1760 AliESDtrack *track = fesd->GetTrack(itrack);
1761 if (!track) continue;
1763 AliExternalTrackParam *tpcTrack =
1764 (AliExternalTrackParam *)track->GetTPCInnerParam();
1768 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1771 Int_t n=trkArray.GetEntriesFast();
1772 selectedIdx[n]=track->GetID();
1773 trkArray.AddLast(tpcTrack);
1776 //Tracks refitted by ITS should already be at the SPD vertex
1777 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1780 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1781 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1786 // Improve the reconstructed primary vertex position using the tracks
1788 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1789 if(fesd->GetPrimaryVertexSPD()) {
1790 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1791 if(vtitle.Contains("cosmics")) {
1792 runVertexFinderTracks=kFALSE;
1796 if (runVertexFinderTracks) {
1797 // TPC + ITS primary vertex
1798 ftVertexer->SetITSMode();
1799 ftVertexer->SetConstraintOff();
1800 // get cuts for vertexer from AliGRPRecoParam
1802 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1803 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1804 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1805 ftVertexer->SetCuts(cutsVertexer);
1806 delete [] cutsVertexer; cutsVertexer = NULL;
1807 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1808 ftVertexer->SetVtxStart(fDiamondProfile);
1810 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1812 if (pvtx->GetStatus()) {
1813 fesd->SetPrimaryVertexTracks(pvtx);
1814 for (Int_t i=0; i<ntracks; i++) {
1815 AliESDtrack *t = fesd->GetTrack(i);
1816 t->RelateToVertex(pvtx, kBz, kVeryBig);
1821 // TPC-only primary vertex
1822 ftVertexer->SetTPCMode();
1823 ftVertexer->SetConstraintOff();
1824 // get cuts for vertexer from AliGRPRecoParam
1826 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1827 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1828 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1829 ftVertexer->SetCuts(cutsVertexer);
1830 delete [] cutsVertexer; cutsVertexer = NULL;
1831 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1832 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1834 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1836 if (pvtx->GetStatus()) {
1837 fesd->SetPrimaryVertexTPC(pvtx);
1838 for (Int_t i=0; i<ntracks; i++) {
1839 AliESDtrack *t = fesd->GetTrack(i);
1840 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1846 delete[] selectedIdx;
1848 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1853 AliV0vertexer vtxer;
1854 vtxer.Tracks2V0vertices(fesd);
1856 if (fRunCascadeFinder) {
1858 AliCascadeVertexer cvtxer;
1859 cvtxer.V0sTracks2CascadeVertices(fesd);
1864 if (fCleanESD) CleanESD(fesd);
1867 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1868 AliQAManager::QAManager()->RunOneEvent(fesd) ;
1871 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1872 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1873 if (qadm && fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1874 qadm->Exec(AliQAv1::kESDS, fesd);
1877 if (fWriteESDfriend) {
1878 // fesdf->~AliESDfriend();
1879 // new (fesdf) AliESDfriend(); // Reset...
1880 fesd->GetESDfriend(fesdf);
1884 // Auto-save the ESD tree in case of prompt reco @P2
1885 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1886 ftree->AutoSave("SaveSelf");
1887 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1888 if (friendfile) friendfile->Save();
1895 if (fRunAliEVE) RunAliEVE();
1899 if (fWriteESDfriend) {
1900 fesdf->~AliESDfriend();
1901 new (fesdf) AliESDfriend(); // Reset...
1904 ProcInfo_t procInfo;
1905 gSystem->GetProcInfo(&procInfo);
1906 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1909 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1910 if (fReconstructor[iDet])
1911 fReconstructor[iDet]->SetRecoParam(NULL);
1914 if (fRunQA || fRunGlobalQA)
1915 AliQAManager::QAManager()->Increment() ;
1920 //_____________________________________________________________________________
1921 void AliReconstruction::SlaveTerminate()
1923 // Finalize the run on the slave side
1924 // Called after the exit
1925 // from the event loop
1926 AliCodeTimerAuto("");
1928 if (fIsNewRunLoader) { // galice.root didn't exist
1929 fRunLoader->WriteHeader("OVERWRITE");
1930 fRunLoader->CdGAFile();
1931 fRunLoader->Write(0, TObject::kOverwrite);
1934 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1935 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1937 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1938 cdbMapCopy->SetOwner(1);
1939 cdbMapCopy->SetName("cdbMap");
1940 TIter iter(cdbMap->GetTable());
1943 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1944 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1945 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1946 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1949 TList *cdbListCopy = new TList();
1950 cdbListCopy->SetOwner(1);
1951 cdbListCopy->SetName("cdbList");
1953 TIter iter2(cdbList);
1956 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1957 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1960 ftree->GetUserInfo()->Add(cdbMapCopy);
1961 ftree->GetUserInfo()->Add(cdbListCopy);
1966 if (fWriteESDfriend)
1967 ftree->SetBranchStatus("ESDfriend*",0);
1968 // we want to have only one tree version number
1969 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1972 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1973 if (fRunPlaneEff && !FinishPlaneEff()) {
1974 AliWarning("Finish PlaneEff evaluation failed");
1977 // End of cycle for the in-loop
1979 AliQAManager::QAManager()->EndOfCycle() ;
1982 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1984 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS)))
1985 qadm->EndOfCycle(AliQAv1::kRECPOINTS);
1986 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1987 qadm->EndOfCycle(AliQAv1::kESDS);
1992 if (fRunQA || fRunGlobalQA) {
1994 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1995 TString qaOutputFile = outputFileName->GetTitle();
1996 qaOutputFile.ReplaceAll(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()),
1997 Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
1998 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
1999 qaProofFile->SetOutputFileName(qaOutputFile.Data());
2000 fOutput->Add(qaProofFile);
2001 MergeQA(qaProofFile->GetFileName());
2013 //_____________________________________________________________________________
2014 void AliReconstruction::Terminate()
2016 // Create tags for the events in the ESD tree (the ESD tree is always present)
2017 // In case of empty events the tags will contain dummy values
2018 AliCodeTimerAuto("");
2020 // Do not call the ESD tag creator in case of PROOF-based reconstruction
2022 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
2023 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
2026 // Cleanup of CDB manager: cache and active storages!
2027 AliCDBManager::Instance()->ClearCache();
2030 //_____________________________________________________________________________
2031 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
2033 // run the local reconstruction
2035 static Int_t eventNr=0;
2036 AliCodeTimerAuto("")
2038 TString detStr = detectors;
2039 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2040 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2041 AliReconstructor* reconstructor = GetReconstructor(iDet);
2042 if (!reconstructor) continue;
2043 AliLoader* loader = fLoader[iDet];
2044 // Matthias April 2008: temporary fix to run HLT reconstruction
2045 // although the HLT loader is missing
2046 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2048 reconstructor->Reconstruct(fRawReader, NULL);
2051 reconstructor->Reconstruct(dummy, NULL);
2056 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2059 // conversion of digits
2060 if (fRawReader && reconstructor->HasDigitConversion()) {
2061 AliInfo(Form("converting raw data digits into root objects for %s",
2062 fgkDetectorName[iDet]));
2063 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2064 // fgkDetectorName[iDet]));
2065 loader->LoadDigits("update");
2066 loader->CleanDigits();
2067 loader->MakeDigitsContainer();
2068 TTree* digitsTree = loader->TreeD();
2069 reconstructor->ConvertDigits(fRawReader, digitsTree);
2070 loader->WriteDigits("OVERWRITE");
2071 loader->UnloadDigits();
2073 // local reconstruction
2074 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2075 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2076 loader->LoadRecPoints("update");
2077 loader->CleanRecPoints();
2078 loader->MakeRecPointsContainer();
2079 TTree* clustersTree = loader->TreeR();
2080 if (fRawReader && !reconstructor->HasDigitConversion()) {
2081 reconstructor->Reconstruct(fRawReader, clustersTree);
2083 loader->LoadDigits("read");
2084 TTree* digitsTree = loader->TreeD();
2086 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2087 if (fStopOnError) return kFALSE;
2089 reconstructor->Reconstruct(digitsTree, clustersTree);
2091 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2092 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2095 loader->UnloadDigits();
2098 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2099 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2101 loader->WriteRecPoints("OVERWRITE");
2102 loader->UnloadRecPoints();
2103 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2105 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2106 AliError(Form("the following detectors were not found: %s",
2108 if (fStopOnError) return kFALSE;
2113 //_____________________________________________________________________________
2114 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2116 // run the SPD trackleting (for SPD efficiency purpouses)
2118 AliCodeTimerAuto("")
2120 Double_t vtxPos[3] = {0, 0, 0};
2121 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2123 TArrayF mcVertex(3);
2125 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2126 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2127 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2130 const AliESDVertex *vertex = esd->GetVertex();
2132 AliWarning("Vertex not found");
2135 vertex->GetXYZ(vtxPos);
2136 vertex->GetSigmaXYZ(vtxErr);
2137 if (fSPDTrackleter) {
2138 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2141 fLoader[0]->LoadRecPoints("read");
2142 TTree* tree = fLoader[0]->TreeR();
2144 AliError("Can't get the ITS cluster tree");
2147 fSPDTrackleter->LoadClusters(tree);
2148 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2150 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2151 AliError("AliITSTrackleterSPDEff Clusters2Tracks failed");
2152 // fLoader[0]->UnloadRecPoints();
2155 //fSPDTrackleter->UnloadRecPoints();
2157 AliWarning("SPDTrackleter not available");
2163 //_____________________________________________________________________________
2164 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2166 // run the barrel tracking
2168 AliCodeTimerAuto("")
2170 AliVertexer *vertexer = CreateVertexer();
2171 if (!vertexer) return kFALSE;
2173 AliInfo("running the ITS vertex finder");
2174 AliESDVertex* vertex = NULL;
2176 fLoader[0]->LoadRecPoints();
2177 TTree* cltree = fLoader[0]->TreeR();
2179 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2180 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2183 AliError("Can't get the ITS cluster tree");
2185 fLoader[0]->UnloadRecPoints();
2188 AliError("Can't get the ITS loader");
2191 AliWarning("Vertex not found");
2192 vertex = new AliESDVertex();
2193 vertex->SetName("default");
2196 vertex->SetName("reconstructed");
2201 vertex->GetXYZ(vtxPos);
2202 vertex->GetSigmaXYZ(vtxErr);
2204 esd->SetPrimaryVertexSPD(vertex);
2205 AliESDVertex *vpileup = NULL;
2206 Int_t novertices = 0;
2207 vpileup = vertexer->GetAllVertices(novertices);
2209 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2211 // if SPD multiplicity has been determined, it is stored in the ESD
2212 AliMultiplicity *mult = vertexer->GetMultiplicity();
2213 if(mult)esd->SetMultiplicity(mult);
2215 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2216 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2225 //_____________________________________________________________________________
2226 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2228 // run the HLT barrel tracking
2230 AliCodeTimerAuto("")
2233 AliError("Missing runLoader!");
2237 AliInfo("running HLT tracking");
2239 // Get a pointer to the HLT reconstructor
2240 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2241 if (!reconstructor) return kFALSE;
2244 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2245 TString detName = fgkDetectorName[iDet];
2246 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2247 reconstructor->SetOption(detName.Data());
2248 AliTracker *tracker = reconstructor->CreateTracker();
2250 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2251 if (fStopOnError) return kFALSE;
2255 Double_t vtxErr[3]={0.005,0.005,0.010};
2256 const AliESDVertex *vertex = esd->GetVertex();
2257 vertex->GetXYZ(vtxPos);
2258 tracker->SetVertex(vtxPos,vtxErr);
2260 fLoader[iDet]->LoadRecPoints("read");
2261 TTree* tree = fLoader[iDet]->TreeR();
2263 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2266 tracker->LoadClusters(tree);
2268 if (tracker->Clusters2Tracks(esd) != 0) {
2269 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2273 tracker->UnloadClusters();
2281 //_____________________________________________________________________________
2282 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2284 // run the muon spectrometer tracking
2286 AliCodeTimerAuto("")
2289 AliError("Missing runLoader!");
2292 Int_t iDet = 7; // for MUON
2294 AliInfo("is running...");
2296 // Get a pointer to the MUON reconstructor
2297 AliReconstructor *reconstructor = GetReconstructor(iDet);
2298 if (!reconstructor) return kFALSE;
2301 TString detName = fgkDetectorName[iDet];
2302 AliDebug(1, Form("%s tracking", detName.Data()));
2303 AliTracker *tracker = reconstructor->CreateTracker();
2305 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2310 fLoader[iDet]->LoadRecPoints("read");
2312 tracker->LoadClusters(fLoader[iDet]->TreeR());
2314 Int_t rv = tracker->Clusters2Tracks(esd);
2318 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2322 fLoader[iDet]->UnloadRecPoints();
2324 tracker->UnloadClusters();
2332 //_____________________________________________________________________________
2333 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2335 // run the barrel tracking
2336 static Int_t eventNr=0;
2337 AliCodeTimerAuto("")
2339 AliInfo("running tracking");
2342 //Fill the ESD with the T0 info (will be used by the TOF)
2343 if (fReconstructor[11] && fLoader[11]) {
2344 fLoader[11]->LoadRecPoints("READ");
2345 TTree *treeR = fLoader[11]->TreeR();
2347 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2351 // pass 1: TPC + ITS inwards
2352 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2353 if (!fTracker[iDet]) continue;
2354 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2357 fLoader[iDet]->LoadRecPoints("read");
2358 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2359 TTree* tree = fLoader[iDet]->TreeR();
2361 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2364 fTracker[iDet]->LoadClusters(tree);
2365 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2367 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2368 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2371 // preliminary PID in TPC needed by the ITS tracker
2373 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2374 AliESDpid::MakePID(esd);
2376 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2379 // pass 2: ALL backwards
2381 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2382 if (!fTracker[iDet]) continue;
2383 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2386 if (iDet > 1) { // all except ITS, TPC
2388 fLoader[iDet]->LoadRecPoints("read");
2389 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2390 tree = fLoader[iDet]->TreeR();
2392 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2395 fTracker[iDet]->LoadClusters(tree);
2396 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2400 if (iDet>1) // start filling residuals for the "outer" detectors
2402 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2403 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2405 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2406 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2407 if ( elem && (! elem->At(0)) ) {
2408 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2409 if (qadm) qadm->InitRecPointsForTracker() ;
2413 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2414 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2419 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2420 fTracker[iDet]->UnloadClusters();
2421 fLoader[iDet]->UnloadRecPoints();
2423 // updated PID in TPC needed by the ITS tracker -MI
2425 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2426 AliESDpid::MakePID(esd);
2428 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2430 //stop filling residuals for the "outer" detectors
2431 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2433 // pass 3: TRD + TPC + ITS refit inwards
2435 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2436 if (!fTracker[iDet]) continue;
2437 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2440 if (iDet<2) // start filling residuals for TPC and ITS
2442 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2443 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2445 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2446 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2447 if ( elem && (! elem->At(0)) ) {
2448 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2449 if (qadm) qadm->InitRecPointsForTracker() ;
2454 if (fTracker[iDet]->RefitInward(esd) != 0) {
2455 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2458 // run postprocessing
2459 if (fTracker[iDet]->PostProcess(esd) != 0) {
2460 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2463 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2466 // write space-points to the ESD in case alignment data output
2468 if (fWriteAlignmentData)
2469 WriteAlignmentData(esd);
2471 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2472 if (!fTracker[iDet]) continue;
2474 fTracker[iDet]->UnloadClusters();
2475 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2476 fLoader[iDet]->UnloadRecPoints();
2477 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2479 // stop filling residuals for TPC and ITS
2480 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2486 //_____________________________________________________________________________
2487 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2489 // Remove the data which are not needed for the physics analysis.
2492 Int_t nTracks=esd->GetNumberOfTracks();
2493 Int_t nV0s=esd->GetNumberOfV0s();
2495 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2497 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2498 Bool_t rc=esd->Clean(cleanPars);
2500 nTracks=esd->GetNumberOfTracks();
2501 nV0s=esd->GetNumberOfV0s();
2503 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2508 //_____________________________________________________________________________
2509 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2511 // fill the event summary data
2513 AliCodeTimerAuto("")
2514 static Int_t eventNr=0;
2515 TString detStr = detectors;
2517 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2518 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2519 AliReconstructor* reconstructor = GetReconstructor(iDet);
2520 if (!reconstructor) continue;
2521 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2522 TTree* clustersTree = NULL;
2523 if (fLoader[iDet]) {
2524 fLoader[iDet]->LoadRecPoints("read");
2525 clustersTree = fLoader[iDet]->TreeR();
2526 if (!clustersTree) {
2527 AliError(Form("Can't get the %s clusters tree",
2528 fgkDetectorName[iDet]));
2529 if (fStopOnError) return kFALSE;
2532 if (fRawReader && !reconstructor->HasDigitConversion()) {
2533 reconstructor->FillESD(fRawReader, clustersTree, esd);
2535 TTree* digitsTree = NULL;
2536 if (fLoader[iDet]) {
2537 fLoader[iDet]->LoadDigits("read");
2538 digitsTree = fLoader[iDet]->TreeD();
2540 AliError(Form("Can't get the %s digits tree",
2541 fgkDetectorName[iDet]));
2542 if (fStopOnError) return kFALSE;
2545 reconstructor->FillESD(digitsTree, clustersTree, esd);
2546 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2548 if (fLoader[iDet]) {
2549 fLoader[iDet]->UnloadRecPoints();
2553 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2554 AliError(Form("the following detectors were not found: %s",
2556 if (fStopOnError) return kFALSE;
2558 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2563 //_____________________________________________________________________________
2564 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2566 // Reads the trigger decision which is
2567 // stored in Trigger.root file and fills
2568 // the corresponding esd entries
2570 AliCodeTimerAuto("")
2572 AliInfo("Filling trigger information into the ESD");
2575 AliCTPRawStream input(fRawReader);
2576 if (!input.Next()) {
2577 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2580 if (esd->GetTriggerMask() != input.GetClassMask())
2581 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2582 input.GetClassMask(),esd->GetTriggerMask()));
2583 if (esd->GetOrbitNumber() != input.GetOrbitID())
2584 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2585 input.GetOrbitID(),esd->GetOrbitNumber()));
2586 if (esd->GetBunchCrossNumber() != input.GetBCID())
2587 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2588 input.GetBCID(),esd->GetBunchCrossNumber()));
2589 AliESDHeader* esdheader = esd->GetHeader();
2590 esdheader->SetL0TriggerInputs(input.GetL0Inputs());
2591 esdheader->SetL1TriggerInputs(input.GetL1Inputs());
2592 esdheader->SetL2TriggerInputs(input.GetL2Inputs());
2594 UInt_t orbit=input.GetOrbitID();
2595 for(Int_t i=0 ; i<input.GetNIRs() ; i++ )
2596 if(TMath::Abs(Int_t(orbit-(input.GetIR(i))->GetOrbit()))<=1){
2597 esdheader->AddTriggerIR(input.GetIR(i));
2603 AliTimeStamp* timestamp = new AliTimeStamp(esd->GetOrbitNumber(), esd->GetPeriodNumber(), esd->GetBunchCrossNumber());
2604 //AliTimeStamp* timestamp = new AliTimeStamp(14322992, 5, (ULong64_t)486238);
2605 AliESDHeader* esdheader = fesd->GetHeader();
2606 for(Int_t i=0;i<50;i++){
2607 if((1<<i) & esd->GetTriggerMask()){
2608 AliTriggerScalersESD* scalesd = fRunScalers->GetScalersForEventClass( timestamp, i);
2609 esdheader->SetTriggerScalersRecord(scalesd);
2620 //_____________________________________________________________________________
2621 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2624 // Filling information from RawReader Header
2627 if (!fRawReader) return kFALSE;
2629 AliInfo("Filling information from RawReader Header");
2631 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2632 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2633 esd->SetPeriodNumber(fRawReader->GetPeriod());
2635 esd->SetTimeStamp(fRawReader->GetTimestamp());
2636 esd->SetEventType(fRawReader->GetType());
2642 //_____________________________________________________________________________
2643 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2645 // check whether detName is contained in detectors
2646 // if yes, it is removed from detectors
2648 // check if all detectors are selected
2649 if ((detectors.CompareTo("ALL") == 0) ||
2650 detectors.BeginsWith("ALL ") ||
2651 detectors.EndsWith(" ALL") ||
2652 detectors.Contains(" ALL ")) {
2657 // search for the given detector
2658 Bool_t result = kFALSE;
2659 if ((detectors.CompareTo(detName) == 0) ||
2660 detectors.BeginsWith(detName+" ") ||
2661 detectors.EndsWith(" "+detName) ||
2662 detectors.Contains(" "+detName+" ")) {
2663 detectors.ReplaceAll(detName, "");
2667 // clean up the detectors string
2668 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2669 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2670 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2675 //_____________________________________________________________________________
2676 Bool_t AliReconstruction::InitRunLoader()
2678 // get or create the run loader
2680 if (gAlice) delete gAlice;
2683 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2684 // load all base libraries to get the loader classes
2685 TString libs = gSystem->GetLibraries();
2686 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2687 TString detName = fgkDetectorName[iDet];
2688 if (detName == "HLT") continue;
2689 if (libs.Contains("lib" + detName + "base.so")) continue;
2690 gSystem->Load("lib" + detName + "base.so");
2692 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2694 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2699 fRunLoader->CdGAFile();
2700 fRunLoader->LoadgAlice();
2702 //PH This is a temporary fix to give access to the kinematics
2703 //PH that is needed for the labels of ITS clusters
2704 fRunLoader->LoadHeader();
2705 fRunLoader->LoadKinematics();
2707 } else { // galice.root does not exist
2709 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2711 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2712 AliConfig::GetDefaultEventFolderName(),
2715 AliError(Form("could not create run loader in file %s",
2716 fGAliceFileName.Data()));
2720 fIsNewRunLoader = kTRUE;
2721 fRunLoader->MakeTree("E");
2723 if (fNumberOfEventsPerFile > 0)
2724 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2726 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2732 //_____________________________________________________________________________
2733 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2735 // get the reconstructor object and the loader for a detector
2737 if (fReconstructor[iDet]) {
2738 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2739 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2740 fReconstructor[iDet]->SetRecoParam(par);
2742 return fReconstructor[iDet];
2745 // load the reconstructor object
2746 TPluginManager* pluginManager = gROOT->GetPluginManager();
2747 TString detName = fgkDetectorName[iDet];
2748 TString recName = "Ali" + detName + "Reconstructor";
2750 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2752 AliReconstructor* reconstructor = NULL;
2753 // first check if a plugin is defined for the reconstructor
2754 TPluginHandler* pluginHandler =
2755 pluginManager->FindHandler("AliReconstructor", detName);
2756 // if not, add a plugin for it
2757 if (!pluginHandler) {
2758 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2759 TString libs = gSystem->GetLibraries();
2760 if (libs.Contains("lib" + detName + "base.so") ||
2761 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2762 pluginManager->AddHandler("AliReconstructor", detName,
2763 recName, detName + "rec", recName + "()");
2765 pluginManager->AddHandler("AliReconstructor", detName,
2766 recName, detName, recName + "()");
2768 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2770 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2771 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2773 if (reconstructor) {
2774 TObject* obj = fOptions.FindObject(detName.Data());
2775 if (obj) reconstructor->SetOption(obj->GetTitle());
2776 reconstructor->Init();
2777 fReconstructor[iDet] = reconstructor;
2780 // get or create the loader
2781 if (detName != "HLT") {
2782 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2783 if (!fLoader[iDet]) {
2784 AliConfig::Instance()
2785 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2787 // first check if a plugin is defined for the loader
2789 pluginManager->FindHandler("AliLoader", detName);
2790 // if not, add a plugin for it
2791 if (!pluginHandler) {
2792 TString loaderName = "Ali" + detName + "Loader";
2793 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2794 pluginManager->AddHandler("AliLoader", detName,
2795 loaderName, detName + "base",
2796 loaderName + "(const char*, TFolder*)");
2797 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2799 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2801 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2802 fRunLoader->GetEventFolder());
2804 if (!fLoader[iDet]) { // use default loader
2805 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2807 if (!fLoader[iDet]) {
2808 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2809 if (fStopOnError) return NULL;
2811 fRunLoader->AddLoader(fLoader[iDet]);
2812 fRunLoader->CdGAFile();
2813 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2814 fRunLoader->Write(0, TObject::kOverwrite);
2819 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2820 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2821 reconstructor->SetRecoParam(par);
2823 return reconstructor;
2826 //_____________________________________________________________________________
2827 AliVertexer* AliReconstruction::CreateVertexer()
2829 // create the vertexer
2830 // Please note that the caller is the owner of the
2833 AliVertexer* vertexer = NULL;
2834 AliReconstructor* itsReconstructor = GetReconstructor(0);
2835 if (itsReconstructor && (fRunLocalReconstruction.Contains("ITS"))) {
2836 vertexer = itsReconstructor->CreateVertexer();
2839 AliWarning("couldn't create a vertexer for ITS");
2845 //_____________________________________________________________________________
2846 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2848 // create the trackers
2849 AliInfo("Creating trackers");
2851 TString detStr = detectors;
2852 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2853 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2854 AliReconstructor* reconstructor = GetReconstructor(iDet);
2855 if (!reconstructor) continue;
2856 TString detName = fgkDetectorName[iDet];
2857 if (detName == "HLT") {
2858 fRunHLTTracking = kTRUE;
2861 if (detName == "MUON") {
2862 fRunMuonTracking = kTRUE;
2867 fTracker[iDet] = reconstructor->CreateTracker();
2868 if (!fTracker[iDet] && (iDet < 7)) {
2869 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2870 if (fStopOnError) return kFALSE;
2872 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2878 //_____________________________________________________________________________
2879 void AliReconstruction::CleanUp()
2881 // delete trackers and the run loader and close and delete the file
2883 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2884 delete fReconstructor[iDet];
2885 fReconstructor[iDet] = NULL;
2886 fLoader[iDet] = NULL;
2887 delete fTracker[iDet];
2888 fTracker[iDet] = NULL;
2893 delete fSPDTrackleter;
2894 fSPDTrackleter = NULL;
2903 delete fParentRawReader;
2904 fParentRawReader=NULL;
2912 if (AliQAManager::QAManager())
2913 AliQAManager::QAManager()->ShowQA() ;
2914 AliQAManager::Destroy() ;
2916 TGeoGlobalMagField::Instance()->SetField(NULL);
2919 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2921 // Write space-points which are then used in the alignment procedures
2922 // For the moment only ITS, TPC, TRD and TOF
2924 Int_t ntracks = esd->GetNumberOfTracks();
2925 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2927 AliESDtrack *track = esd->GetTrack(itrack);
2930 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
2931 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2932 nsp += track->GetNcls(iDet);
2934 if (iDet==0) { // ITS "extra" clusters
2935 track->GetClusters(iDet,idx);
2936 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2941 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2942 track->SetTrackPointArray(sp);
2944 for (Int_t iDet = 5; iDet >= 0; iDet--) {
2945 AliTracker *tracker = fTracker[iDet];
2946 if (!tracker) continue;
2947 Int_t nspdet = track->GetClusters(iDet,idx);
2949 if (iDet==0) // ITS "extra" clusters
2950 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2952 if (nspdet <= 0) continue;
2956 while (isp2 < nspdet) {
2957 Bool_t isvalid=kTRUE;
2959 Int_t index=idx[isp++];
2960 if (index < 0) continue;
2962 TString dets = fgkDetectorName[iDet];
2963 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2964 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2965 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2966 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2967 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2969 isvalid = tracker->GetTrackPoint(index,p);
2972 if (!isvalid) continue;
2973 if (iDet==0 && (isp-1)>=6) p.SetExtra();
2974 sp->AddPoint(isptrack,&p); isptrack++;
2981 //_____________________________________________________________________________
2982 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2984 // The method reads the raw-data error log
2985 // accumulated within the rawReader.
2986 // It extracts the raw-data errors related to
2987 // the current event and stores them into
2988 // a TClonesArray inside the esd object.
2990 if (!fRawReader) return;
2992 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2994 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2996 if (iEvent != log->GetEventNumber()) continue;
2998 esd->AddRawDataErrorLog(log);
3003 //_____________________________________________________________________________
3004 void AliReconstruction::CheckQA()
3006 // check the QA of SIM for this run and remove the detectors
3007 // with status Fatal
3009 // TString newRunLocalReconstruction ;
3010 // TString newRunTracking ;
3011 // TString newFillESD ;
3013 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
3014 // TString detName(AliQAv1::GetDetName(iDet)) ;
3015 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
3016 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
3017 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
3018 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
3020 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
3021 // fRunLocalReconstruction.Contains("ALL") ) {
3022 // newRunLocalReconstruction += detName ;
3023 // newRunLocalReconstruction += " " ;
3025 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
3026 // fRunTracking.Contains("ALL") ) {
3027 // newRunTracking += detName ;
3028 // newRunTracking += " " ;
3030 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
3031 // fFillESD.Contains("ALL") ) {
3032 // newFillESD += detName ;
3033 // newFillESD += " " ;
3037 // fRunLocalReconstruction = newRunLocalReconstruction ;
3038 // fRunTracking = newRunTracking ;
3039 // fFillESD = newFillESD ;
3042 //_____________________________________________________________________________
3043 Int_t AliReconstruction::GetDetIndex(const char* detector)
3045 // return the detector index corresponding to detector
3047 for (index = 0; index < kNDetectors ; index++) {
3048 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3053 //_____________________________________________________________________________
3054 Bool_t AliReconstruction::FinishPlaneEff() {
3056 // Here execute all the necessary operationis, at the end of the tracking phase,
3057 // in case that evaluation of PlaneEfficiencies was required for some detector.
3058 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3060 // This Preliminary version works only FOR ITS !!!!!
3061 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3064 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3067 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3068 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3069 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3070 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
3071 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3072 TString name=planeeff->GetName();
3074 TFile* pefile = TFile::Open(name, "RECREATE");
3075 ret=(Bool_t)planeeff->Write();
3077 if(planeeff->GetCreateHistos()) {
3078 TString hname=planeeff->GetName();
3079 hname+="Histo.root";
3080 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3083 if(fSPDTrackleter) {
3084 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3085 TString name="AliITSPlaneEffSPDtracklet.root";
3086 TFile* pefile = TFile::Open(name, "RECREATE");
3087 ret=(Bool_t)planeeff->Write();
3089 AliESDEvent *dummy=NULL;
3090 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3095 //_____________________________________________________________________________
3096 Bool_t AliReconstruction::InitPlaneEff() {
3098 // Here execute all the necessary operations, before of the tracking phase,
3099 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3100 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3101 // which should be updated/recalculated.
3103 // This Preliminary version will work only FOR ITS !!!!!
3104 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3107 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3109 AliWarning(Form("Implementation of this method not yet completed !! Method return kTRUE"));
3111 fSPDTrackleter = NULL;
3112 AliReconstructor* itsReconstructor = GetReconstructor(0);
3113 if (itsReconstructor) {
3114 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3116 if (fSPDTrackleter) {
3117 AliInfo("Trackleter for SPD has been created");
3123 //_____________________________________________________________________________
3124 Bool_t AliReconstruction::InitAliEVE()
3126 // This method should be called only in case
3127 // AliReconstruction is run
3128 // within the alieve environment.
3129 // It will initialize AliEVE in a way
3130 // so that it can visualize event processed
3131 // by AliReconstruction.
3132 // The return flag shows whenever the
3133 // AliEVE initialization was successful or not.
3136 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3137 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3138 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3140 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3141 gROOT->ProcessLine("alieve_online_init()");
3146 //_____________________________________________________________________________
3147 void AliReconstruction::RunAliEVE()
3149 // Runs AliEVE visualisation of
3150 // the current event.
3151 // Should be executed only after
3152 // successful initialization of AliEVE.
3154 AliInfo("Running AliEVE...");
3155 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3159 //_____________________________________________________________________________
3160 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3162 // Allows to run QA for a selected set of detectors
3163 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3164 // all selected detectors run the same selected tasks
3166 if (!detAndAction.Contains(":")) {
3167 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3171 Int_t colon = detAndAction.Index(":") ;
3172 fQADetectors = detAndAction(0, colon) ;
3173 if (fQADetectors.Contains("ALL") )
3174 fQADetectors = fFillESD ;
3175 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3176 if (fQATasks.Contains("ALL") ) {
3177 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3179 fQATasks.ToUpper() ;
3181 if ( fQATasks.Contains("RAW") )
3182 tempo = Form("%d ", AliQAv1::kRAWS) ;
3183 if ( fQATasks.Contains("DIGIT") )
3184 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3185 if ( fQATasks.Contains("RECPOINT") )
3186 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3187 if ( fQATasks.Contains("ESD") )
3188 tempo += Form("%d ", AliQAv1::kESDS) ;
3190 if (fQATasks.IsNull()) {
3191 AliInfo("No QA requested\n") ;
3196 TString tempo(fQATasks) ;
3197 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3198 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3199 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3200 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3201 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3206 //_____________________________________________________________________________
3207 Bool_t AliReconstruction::InitRecoParams()
3209 // The method accesses OCDB and retrieves all
3210 // the available reco-param objects from there.
3212 Bool_t isOK = kTRUE;
3214 TString detStr = fLoadCDB;
3215 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3217 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3219 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3220 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3224 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3226 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3227 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3229 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3233 TObject *recoParamObj = entry->GetObject();
3234 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3235 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3236 // Registering them in AliRecoParam
3237 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3239 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3240 // The detector has only onse set of reco parameters
3241 // Registering it in AliRecoParam
3242 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3243 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3244 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3247 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3251 AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3255 if (AliDebugLevel() > 0) fRecoParam.Print();
3260 //_____________________________________________________________________________
3261 Bool_t AliReconstruction::GetEventInfo()
3263 // Fill the event info object
3265 AliCodeTimerAuto("")
3267 AliCentralTrigger *aCTP = NULL;
3269 fEventInfo.SetEventType(fRawReader->GetType());
3271 ULong64_t mask = fRawReader->GetClassMask();
3272 fEventInfo.SetTriggerMask(mask);
3273 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3274 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3276 aCTP = new AliCentralTrigger();
3277 TString configstr("");
3278 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3279 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3283 aCTP->SetClassMask(mask);
3284 aCTP->SetClusterMask(clmask);
3287 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3289 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3290 aCTP = fRunLoader->GetTrigger();
3291 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3292 // get inputs from actp - just get
3293 AliESDHeader* esdheader = fesd->GetHeader();
3294 esdheader->SetL0TriggerInputs(aCTP->GetL0TriggerInputs());
3295 esdheader->SetL1TriggerInputs(aCTP->GetL1TriggerInputs());
3296 esdheader->SetL2TriggerInputs(aCTP->GetL2TriggerInputs());
3297 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3300 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3305 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3307 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3308 if (fRawReader) delete aCTP;
3312 UChar_t clustmask = 0;
3314 ULong64_t trmask = fEventInfo.GetTriggerMask();
3315 const TObjArray& classesArray = config->GetClasses();
3316 Int_t nclasses = classesArray.GetEntriesFast();
3317 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3318 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3320 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3321 fesd->SetTriggerClass(trclass->GetName(),trindex);
3322 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3323 if (trmask & (1ull << trindex)) {
3325 trclasses += trclass->GetName();
3327 clustmask |= trclass->GetCluster()->GetClusterMask();
3331 fEventInfo.SetTriggerClasses(trclasses);
3333 // Set the information in ESD
3334 fesd->SetTriggerMask(trmask);
3335 fesd->SetTriggerCluster(clustmask);
3337 if (!aCTP->CheckTriggeredDetectors()) {
3338 if (fRawReader) delete aCTP;
3342 if (fRawReader) delete aCTP;
3344 // We have to fill also the HLT decision here!!
3350 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3352 // Match the detector list found in the rec.C or the default 'ALL'
3353 // to the list found in the GRP (stored there by the shuttle PP which
3354 // gets the information from ECS)
3355 static TString resultList;
3356 TString detList = detectorList;
3360 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3361 if ((detectorMask >> iDet) & 0x1) {
3362 TString det = AliDAQ::OfflineModuleName(iDet);
3363 if ((detList.CompareTo("ALL") == 0) ||
3364 ((detList.BeginsWith("ALL ") ||
3365 detList.EndsWith(" ALL") ||
3366 detList.Contains(" ALL ")) &&
3367 !(detList.BeginsWith("-"+det+" ") ||
3368 detList.EndsWith(" -"+det) ||
3369 detList.Contains(" -"+det+" "))) ||
3370 (detList.CompareTo(det) == 0) ||
3371 detList.BeginsWith(det+" ") ||
3372 detList.EndsWith(" "+det) ||
3373 detList.Contains( " "+det+" " )) {
3374 if (!resultList.EndsWith(det + " ")) {
3383 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3384 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3385 if ((detList.CompareTo("ALL") == 0) ||
3386 ((detList.BeginsWith("ALL ") ||
3387 detList.EndsWith(" ALL") ||
3388 detList.Contains(" ALL ")) &&
3389 !(detList.BeginsWith("-"+hltDet+" ") ||
3390 detList.EndsWith(" -"+hltDet) ||
3391 detList.Contains(" -"+hltDet+" "))) ||
3392 (detList.CompareTo(hltDet) == 0) ||
3393 detList.BeginsWith(hltDet+" ") ||
3394 detList.EndsWith(" "+hltDet) ||
3395 detList.Contains( " "+hltDet+" " )) {
3396 resultList += hltDet;
3400 return resultList.Data();
3404 //______________________________________________________________________________
3405 void AliReconstruction::Abort(const char *method, EAbort what)
3407 // Abort processing. If what = kAbortProcess, the Process() loop will be
3408 // aborted. If what = kAbortFile, the current file in a chain will be
3409 // aborted and the processing will continue with the next file, if there
3410 // is no next file then Process() will be aborted. Abort() can also be
3411 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3412 // the SlaveTerminate() and Terminate() are always called. The abort flag
3413 // can be checked in these methods using GetAbort().
3415 // The method is overwritten in AliReconstruction for better handling of
3416 // reco specific errors
3418 if (!fStopOnError) return;
3422 TString whyMess = method;
3423 whyMess += " failed! Aborting...";
3425 AliError(whyMess.Data());
3428 TString mess = "Abort";
3429 if (fAbort == kAbortProcess)
3430 mess = "AbortProcess";
3431 else if (fAbort == kAbortFile)
3434 Info(mess, whyMess.Data());