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(AliQAv1::kRECMODE) ;
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, Int_t convention, Bool_t uniform,
911 Float_t beamenergy, 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 match the convention (LHC or DCS2008)
919 // unless the special uniform map was used for MC
920 //------------------------------------------------
921 const Float_t l3NominalCurrent1=30000.; // (A)
922 const Float_t l3NominalCurrent2=12000.; // (A)
923 const Float_t diNominalCurrent =6000. ; // (A)
925 const Float_t tolerance=0.03; // relative current tolerance
926 const Float_t zero=77.; // "zero" current (A)
931 l3Cur = TMath::Abs(l3Cur);
932 diCur = TMath::Abs(diCur);
934 if (TMath::Abs((sclDip=diCur/diNominalCurrent)-1.) > tolerance && !uniform) {
935 if (diCur <= zero) sclDip = 0.; // some small current.. -> Dipole OFF
937 AliError(Form("Wrong dipole current (%f A)!",diCur));
943 // special treatment of special MC with uniform mag field (normalized to 0.5 T)
944 // no check for scaling/polarities are done
945 map = AliMagF::k5kGUniform;
946 sclL3 = l3Cur/l3NominalCurrent1;
949 if (TMath::Abs((sclL3=l3Cur/l3NominalCurrent1)-1.) < tolerance) map = AliMagF::k5kG;
950 else if (TMath::Abs((sclL3=l3Cur/l3NominalCurrent2)-1.) < tolerance) map = AliMagF::k2kG;
951 else if (l3Cur <= zero) { sclL3 = 0; map = AliMagF::k5kGUniform;}
953 AliError(Form("Wrong L3 current (%f A)!",l3Cur));
958 if (sclDip!=0 && (map==AliMagF::k5kG || map==AliMagF::k2kG) &&
959 ((convention==AliMagF::kConvLHC && l3Pol!=diPol) ||
960 (convention==AliMagF::kConvDCS2008 && l3Pol==diPol)) ) {
961 AliError(Form("Wrong combination for L3/Dipole polarities (%c/%c) for convention %d",
962 l3Pol>0?'+':'-',diPol>0?'+':'-',AliMagF::GetPolarityConvention()));
966 if (l3Pol<0) sclL3 = -sclL3;
967 if (diPol<0) sclDip = -sclDip;
969 AliMagF::BeamType_t btype = AliMagF::kNoBeamField;
970 TString btypestr = beamtype;
972 TPRegexp protonBeam("(proton|p)\\s*-?\\s*\\1");
973 TPRegexp ionBeam("(lead|pb|ion|a)\\s*-?\\s*\\1");
974 if (btypestr.Contains(ionBeam)) btype = AliMagF::kBeamTypeAA;
975 else if (btypestr.Contains(protonBeam)) btype = AliMagF::kBeamTypepp;
977 AliInfo(Form("Cannot determine the beam type from %s, assume no LHC magnet field",beamtype));
980 sprintf(ttl,"L3: %+5d Dip: %+4d kA; %s",(int)TMath::Sign(l3Cur,float(sclL3)),
981 (int)TMath::Sign(diCur,float(sclDip)),uniform ? " Constant":"");
982 AliMagF* fld = new AliMagF("MagneticFieldMap", ttl, 2, sclL3, sclDip, 10., map, path,
984 TGeoGlobalMagField::Instance()->SetField( fld );
985 TGeoGlobalMagField::Instance()->Lock();
990 Bool_t AliReconstruction::InitGRP() {
991 //------------------------------------
992 // Initialization of the GRP entry
993 //------------------------------------
994 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
998 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
1001 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
1003 fGRPData = new AliGRPObject();
1004 fGRPData->ReadValuesFromMap(m);
1008 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
1009 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
1013 // FIX ME: The unloading of GRP entry is temporarily disabled
1014 // because ZDC and VZERO are using it in order to initialize
1015 // their reconstructor objects. In the future one has to think
1016 // of propagating AliRunInfo to the reconstructors.
1017 // AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
1021 AliError("No GRP entry found in OCDB!");
1025 TString lhcState = fGRPData->GetLHCState();
1026 if (lhcState==AliGRPObject::GetInvalidString()) {
1027 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
1028 lhcState = "UNKNOWN";
1031 TString beamType = fGRPData->GetBeamType();
1032 if (beamType==AliGRPObject::GetInvalidString()) {
1033 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
1034 beamType = "UNKNOWN";
1037 Float_t beamEnergy = fGRPData->GetBeamEnergy();
1038 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
1039 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
1042 // energy is provided in MeV*120
1043 beamEnergy /= 120E3;
1045 TString runType = fGRPData->GetRunType();
1046 if (runType==AliGRPObject::GetInvalidString()) {
1047 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
1048 runType = "UNKNOWN";
1051 Int_t activeDetectors = fGRPData->GetDetectorMask();
1052 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
1053 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
1054 activeDetectors = 1074790399;
1057 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1061 // Process the list of active detectors
1062 if (activeDetectors) {
1063 UInt_t detMask = activeDetectors;
1064 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1065 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1066 fFillESD = MatchDetectorList(fFillESD,detMask);
1067 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1068 fLoadCDB.Form("%s %s %s %s",
1069 fRunLocalReconstruction.Data(),
1070 fRunTracking.Data(),
1072 fQADetectors.Data());
1073 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1074 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1)) {
1075 // switch off the vertexer
1076 AliInfo("SPD is not in the list of active detectors. Vertexer switched off.");
1077 fRunVertexFinder = kFALSE;
1079 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1080 // switch off the reading of CTP raw-data payload
1081 if (fFillTriggerESD) {
1082 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1083 fFillTriggerESD = kFALSE;
1088 AliInfo("===================================================================================");
1089 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1090 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1091 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1092 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1093 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1094 AliInfo("===================================================================================");
1096 //*** Dealing with the magnetic field map
1097 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {AliInfo("Running with the externally locked B field !");}
1099 // Construct the field map out of the information retrieved from GRP.
1102 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1103 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1104 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1108 Char_t l3Polarity = fGRPData->GetL3Polarity();
1109 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1110 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1115 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1116 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1117 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1121 Char_t diPolarity = fGRPData->GetDipolePolarity();
1122 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1123 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1127 // read special bits for the polarity convention and map type
1128 Int_t polConvention = fGRPData->IsPolarityConventionLHC() ? AliMagF::kConvLHC : AliMagF::kConvDCS2008;
1129 Bool_t uniformB = fGRPData->IsUniformBMap();
1132 if ( !SetFieldMap(l3Current, diCurrent, l3Polarity ? -1:1, diPolarity ? -1:1,
1133 polConvention,uniformB,beamEnergy, beamType.Data()))
1134 AliFatal("Failed to creat a B field map ! Exiting...");
1135 AliInfo("Running with the B field constructed out of GRP !");
1137 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1140 //*** Get the diamond profiles from OCDB
1141 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1143 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1145 AliError("No SPD diamond profile found in OCDB!");
1148 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1150 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1152 AliError("No diamond profile found in OCDB!");
1155 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1157 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1159 AliError("No TPC diamond profile found in OCDB!");
1162 entry = AliCDBManager::Instance()->Get("GRP/Calib/CosmicTriggers");
1164 fListOfCosmicTriggers = dynamic_cast<THashTable*>(entry->GetObject());
1166 AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
1169 if (!fListOfCosmicTriggers) {
1170 AliWarning("Can not get list of cosmic triggers from OCDB! Cosmic event specie will be effectively disabled!");
1176 //_____________________________________________________________________________
1177 Bool_t AliReconstruction::LoadCDB()
1179 AliCodeTimerAuto("");
1181 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1183 TString detStr = fLoadCDB;
1184 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1185 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1186 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1190 //_____________________________________________________________________________
1191 Bool_t AliReconstruction::LoadTriggerScalersCDB()
1193 AliCodeTimerAuto("");
1195 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
1199 AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
1200 fRunScalers = dynamic_cast<AliTriggerRunScalers*> (entry->GetObject());
1202 if (fRunScalers->CorrectScalersOverflow() == 0) AliInfo("32bit Trigger counters corrected for overflow");
1207 //_____________________________________________________________________________
1208 Bool_t AliReconstruction::Run(const char* input)
1211 AliCodeTimerAuto("");
1214 if (GetAbort() != TSelector::kContinue) return kFALSE;
1216 TChain *chain = NULL;
1217 if (fRawReader && (chain = fRawReader->GetChain())) {
1220 gProof->AddInput(this);
1221 if (!fESDOutput.IsValid()) {
1222 fESDOutput.SetProtocol("root",kTRUE);
1223 fESDOutput.SetHost(gSystem->HostName());
1224 fESDOutput.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1226 AliInfo(Form("Output file with ESDs is %s",fESDOutput.GetUrl()));
1227 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",fESDOutput.GetUrl()));
1228 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1230 chain->Process("AliReconstruction");
1233 chain->Process(this);
1238 if (GetAbort() != TSelector::kContinue) return kFALSE;
1240 if (GetAbort() != TSelector::kContinue) return kFALSE;
1241 //******* The loop over events
1242 AliInfo("Starting looping over events");
1244 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1245 (fRawReader && fRawReader->NextEvent())) {
1246 if (!ProcessEvent(iEvent)) {
1247 Abort("ProcessEvent",TSelector::kAbortFile);
1253 if (GetAbort() != TSelector::kContinue) return kFALSE;
1255 if (GetAbort() != TSelector::kContinue) return kFALSE;
1261 //_____________________________________________________________________________
1262 void AliReconstruction::InitRawReader(const char* input)
1264 AliCodeTimerAuto("");
1266 // Init raw-reader and
1267 // set the input in case of raw data
1268 if (input) fRawInput = input;
1269 fRawReader = AliRawReader::Create(fRawInput.Data());
1271 AliInfo("Reconstruction will run over digits");
1273 if (!fEquipIdMap.IsNull() && fRawReader)
1274 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1276 if (!fUseHLTData.IsNull()) {
1277 // create the RawReaderHLT which performs redirection of HLT input data for
1278 // the specified detectors
1279 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1281 fParentRawReader=fRawReader;
1282 fRawReader=pRawReader;
1284 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1287 AliSysInfo::AddStamp("CreateRawReader");
1290 //_____________________________________________________________________________
1291 void AliReconstruction::InitRun(const char* input)
1293 // Initialization of raw-reader,
1294 // run number, CDB etc.
1295 AliCodeTimerAuto("");
1296 AliSysInfo::AddStamp("Start");
1298 // Initialize raw-reader if any
1299 InitRawReader(input);
1301 // Initialize the CDB storage
1304 // Set run number in CDBManager (if it is not already set by the user)
1305 if (!SetRunNumberFromData()) {
1306 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1310 // Set CDB lock: from now on it is forbidden to reset the run number
1311 // or the default storage or to activate any further storage!
1316 //_____________________________________________________________________________
1317 void AliReconstruction::Begin(TTree *)
1319 // Initialize AlReconstruction before
1320 // going into the event loop
1321 // Should follow the TSelector convention
1322 // i.e. initialize only the object on the client side
1323 AliCodeTimerAuto("");
1325 AliReconstruction *reco = NULL;
1327 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1330 AliSysInfo::AddStamp("ReadInputInBegin");
1333 // Import ideal TGeo geometry and apply misalignment
1335 TString geom(gSystem->DirName(fGAliceFileName));
1336 geom += "/geometry.root";
1337 AliGeomManager::LoadGeometry(geom.Data());
1339 Abort("LoadGeometry", TSelector::kAbortProcess);
1342 AliSysInfo::AddStamp("LoadGeom");
1343 TString detsToCheck=fRunLocalReconstruction;
1344 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1345 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1348 AliSysInfo::AddStamp("CheckGeom");
1351 if (!MisalignGeometry(fLoadAlignData)) {
1352 Abort("MisalignGeometry", TSelector::kAbortProcess);
1355 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1356 AliSysInfo::AddStamp("MisalignGeom");
1359 Abort("InitGRP", TSelector::kAbortProcess);
1362 AliSysInfo::AddStamp("InitGRP");
1365 Abort("LoadCDB", TSelector::kAbortProcess);
1368 AliSysInfo::AddStamp("LoadCDB");
1370 if (!LoadTriggerScalersCDB()) {
1371 Abort("LoadTriggerScalersCDB", TSelector::kAbortProcess);
1374 AliSysInfo::AddStamp("LoadTriggerScalersCDB");
1377 // Read the reconstruction parameters from OCDB
1378 if (!InitRecoParams()) {
1379 AliWarning("Not all detectors have correct RecoParam objects initialized");
1381 AliSysInfo::AddStamp("InitRecoParams");
1383 if (fInput && gProof) {
1384 if (reco) *reco = *this;
1386 gProof->AddInputData(gGeoManager,kTRUE);
1388 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1389 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1390 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1391 magFieldMap->SetName("MagneticFieldMap");
1392 gProof->AddInputData(magFieldMap,kTRUE);
1397 //_____________________________________________________________________________
1398 void AliReconstruction::SlaveBegin(TTree*)
1400 // Initialization related to run-loader,
1401 // vertexer, trackers, recontructors
1402 // In proof mode it is executed on the slave
1403 AliCodeTimerAuto("");
1405 TProofOutputFile *outProofFile = NULL;
1407 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1410 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1412 AliGeomManager::SetGeometry(tgeo);
1414 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1415 Int_t runNumber = -1;
1416 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1417 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1418 man->SetCacheFlag(kTRUE);
1419 man->SetLock(kTRUE);
1423 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1424 TGeoGlobalMagField::Instance()->SetField(map);
1426 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1427 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1428 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1429 fOutput->Add(outProofFile);
1431 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1434 // get the run loader
1435 if (!InitRunLoader()) {
1436 Abort("InitRunLoader", TSelector::kAbortProcess);
1439 AliSysInfo::AddStamp("LoadLoader");
1441 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1444 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1445 Abort("CreateTrackers", TSelector::kAbortProcess);
1448 AliSysInfo::AddStamp("CreateTrackers");
1450 // create the ESD output file and tree
1451 if (!outProofFile) {
1452 ffile = TFile::Open("AliESDs.root", "RECREATE");
1453 ffile->SetCompressionLevel(2);
1454 if (!ffile->IsOpen()) {
1455 Abort("OpenESDFile", TSelector::kAbortProcess);
1460 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1461 Abort(Form("Problems opening output PROOF file: %s/%s",
1462 outProofFile->GetDir(), outProofFile->GetFileName()),
1463 TSelector::kAbortProcess);
1468 ftree = new TTree("esdTree", "Tree with ESD objects");
1469 fesd = new AliESDEvent();
1470 fesd->CreateStdContent();
1472 fesd->WriteToTree(ftree);
1473 if (fWriteESDfriend) {
1475 // Since we add the branch manually we must
1476 // book and add it after WriteToTree
1477 // otherwise it is created twice,
1478 // once via writetotree and once here.
1479 // The case for AliESDfriend is now
1480 // caught also in AlIESDEvent::WriteToTree but
1481 // be careful when changing the name (AliESDfriend is not
1482 // a TNamed so we had to hardwire it)
1483 fesdf = new AliESDfriend();
1484 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1485 br->SetFile("AliESDfriends.root");
1486 fesd->AddObject(fesdf);
1488 ftree->GetUserInfo()->Add(fesd);
1490 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1491 fhltesd = new AliESDEvent();
1492 fhltesd->CreateStdContent();
1494 // read the ESD template from CDB
1495 // HLT is allowed to put non-std content to its ESD, the non-std
1496 // objects need to be created before invocation of WriteToTree in
1497 // order to create all branches. Initialization is done from an
1498 // ESD layout template in CDB
1499 AliCDBManager* man = AliCDBManager::Instance();
1500 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1501 AliCDBEntry* hltESDConfig=NULL;
1502 if (man->GetId(hltESDConfigPath)!=NULL &&
1503 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1504 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1506 // init all internal variables from the list of objects
1507 pESDLayout->GetStdContent();
1509 // copy content and create non-std objects
1510 *fhltesd=*pESDLayout;
1513 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1514 hltESDConfigPath.GetPath().Data()));
1518 fhltesd->WriteToTree(fhlttree);
1519 fhlttree->GetUserInfo()->Add(fhltesd);
1521 ProcInfo_t procInfo;
1522 gSystem->GetProcInfo(&procInfo);
1523 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1526 //Initialize the QA and start of cycle
1527 if (fRunQA || fRunGlobalQA)
1530 //Initialize the Plane Efficiency framework
1531 if (fRunPlaneEff && !InitPlaneEff()) {
1532 Abort("InitPlaneEff", TSelector::kAbortProcess);
1536 if (strcmp(gProgName,"alieve") == 0)
1537 fRunAliEVE = InitAliEVE();
1542 //_____________________________________________________________________________
1543 Bool_t AliReconstruction::Process(Long64_t entry)
1545 // run the reconstruction over a single entry
1546 // from the chain with raw data
1547 AliCodeTimerAuto("");
1549 TTree *currTree = fChain->GetTree();
1550 AliRawVEvent *event = NULL;
1551 currTree->SetBranchAddress("rawevent",&event);
1552 currTree->GetEntry(entry);
1553 fRawReader = new AliRawReaderRoot(event);
1554 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1562 //_____________________________________________________________________________
1563 void AliReconstruction::Init(TTree *tree)
1566 AliError("The input tree is not found!");
1572 //_____________________________________________________________________________
1573 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1575 // run the reconstruction over a single event
1576 // The event loop is steered in Run method
1578 AliCodeTimerAuto("");
1580 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1581 fRunLoader->SetEventNumber(iEvent);
1582 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1584 fRunLoader->TreeE()->Fill();
1585 if (fRawReader && fRawReader->UseAutoSaveESD())
1586 fRunLoader->TreeE()->AutoSave("SaveSelf");
1589 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1593 AliInfo(Form("processing event %d", iEvent));
1595 fRunLoader->GetEvent(iEvent);
1597 // Fill Event-info object
1599 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo,fListOfCosmicTriggers);
1600 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1602 // Set the reco-params
1604 TString detStr = fLoadCDB;
1605 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1606 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1607 AliReconstructor *reconstructor = GetReconstructor(iDet);
1608 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1609 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1610 reconstructor->SetRecoParam(par);
1611 reconstructor->SetEventInfo(&fEventInfo);
1613 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1614 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(par->GetEventSpecie())) ;
1622 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1623 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1625 // local single event reconstruction
1626 if (!fRunLocalReconstruction.IsNull()) {
1627 TString detectors=fRunLocalReconstruction;
1628 // run HLT event reconstruction first
1629 // ;-( IsSelected changes the string
1630 if (IsSelected("HLT", detectors) &&
1631 !RunLocalEventReconstruction("HLT")) {
1632 if (fStopOnError) {CleanUp(); return kFALSE;}
1634 detectors=fRunLocalReconstruction;
1635 detectors.ReplaceAll("HLT", "");
1636 if (!RunLocalEventReconstruction(detectors)) {
1637 if (fStopOnError) {CleanUp(); return kFALSE;}
1641 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1642 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1643 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1644 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1646 // Set magnetic field from the tracker
1647 fesd->SetMagneticField(AliTracker::GetBz());
1648 fhltesd->SetMagneticField(AliTracker::GetBz());
1650 // Set most probable pt, for B=0 tracking
1651 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1652 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1653 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1655 // Fill raw-data error log into the ESD
1656 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1659 if (fRunVertexFinder) {
1660 if (!RunVertexFinder(fesd)) {
1661 if (fStopOnError) {CleanUp(); return kFALSE;}
1665 // For Plane Efficiency: run the SPD trackleter
1666 if (fRunPlaneEff && fSPDTrackleter) {
1667 if (!RunSPDTrackleting(fesd)) {
1668 if (fStopOnError) {CleanUp(); return kFALSE;}
1673 if (!fRunTracking.IsNull()) {
1674 if (fRunMuonTracking) {
1675 if (!RunMuonTracking(fesd)) {
1676 if (fStopOnError) {CleanUp(); return kFALSE;}
1682 if (!fRunTracking.IsNull()) {
1683 if (!RunTracking(fesd)) {
1684 if (fStopOnError) {CleanUp(); return kFALSE;}
1689 if (!fFillESD.IsNull()) {
1690 TString detectors=fFillESD;
1691 // run HLT first and on hltesd
1692 // ;-( IsSelected changes the string
1693 if (IsSelected("HLT", detectors) &&
1694 !FillESD(fhltesd, "HLT")) {
1695 if (fStopOnError) {CleanUp(); return kFALSE;}
1698 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1699 if (detectors.Contains("ALL")) {
1701 for (Int_t idet=0; idet<kNDetectors; ++idet){
1702 detectors += fgkDetectorName[idet];
1706 detectors.ReplaceAll("HLT", "");
1707 if (!FillESD(fesd, detectors)) {
1708 if (fStopOnError) {CleanUp(); return kFALSE;}
1712 // fill Event header information from the RawEventHeader
1713 if (fRawReader){FillRawEventHeaderESD(fesd);}
1714 if (fRawReader){FillRawEventHeaderESD(fhltesd);}
1717 AliESDpid::MakePID(fesd);
1719 if (fFillTriggerESD) {
1720 if (!FillTriggerESD(fesd)) {
1721 if (fStopOnError) {CleanUp(); return kFALSE;}
1728 // Propagate track to the beam pipe (if not already done by ITS)
1730 const Int_t ntracks = fesd->GetNumberOfTracks();
1731 const Double_t kBz = fesd->GetMagneticField();
1732 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1735 UShort_t *selectedIdx=new UShort_t[ntracks];
1737 for (Int_t itrack=0; itrack<ntracks; itrack++){
1738 const Double_t kMaxStep = 1; //max step over the material
1741 AliESDtrack *track = fesd->GetTrack(itrack);
1742 if (!track) continue;
1744 AliExternalTrackParam *tpcTrack =
1745 (AliExternalTrackParam *)track->GetTPCInnerParam();
1749 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1752 Int_t n=trkArray.GetEntriesFast();
1753 selectedIdx[n]=track->GetID();
1754 trkArray.AddLast(tpcTrack);
1757 //Tracks refitted by ITS should already be at the SPD vertex
1758 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1761 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1762 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1767 // Improve the reconstructed primary vertex position using the tracks
1769 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1770 if(fesd->GetPrimaryVertexSPD()) {
1771 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1772 if(vtitle.Contains("cosmics")) {
1773 runVertexFinderTracks=kFALSE;
1777 if (runVertexFinderTracks) {
1778 // TPC + ITS primary vertex
1779 ftVertexer->SetITSMode();
1780 ftVertexer->SetConstraintOff();
1781 // get cuts for vertexer from AliGRPRecoParam
1783 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1784 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1785 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1786 ftVertexer->SetCuts(cutsVertexer);
1787 delete [] cutsVertexer; cutsVertexer = NULL;
1788 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1789 ftVertexer->SetVtxStart(fDiamondProfile);
1791 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1793 if (pvtx->GetStatus()) {
1794 fesd->SetPrimaryVertexTracks(pvtx);
1795 for (Int_t i=0; i<ntracks; i++) {
1796 AliESDtrack *t = fesd->GetTrack(i);
1797 t->RelateToVertex(pvtx, kBz, kVeryBig);
1802 // TPC-only primary vertex
1803 ftVertexer->SetTPCMode();
1804 ftVertexer->SetConstraintOff();
1805 // get cuts for vertexer from AliGRPRecoParam
1807 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1808 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1809 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1810 ftVertexer->SetCuts(cutsVertexer);
1811 delete [] cutsVertexer; cutsVertexer = NULL;
1812 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1813 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1815 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1817 if (pvtx->GetStatus()) {
1818 fesd->SetPrimaryVertexTPC(pvtx);
1819 for (Int_t i=0; i<ntracks; i++) {
1820 AliESDtrack *t = fesd->GetTrack(i);
1821 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1827 delete[] selectedIdx;
1829 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1834 AliV0vertexer vtxer;
1835 vtxer.Tracks2V0vertices(fesd);
1837 if (fRunCascadeFinder) {
1839 AliCascadeVertexer cvtxer;
1840 cvtxer.V0sTracks2CascadeVertices(fesd);
1845 if (fCleanESD) CleanESD(fesd);
1848 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1849 AliQAManager::QAManager()->RunOneEvent(fesd) ;
1852 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1853 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1854 if (qadm && fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1855 qadm->Exec(AliQAv1::kESDS, fesd);
1858 if (fWriteESDfriend) {
1859 // fesdf->~AliESDfriend();
1860 // new (fesdf) AliESDfriend(); // Reset...
1861 fesd->GetESDfriend(fesdf);
1865 // Auto-save the ESD tree in case of prompt reco @P2
1866 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1867 ftree->AutoSave("SaveSelf");
1868 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1869 if (friendfile) friendfile->Save();
1876 if (fRunAliEVE) RunAliEVE();
1880 if (fWriteESDfriend) {
1881 fesdf->~AliESDfriend();
1882 new (fesdf) AliESDfriend(); // Reset...
1885 ProcInfo_t procInfo;
1886 gSystem->GetProcInfo(&procInfo);
1887 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1890 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1891 if (fReconstructor[iDet]) {
1892 fReconstructor[iDet]->SetRecoParam(NULL);
1893 fReconstructor[iDet]->SetEventInfo(NULL);
1895 if (fTracker[iDet]) fTracker[iDet]->SetEventInfo(NULL);
1898 if (fRunQA || fRunGlobalQA)
1899 AliQAManager::QAManager()->Increment() ;
1904 //_____________________________________________________________________________
1905 void AliReconstruction::SlaveTerminate()
1907 // Finalize the run on the slave side
1908 // Called after the exit
1909 // from the event loop
1910 AliCodeTimerAuto("");
1912 if (fIsNewRunLoader) { // galice.root didn't exist
1913 fRunLoader->WriteHeader("OVERWRITE");
1914 fRunLoader->CdGAFile();
1915 fRunLoader->Write(0, TObject::kOverwrite);
1918 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1919 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1921 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1922 cdbMapCopy->SetOwner(1);
1923 cdbMapCopy->SetName("cdbMap");
1924 TIter iter(cdbMap->GetTable());
1927 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1928 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1929 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1930 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1933 TList *cdbListCopy = new TList();
1934 cdbListCopy->SetOwner(1);
1935 cdbListCopy->SetName("cdbList");
1937 TIter iter2(cdbList);
1940 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1941 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1944 ftree->GetUserInfo()->Add(cdbMapCopy);
1945 ftree->GetUserInfo()->Add(cdbListCopy);
1950 if (fWriteESDfriend)
1951 ftree->SetBranchStatus("ESDfriend*",0);
1952 // we want to have only one tree version number
1953 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1954 fhlttree->Write(fhlttree->GetName(),TObject::kOverwrite);
1956 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1957 if (fRunPlaneEff && !FinishPlaneEff()) {
1958 AliWarning("Finish PlaneEff evaluation failed");
1961 // End of cycle for the in-loop
1963 AliQAManager::QAManager()->EndOfCycle() ;
1966 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1968 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS)))
1969 qadm->EndOfCycle(AliQAv1::kRECPOINTS);
1970 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1971 qadm->EndOfCycle(AliQAv1::kESDS);
1976 if (fRunQA || fRunGlobalQA) {
1978 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1979 TString qaOutputFile = outputFileName->GetTitle();
1980 qaOutputFile.ReplaceAll(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()),
1981 Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
1982 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
1983 qaProofFile->SetOutputFileName(qaOutputFile.Data());
1984 fOutput->Add(qaProofFile);
1985 MergeQA(qaProofFile->GetFileName());
1997 //_____________________________________________________________________________
1998 void AliReconstruction::Terminate()
2000 // Create tags for the events in the ESD tree (the ESD tree is always present)
2001 // In case of empty events the tags will contain dummy values
2002 AliCodeTimerAuto("");
2004 // Do not call the ESD tag creator in case of PROOF-based reconstruction
2006 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
2007 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
2010 // Cleanup of CDB manager: cache and active storages!
2011 AliCDBManager::Instance()->ClearCache();
2014 //_____________________________________________________________________________
2015 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
2017 // run the local reconstruction
2019 static Int_t eventNr=0;
2020 AliCodeTimerAuto("")
2022 TString detStr = detectors;
2023 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2024 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2025 AliReconstructor* reconstructor = GetReconstructor(iDet);
2026 if (!reconstructor) continue;
2027 AliLoader* loader = fLoader[iDet];
2028 // Matthias April 2008: temporary fix to run HLT reconstruction
2029 // although the HLT loader is missing
2030 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2032 reconstructor->Reconstruct(fRawReader, NULL);
2035 reconstructor->Reconstruct(dummy, NULL);
2040 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2043 // conversion of digits
2044 if (fRawReader && reconstructor->HasDigitConversion()) {
2045 AliInfo(Form("converting raw data digits into root objects for %s",
2046 fgkDetectorName[iDet]));
2047 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2048 // fgkDetectorName[iDet]));
2049 loader->LoadDigits("update");
2050 loader->CleanDigits();
2051 loader->MakeDigitsContainer();
2052 TTree* digitsTree = loader->TreeD();
2053 reconstructor->ConvertDigits(fRawReader, digitsTree);
2054 loader->WriteDigits("OVERWRITE");
2055 loader->UnloadDigits();
2057 // local reconstruction
2058 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2059 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2060 loader->LoadRecPoints("update");
2061 loader->CleanRecPoints();
2062 loader->MakeRecPointsContainer();
2063 TTree* clustersTree = loader->TreeR();
2064 if (fRawReader && !reconstructor->HasDigitConversion()) {
2065 reconstructor->Reconstruct(fRawReader, clustersTree);
2067 loader->LoadDigits("read");
2068 TTree* digitsTree = loader->TreeD();
2070 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2071 if (fStopOnError) return kFALSE;
2073 reconstructor->Reconstruct(digitsTree, clustersTree);
2075 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2076 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2079 loader->UnloadDigits();
2082 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2083 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2085 loader->WriteRecPoints("OVERWRITE");
2086 loader->UnloadRecPoints();
2087 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2089 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2090 AliError(Form("the following detectors were not found: %s",
2092 if (fStopOnError) return kFALSE;
2097 //_____________________________________________________________________________
2098 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2100 // run the SPD trackleting (for SPD efficiency purpouses)
2102 AliCodeTimerAuto("")
2104 Double_t vtxPos[3] = {0, 0, 0};
2105 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2107 TArrayF mcVertex(3);
2109 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2110 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2111 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2114 const AliESDVertex *vertex = esd->GetVertex();
2116 AliWarning("Vertex not found");
2119 vertex->GetXYZ(vtxPos);
2120 vertex->GetSigmaXYZ(vtxErr);
2121 if (fSPDTrackleter) {
2122 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2125 fLoader[0]->LoadRecPoints("read");
2126 TTree* tree = fLoader[0]->TreeR();
2128 AliError("Can't get the ITS cluster tree");
2131 fSPDTrackleter->LoadClusters(tree);
2132 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2134 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2135 AliError("AliITSTrackleterSPDEff Clusters2Tracks failed");
2136 // fLoader[0]->UnloadRecPoints();
2139 //fSPDTrackleter->UnloadRecPoints();
2141 AliWarning("SPDTrackleter not available");
2147 //_____________________________________________________________________________
2148 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2150 // run the barrel tracking
2152 AliCodeTimerAuto("")
2154 AliVertexer *vertexer = CreateVertexer();
2155 if (!vertexer) return kFALSE;
2157 AliInfo("running the ITS vertex finder");
2158 AliESDVertex* vertex = NULL;
2160 fLoader[0]->LoadRecPoints();
2161 TTree* cltree = fLoader[0]->TreeR();
2163 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2164 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2167 AliError("Can't get the ITS cluster tree");
2169 fLoader[0]->UnloadRecPoints();
2172 AliError("Can't get the ITS loader");
2175 AliWarning("Vertex not found");
2176 vertex = new AliESDVertex();
2177 vertex->SetName("default");
2180 vertex->SetName("reconstructed");
2185 vertex->GetXYZ(vtxPos);
2186 vertex->GetSigmaXYZ(vtxErr);
2188 esd->SetPrimaryVertexSPD(vertex);
2189 AliESDVertex *vpileup = NULL;
2190 Int_t novertices = 0;
2191 vpileup = vertexer->GetAllVertices(novertices);
2193 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2195 // if SPD multiplicity has been determined, it is stored in the ESD
2196 AliMultiplicity *mult = vertexer->GetMultiplicity();
2197 if(mult)esd->SetMultiplicity(mult);
2199 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2200 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2209 //_____________________________________________________________________________
2210 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2212 // run the HLT barrel tracking
2214 AliCodeTimerAuto("")
2217 AliError("Missing runLoader!");
2221 AliInfo("running HLT tracking");
2223 // Get a pointer to the HLT reconstructor
2224 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2225 if (!reconstructor) return kFALSE;
2228 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2229 TString detName = fgkDetectorName[iDet];
2230 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2231 reconstructor->SetOption(detName.Data());
2232 AliTracker *tracker = reconstructor->CreateTracker();
2234 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2235 if (fStopOnError) return kFALSE;
2239 Double_t vtxErr[3]={0.005,0.005,0.010};
2240 const AliESDVertex *vertex = esd->GetVertex();
2241 vertex->GetXYZ(vtxPos);
2242 tracker->SetVertex(vtxPos,vtxErr);
2244 fLoader[iDet]->LoadRecPoints("read");
2245 TTree* tree = fLoader[iDet]->TreeR();
2247 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2250 tracker->LoadClusters(tree);
2252 if (tracker->Clusters2Tracks(esd) != 0) {
2253 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2257 tracker->UnloadClusters();
2265 //_____________________________________________________________________________
2266 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2268 // run the muon spectrometer tracking
2270 AliCodeTimerAuto("")
2273 AliError("Missing runLoader!");
2276 Int_t iDet = 7; // for MUON
2278 AliInfo("is running...");
2280 // Get a pointer to the MUON reconstructor
2281 AliReconstructor *reconstructor = GetReconstructor(iDet);
2282 if (!reconstructor) return kFALSE;
2285 TString detName = fgkDetectorName[iDet];
2286 AliDebug(1, Form("%s tracking", detName.Data()));
2287 AliTracker *tracker = reconstructor->CreateTracker();
2289 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2294 fLoader[iDet]->LoadRecPoints("read");
2296 tracker->LoadClusters(fLoader[iDet]->TreeR());
2298 Int_t rv = tracker->Clusters2Tracks(esd);
2302 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2306 fLoader[iDet]->UnloadRecPoints();
2308 tracker->UnloadClusters();
2316 //_____________________________________________________________________________
2317 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2319 // run the barrel tracking
2320 static Int_t eventNr=0;
2321 AliCodeTimerAuto("")
2323 AliInfo("running tracking");
2325 // Set the event info which is used
2326 // by the trackers in order to obtain
2327 // information about read-out detectors,
2329 AliDebug(1, "Setting event info");
2330 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2331 if (!fTracker[iDet]) continue;
2332 fTracker[iDet]->SetEventInfo(&fEventInfo);
2335 //Fill the ESD with the T0 info (will be used by the TOF)
2336 if (fReconstructor[11] && fLoader[11]) {
2337 fLoader[11]->LoadRecPoints("READ");
2338 TTree *treeR = fLoader[11]->TreeR();
2340 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2344 // pass 1: TPC + ITS inwards
2345 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2346 if (!fTracker[iDet]) continue;
2347 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2350 fLoader[iDet]->LoadRecPoints("read");
2351 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2352 TTree* tree = fLoader[iDet]->TreeR();
2354 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2357 fTracker[iDet]->LoadClusters(tree);
2358 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2360 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2361 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2364 // preliminary PID in TPC needed by the ITS tracker
2366 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2367 AliESDpid::MakePID(esd);
2369 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2372 // pass 2: ALL backwards
2374 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2375 if (!fTracker[iDet]) continue;
2376 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2379 if (iDet > 1) { // all except ITS, TPC
2381 fLoader[iDet]->LoadRecPoints("read");
2382 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2383 tree = fLoader[iDet]->TreeR();
2385 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2388 fTracker[iDet]->LoadClusters(tree);
2389 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2393 if (iDet>1) // start filling residuals for the "outer" detectors
2395 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2396 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2398 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2399 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2400 if ( elem && (! elem->At(0)) ) {
2401 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2402 if (qadm) qadm->InitRecPointsForTracker() ;
2406 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2407 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2412 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2413 fTracker[iDet]->UnloadClusters();
2414 fLoader[iDet]->UnloadRecPoints();
2416 // updated PID in TPC needed by the ITS tracker -MI
2418 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2419 AliESDpid::MakePID(esd);
2421 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2423 //stop filling residuals for the "outer" detectors
2424 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2426 // pass 3: TRD + TPC + ITS refit inwards
2428 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2429 if (!fTracker[iDet]) continue;
2430 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2433 if (iDet<2) // start filling residuals for TPC and ITS
2435 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2436 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2438 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2439 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2440 if ( elem && (! elem->At(0)) ) {
2441 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2442 if (qadm) qadm->InitRecPointsForTracker() ;
2447 if (fTracker[iDet]->RefitInward(esd) != 0) {
2448 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2451 // run postprocessing
2452 if (fTracker[iDet]->PostProcess(esd) != 0) {
2453 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2456 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2459 // write space-points to the ESD in case alignment data output
2461 if (fWriteAlignmentData)
2462 WriteAlignmentData(esd);
2464 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2465 if (!fTracker[iDet]) continue;
2467 fTracker[iDet]->UnloadClusters();
2468 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2469 fLoader[iDet]->UnloadRecPoints();
2470 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2472 // stop filling residuals for TPC and ITS
2473 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2479 //_____________________________________________________________________________
2480 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2482 // Remove the data which are not needed for the physics analysis.
2485 Int_t nTracks=esd->GetNumberOfTracks();
2486 Int_t nV0s=esd->GetNumberOfV0s();
2488 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2490 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2491 Bool_t rc=esd->Clean(cleanPars);
2493 nTracks=esd->GetNumberOfTracks();
2494 nV0s=esd->GetNumberOfV0s();
2496 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2501 //_____________________________________________________________________________
2502 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2504 // fill the event summary data
2506 AliCodeTimerAuto("")
2507 static Int_t eventNr=0;
2508 TString detStr = detectors;
2510 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2511 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2512 AliReconstructor* reconstructor = GetReconstructor(iDet);
2513 if (!reconstructor) continue;
2514 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2515 TTree* clustersTree = NULL;
2516 if (fLoader[iDet]) {
2517 fLoader[iDet]->LoadRecPoints("read");
2518 clustersTree = fLoader[iDet]->TreeR();
2519 if (!clustersTree) {
2520 AliError(Form("Can't get the %s clusters tree",
2521 fgkDetectorName[iDet]));
2522 if (fStopOnError) return kFALSE;
2525 if (fRawReader && !reconstructor->HasDigitConversion()) {
2526 reconstructor->FillESD(fRawReader, clustersTree, esd);
2528 TTree* digitsTree = NULL;
2529 if (fLoader[iDet]) {
2530 fLoader[iDet]->LoadDigits("read");
2531 digitsTree = fLoader[iDet]->TreeD();
2533 AliError(Form("Can't get the %s digits tree",
2534 fgkDetectorName[iDet]));
2535 if (fStopOnError) return kFALSE;
2538 reconstructor->FillESD(digitsTree, clustersTree, esd);
2539 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2541 if (fLoader[iDet]) {
2542 fLoader[iDet]->UnloadRecPoints();
2546 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2547 AliError(Form("the following detectors were not found: %s",
2549 if (fStopOnError) return kFALSE;
2551 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2556 //_____________________________________________________________________________
2557 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2559 // Reads the trigger decision which is
2560 // stored in Trigger.root file and fills
2561 // the corresponding esd entries
2563 AliCodeTimerAuto("")
2565 AliInfo("Filling trigger information into the ESD");
2568 AliCTPRawStream input(fRawReader);
2569 if (!input.Next()) {
2570 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2573 if (esd->GetTriggerMask() != input.GetClassMask())
2574 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2575 input.GetClassMask(),esd->GetTriggerMask()));
2576 if (esd->GetOrbitNumber() != input.GetOrbitID())
2577 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2578 input.GetOrbitID(),esd->GetOrbitNumber()));
2579 if (esd->GetBunchCrossNumber() != input.GetBCID())
2580 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2581 input.GetBCID(),esd->GetBunchCrossNumber()));
2582 AliESDHeader* esdheader = esd->GetHeader();
2583 esdheader->SetL0TriggerInputs(input.GetL0Inputs());
2584 esdheader->SetL1TriggerInputs(input.GetL1Inputs());
2585 esdheader->SetL2TriggerInputs(input.GetL2Inputs());
2587 UInt_t orbit=input.GetOrbitID();
2588 for(Int_t i=0 ; i<input.GetNIRs() ; i++ )
2589 if(TMath::Abs(Int_t(orbit-(input.GetIR(i))->GetOrbit()))<=1){
2590 esdheader->AddTriggerIR(input.GetIR(i));
2595 //fRunScalers->Print();
2596 if(fRunScalers && fRunScalers->CheckRunScalers()){
2597 AliTimeStamp* timestamp = new AliTimeStamp(esd->GetOrbitNumber(), esd->GetPeriodNumber(), esd->GetBunchCrossNumber());
2598 //AliTimeStamp* timestamp = new AliTimeStamp(10308000, 0, (ULong64_t)486238);
2599 AliESDHeader* esdheader = fesd->GetHeader();
2600 for(Int_t i=0;i<50;i++){
2601 if((1<<i) & esd->GetTriggerMask()){
2602 AliTriggerScalersESD* scalesd = fRunScalers->GetScalersForEventClass( timestamp, i+1);
2603 if(scalesd)esdheader->SetTriggerScalersRecord(scalesd);
2609 //_____________________________________________________________________________
2610 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2613 // Filling information from RawReader Header
2616 if (!fRawReader) return kFALSE;
2618 AliInfo("Filling information from RawReader Header");
2620 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2621 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2622 esd->SetPeriodNumber(fRawReader->GetPeriod());
2624 esd->SetTimeStamp(fRawReader->GetTimestamp());
2625 esd->SetEventType(fRawReader->GetType());
2631 //_____________________________________________________________________________
2632 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2634 // check whether detName is contained in detectors
2635 // if yes, it is removed from detectors
2637 // check if all detectors are selected
2638 if ((detectors.CompareTo("ALL") == 0) ||
2639 detectors.BeginsWith("ALL ") ||
2640 detectors.EndsWith(" ALL") ||
2641 detectors.Contains(" ALL ")) {
2646 // search for the given detector
2647 Bool_t result = kFALSE;
2648 if ((detectors.CompareTo(detName) == 0) ||
2649 detectors.BeginsWith(detName+" ") ||
2650 detectors.EndsWith(" "+detName) ||
2651 detectors.Contains(" "+detName+" ")) {
2652 detectors.ReplaceAll(detName, "");
2656 // clean up the detectors string
2657 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2658 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2659 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2664 //_____________________________________________________________________________
2665 Bool_t AliReconstruction::InitRunLoader()
2667 // get or create the run loader
2669 if (gAlice) delete gAlice;
2672 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2673 // load all base libraries to get the loader classes
2674 TString libs = gSystem->GetLibraries();
2675 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2676 TString detName = fgkDetectorName[iDet];
2677 if (detName == "HLT") continue;
2678 if (libs.Contains("lib" + detName + "base.so")) continue;
2679 gSystem->Load("lib" + detName + "base.so");
2681 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2683 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2688 fRunLoader->CdGAFile();
2689 fRunLoader->LoadgAlice();
2691 //PH This is a temporary fix to give access to the kinematics
2692 //PH that is needed for the labels of ITS clusters
2693 fRunLoader->LoadHeader();
2694 fRunLoader->LoadKinematics();
2696 } else { // galice.root does not exist
2698 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2700 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2701 AliConfig::GetDefaultEventFolderName(),
2704 AliError(Form("could not create run loader in file %s",
2705 fGAliceFileName.Data()));
2709 fIsNewRunLoader = kTRUE;
2710 fRunLoader->MakeTree("E");
2712 if (fNumberOfEventsPerFile > 0)
2713 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2715 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2721 //_____________________________________________________________________________
2722 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2724 // get the reconstructor object and the loader for a detector
2726 if (fReconstructor[iDet]) {
2727 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2728 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2729 fReconstructor[iDet]->SetRecoParam(par);
2730 fReconstructor[iDet]->SetRunInfo(fRunInfo);
2732 return fReconstructor[iDet];
2735 // load the reconstructor object
2736 TPluginManager* pluginManager = gROOT->GetPluginManager();
2737 TString detName = fgkDetectorName[iDet];
2738 TString recName = "Ali" + detName + "Reconstructor";
2740 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2742 AliReconstructor* reconstructor = NULL;
2743 // first check if a plugin is defined for the reconstructor
2744 TPluginHandler* pluginHandler =
2745 pluginManager->FindHandler("AliReconstructor", detName);
2746 // if not, add a plugin for it
2747 if (!pluginHandler) {
2748 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2749 TString libs = gSystem->GetLibraries();
2750 if (libs.Contains("lib" + detName + "base.so") ||
2751 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2752 pluginManager->AddHandler("AliReconstructor", detName,
2753 recName, detName + "rec", recName + "()");
2755 pluginManager->AddHandler("AliReconstructor", detName,
2756 recName, detName, recName + "()");
2758 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2760 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2761 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2763 if (reconstructor) {
2764 TObject* obj = fOptions.FindObject(detName.Data());
2765 if (obj) reconstructor->SetOption(obj->GetTitle());
2766 reconstructor->Init();
2767 fReconstructor[iDet] = reconstructor;
2770 // get or create the loader
2771 if (detName != "HLT") {
2772 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2773 if (!fLoader[iDet]) {
2774 AliConfig::Instance()
2775 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2777 // first check if a plugin is defined for the loader
2779 pluginManager->FindHandler("AliLoader", detName);
2780 // if not, add a plugin for it
2781 if (!pluginHandler) {
2782 TString loaderName = "Ali" + detName + "Loader";
2783 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2784 pluginManager->AddHandler("AliLoader", detName,
2785 loaderName, detName + "base",
2786 loaderName + "(const char*, TFolder*)");
2787 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2789 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2791 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2792 fRunLoader->GetEventFolder());
2794 if (!fLoader[iDet]) { // use default loader
2795 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2797 if (!fLoader[iDet]) {
2798 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2799 if (fStopOnError) return NULL;
2801 fRunLoader->AddLoader(fLoader[iDet]);
2802 fRunLoader->CdGAFile();
2803 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2804 fRunLoader->Write(0, TObject::kOverwrite);
2809 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2810 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2811 reconstructor->SetRecoParam(par);
2812 reconstructor->SetRunInfo(fRunInfo);
2814 return reconstructor;
2817 //_____________________________________________________________________________
2818 AliVertexer* AliReconstruction::CreateVertexer()
2820 // create the vertexer
2821 // Please note that the caller is the owner of the
2824 AliVertexer* vertexer = NULL;
2825 AliReconstructor* itsReconstructor = GetReconstructor(0);
2826 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) || fRunTracking.Contains("ITS"))) {
2827 vertexer = itsReconstructor->CreateVertexer();
2830 AliWarning("couldn't create a vertexer for ITS");
2836 //_____________________________________________________________________________
2837 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2839 // create the trackers
2840 AliInfo("Creating trackers");
2842 TString detStr = detectors;
2843 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2844 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2845 AliReconstructor* reconstructor = GetReconstructor(iDet);
2846 if (!reconstructor) continue;
2847 TString detName = fgkDetectorName[iDet];
2848 if (detName == "HLT") {
2849 fRunHLTTracking = kTRUE;
2852 if (detName == "MUON") {
2853 fRunMuonTracking = kTRUE;
2858 fTracker[iDet] = reconstructor->CreateTracker();
2859 if (!fTracker[iDet] && (iDet < 7)) {
2860 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2861 if (fStopOnError) return kFALSE;
2863 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2869 //_____________________________________________________________________________
2870 void AliReconstruction::CleanUp()
2872 // delete trackers and the run loader and close and delete the file
2874 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2875 delete fReconstructor[iDet];
2876 fReconstructor[iDet] = NULL;
2877 fLoader[iDet] = NULL;
2878 delete fTracker[iDet];
2879 fTracker[iDet] = NULL;
2884 delete fSPDTrackleter;
2885 fSPDTrackleter = NULL;
2894 delete fParentRawReader;
2895 fParentRawReader=NULL;
2903 if (AliQAManager::QAManager())
2904 AliQAManager::QAManager()->ShowQA() ;
2905 AliQAManager::Destroy() ;
2907 TGeoGlobalMagField::Instance()->SetField(NULL);
2910 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2912 // Write space-points which are then used in the alignment procedures
2913 // For the moment only ITS, TPC, TRD and TOF
2915 Int_t ntracks = esd->GetNumberOfTracks();
2916 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2918 AliESDtrack *track = esd->GetTrack(itrack);
2921 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
2922 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2923 nsp += track->GetNcls(iDet);
2925 if (iDet==0) { // ITS "extra" clusters
2926 track->GetClusters(iDet,idx);
2927 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2932 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2933 track->SetTrackPointArray(sp);
2935 for (Int_t iDet = 5; iDet >= 0; iDet--) {
2936 AliTracker *tracker = fTracker[iDet];
2937 if (!tracker) continue;
2938 Int_t nspdet = track->GetClusters(iDet,idx);
2940 if (iDet==0) // ITS "extra" clusters
2941 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2943 if (nspdet <= 0) continue;
2947 while (isp2 < nspdet) {
2948 Bool_t isvalid=kTRUE;
2950 Int_t index=idx[isp++];
2951 if (index < 0) continue;
2953 TString dets = fgkDetectorName[iDet];
2954 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2955 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2956 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2957 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2958 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2960 isvalid = tracker->GetTrackPoint(index,p);
2963 if (!isvalid) continue;
2964 if (iDet==0 && (isp-1)>=6) p.SetExtra();
2965 sp->AddPoint(isptrack,&p); isptrack++;
2972 //_____________________________________________________________________________
2973 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2975 // The method reads the raw-data error log
2976 // accumulated within the rawReader.
2977 // It extracts the raw-data errors related to
2978 // the current event and stores them into
2979 // a TClonesArray inside the esd object.
2981 if (!fRawReader) return;
2983 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2985 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2987 if (iEvent != log->GetEventNumber()) continue;
2989 esd->AddRawDataErrorLog(log);
2994 //_____________________________________________________________________________
2995 void AliReconstruction::CheckQA()
2997 // check the QA of SIM for this run and remove the detectors
2998 // with status Fatal
3000 // TString newRunLocalReconstruction ;
3001 // TString newRunTracking ;
3002 // TString newFillESD ;
3004 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
3005 // TString detName(AliQAv1::GetDetName(iDet)) ;
3006 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
3007 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
3008 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
3009 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
3011 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
3012 // fRunLocalReconstruction.Contains("ALL") ) {
3013 // newRunLocalReconstruction += detName ;
3014 // newRunLocalReconstruction += " " ;
3016 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
3017 // fRunTracking.Contains("ALL") ) {
3018 // newRunTracking += detName ;
3019 // newRunTracking += " " ;
3021 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
3022 // fFillESD.Contains("ALL") ) {
3023 // newFillESD += detName ;
3024 // newFillESD += " " ;
3028 // fRunLocalReconstruction = newRunLocalReconstruction ;
3029 // fRunTracking = newRunTracking ;
3030 // fFillESD = newFillESD ;
3033 //_____________________________________________________________________________
3034 Int_t AliReconstruction::GetDetIndex(const char* detector)
3036 // return the detector index corresponding to detector
3038 for (index = 0; index < kNDetectors ; index++) {
3039 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3044 //_____________________________________________________________________________
3045 Bool_t AliReconstruction::FinishPlaneEff() {
3047 // Here execute all the necessary operationis, at the end of the tracking phase,
3048 // in case that evaluation of PlaneEfficiencies was required for some detector.
3049 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3051 // This Preliminary version works only FOR ITS !!!!!
3052 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3055 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3058 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3059 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3060 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3061 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
3062 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3063 TString name=planeeff->GetName();
3065 TFile* pefile = TFile::Open(name, "RECREATE");
3066 ret=(Bool_t)planeeff->Write();
3068 if(planeeff->GetCreateHistos()) {
3069 TString hname=planeeff->GetName();
3070 hname+="Histo.root";
3071 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3074 if(fSPDTrackleter) {
3075 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3076 TString name="AliITSPlaneEffSPDtracklet.root";
3077 TFile* pefile = TFile::Open(name, "RECREATE");
3078 ret=(Bool_t)planeeff->Write();
3080 AliESDEvent *dummy=NULL;
3081 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3086 //_____________________________________________________________________________
3087 Bool_t AliReconstruction::InitPlaneEff() {
3089 // Here execute all the necessary operations, before of the tracking phase,
3090 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3091 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3092 // which should be updated/recalculated.
3094 // This Preliminary version will work only FOR ITS !!!!!
3095 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3098 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3100 AliWarning(Form("Implementation of this method not yet completed !! Method return kTRUE"));
3102 fSPDTrackleter = NULL;
3103 AliReconstructor* itsReconstructor = GetReconstructor(0);
3104 if (itsReconstructor) {
3105 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3107 if (fSPDTrackleter) {
3108 AliInfo("Trackleter for SPD has been created");
3114 //_____________________________________________________________________________
3115 Bool_t AliReconstruction::InitAliEVE()
3117 // This method should be called only in case
3118 // AliReconstruction is run
3119 // within the alieve environment.
3120 // It will initialize AliEVE in a way
3121 // so that it can visualize event processed
3122 // by AliReconstruction.
3123 // The return flag shows whenever the
3124 // AliEVE initialization was successful or not.
3127 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3128 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3129 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3131 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3132 gROOT->ProcessLine("alieve_online_init()");
3137 //_____________________________________________________________________________
3138 void AliReconstruction::RunAliEVE()
3140 // Runs AliEVE visualisation of
3141 // the current event.
3142 // Should be executed only after
3143 // successful initialization of AliEVE.
3145 AliInfo("Running AliEVE...");
3146 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3150 //_____________________________________________________________________________
3151 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3153 // Allows to run QA for a selected set of detectors
3154 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3155 // all selected detectors run the same selected tasks
3157 if (!detAndAction.Contains(":")) {
3158 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3162 Int_t colon = detAndAction.Index(":") ;
3163 fQADetectors = detAndAction(0, colon) ;
3164 if (fQADetectors.Contains("ALL") )
3165 fQADetectors = fFillESD ;
3166 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3167 if (fQATasks.Contains("ALL") ) {
3168 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3170 fQATasks.ToUpper() ;
3172 if ( fQATasks.Contains("RAW") )
3173 tempo = Form("%d ", AliQAv1::kRAWS) ;
3174 if ( fQATasks.Contains("DIGIT") )
3175 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3176 if ( fQATasks.Contains("RECPOINT") )
3177 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3178 if ( fQATasks.Contains("ESD") )
3179 tempo += Form("%d ", AliQAv1::kESDS) ;
3181 if (fQATasks.IsNull()) {
3182 AliInfo("No QA requested\n") ;
3187 TString tempo(fQATasks) ;
3188 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3189 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3190 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3191 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3192 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3197 //_____________________________________________________________________________
3198 Bool_t AliReconstruction::InitRecoParams()
3200 // The method accesses OCDB and retrieves all
3201 // the available reco-param objects from there.
3203 Bool_t isOK = kTRUE;
3205 if (fRecoParam.GetDetRecoParamArray(kNDetectors)) {
3206 AliInfo("Using custom GRP reconstruction parameters");
3209 AliInfo("Loading GRP reconstruction parameter objects");
3211 AliCDBPath path("GRP","Calib","RecoParam");
3212 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3214 AliWarning("Couldn't find GRP RecoParam entry in OCDB");
3218 TObject *recoParamObj = entry->GetObject();
3219 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3220 // GRP has a normal TobjArray of AliDetectorRecoParam objects
3221 // Registering them in AliRecoParam
3222 fRecoParam.AddDetRecoParamArray(kNDetectors,dynamic_cast<TObjArray*>(recoParamObj));
3224 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3225 // GRP has only onse set of reco parameters
3226 // Registering it in AliRecoParam
3227 AliInfo("Single set of GRP reconstruction parameters found");
3228 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3229 fRecoParam.AddDetRecoParam(kNDetectors,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3232 AliError("No valid GRP RecoParam object found in the OCDB");
3239 TString detStr = fLoadCDB;
3240 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3242 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3244 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3245 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3249 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3251 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3252 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3254 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3258 TObject *recoParamObj = entry->GetObject();
3259 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3260 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3261 // Registering them in AliRecoParam
3262 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3264 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3265 // The detector has only onse set of reco parameters
3266 // Registering it in AliRecoParam
3267 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3268 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3269 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3272 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3276 // FIX ME: We have to disable the unloading of reco-param CDB
3277 // entries because QA framework is using them. Has to be fix in
3278 // a way that the QA takes the objects already constructed in
3280 // AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3284 if (AliDebugLevel() > 0) fRecoParam.Print();
3289 //_____________________________________________________________________________
3290 Bool_t AliReconstruction::GetEventInfo()
3292 // Fill the event info object
3294 AliCodeTimerAuto("")
3296 AliCentralTrigger *aCTP = NULL;
3298 fEventInfo.SetEventType(fRawReader->GetType());
3300 ULong64_t mask = fRawReader->GetClassMask();
3301 fEventInfo.SetTriggerMask(mask);
3302 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3303 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3305 aCTP = new AliCentralTrigger();
3306 TString configstr("");
3307 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3308 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3312 aCTP->SetClassMask(mask);
3313 aCTP->SetClusterMask(clmask);
3316 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3318 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3319 aCTP = fRunLoader->GetTrigger();
3320 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3321 // get inputs from actp - just get
3322 AliESDHeader* esdheader = fesd->GetHeader();
3323 esdheader->SetL0TriggerInputs(aCTP->GetL0TriggerInputs());
3324 esdheader->SetL1TriggerInputs(aCTP->GetL1TriggerInputs());
3325 esdheader->SetL2TriggerInputs(aCTP->GetL2TriggerInputs());
3326 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3329 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3334 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3336 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3337 if (fRawReader) delete aCTP;
3341 UChar_t clustmask = 0;
3343 ULong64_t trmask = fEventInfo.GetTriggerMask();
3344 const TObjArray& classesArray = config->GetClasses();
3345 Int_t nclasses = classesArray.GetEntriesFast();
3346 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3347 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3349 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3350 fesd->SetTriggerClass(trclass->GetName(),trindex);
3351 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3352 if (trmask & (1ull << trindex)) {
3354 trclasses += trclass->GetName();
3356 clustmask |= trclass->GetCluster()->GetClusterMask();
3360 fEventInfo.SetTriggerClasses(trclasses);
3362 // Set the information in ESD
3363 fesd->SetTriggerMask(trmask);
3364 fesd->SetTriggerCluster(clustmask);
3366 if (!aCTP->CheckTriggeredDetectors()) {
3367 if (fRawReader) delete aCTP;
3371 if (fRawReader) delete aCTP;
3373 // We have to fill also the HLT decision here!!
3379 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3381 // Match the detector list found in the rec.C or the default 'ALL'
3382 // to the list found in the GRP (stored there by the shuttle PP which
3383 // gets the information from ECS)
3384 static TString resultList;
3385 TString detList = detectorList;
3389 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3390 if ((detectorMask >> iDet) & 0x1) {
3391 TString det = AliDAQ::OfflineModuleName(iDet);
3392 if ((detList.CompareTo("ALL") == 0) ||
3393 ((detList.BeginsWith("ALL ") ||
3394 detList.EndsWith(" ALL") ||
3395 detList.Contains(" ALL ")) &&
3396 !(detList.BeginsWith("-"+det+" ") ||
3397 detList.EndsWith(" -"+det) ||
3398 detList.Contains(" -"+det+" "))) ||
3399 (detList.CompareTo(det) == 0) ||
3400 detList.BeginsWith(det+" ") ||
3401 detList.EndsWith(" "+det) ||
3402 detList.Contains( " "+det+" " )) {
3403 if (!resultList.EndsWith(det + " ")) {
3412 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3413 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3414 if ((detList.CompareTo("ALL") == 0) ||
3415 ((detList.BeginsWith("ALL ") ||
3416 detList.EndsWith(" ALL") ||
3417 detList.Contains(" ALL ")) &&
3418 !(detList.BeginsWith("-"+hltDet+" ") ||
3419 detList.EndsWith(" -"+hltDet) ||
3420 detList.Contains(" -"+hltDet+" "))) ||
3421 (detList.CompareTo(hltDet) == 0) ||
3422 detList.BeginsWith(hltDet+" ") ||
3423 detList.EndsWith(" "+hltDet) ||
3424 detList.Contains( " "+hltDet+" " )) {
3425 resultList += hltDet;
3429 return resultList.Data();
3433 //______________________________________________________________________________
3434 void AliReconstruction::Abort(const char *method, EAbort what)
3436 // Abort processing. If what = kAbortProcess, the Process() loop will be
3437 // aborted. If what = kAbortFile, the current file in a chain will be
3438 // aborted and the processing will continue with the next file, if there
3439 // is no next file then Process() will be aborted. Abort() can also be
3440 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3441 // the SlaveTerminate() and Terminate() are always called. The abort flag
3442 // can be checked in these methods using GetAbort().
3444 // The method is overwritten in AliReconstruction for better handling of
3445 // reco specific errors
3447 if (!fStopOnError) return;
3451 TString whyMess = method;
3452 whyMess += " failed! Aborting...";
3454 AliError(whyMess.Data());
3457 TString mess = "Abort";
3458 if (fAbort == kAbortProcess)
3459 mess = "AbortProcess";
3460 else if (fAbort == kAbortFile)
3463 Info(mess, whyMess.Data());