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) in %s convention",
962 l3Pol>0?'+':'-',diPol>0?'+':'-',convention==AliMagF::kConvDCS2008?"DCS2008":"LHC"));
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;
976 else AliInfo(Form("Assume no LHC magnet field for the beam type %s, ",beamtype));
979 sprintf(ttl,"L3: %+5d Dip: %+4d kA; %s | Polarities in %s convention",(int)TMath::Sign(l3Cur,float(sclL3)),
980 (int)TMath::Sign(diCur,float(sclDip)),uniform ? " Constant":"",
981 convention==AliMagF::kConvLHC ? "LHC":"DCS2008");
982 // LHC and DCS08 conventions have opposite dipole polarities
983 if ( AliMagF::GetPolarityConvention() != convention) sclDip = -sclDip;
985 AliMagF* fld = new AliMagF("MagneticFieldMap", ttl, 2, sclL3, sclDip, 10., map, path,
987 TGeoGlobalMagField::Instance()->SetField( fld );
988 TGeoGlobalMagField::Instance()->Lock();
993 Bool_t AliReconstruction::InitGRP() {
994 //------------------------------------
995 // Initialization of the GRP entry
996 //------------------------------------
997 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
1001 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
1004 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
1006 fGRPData = new AliGRPObject();
1007 fGRPData->ReadValuesFromMap(m);
1011 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
1012 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
1016 // FIX ME: The unloading of GRP entry is temporarily disabled
1017 // because ZDC and VZERO are using it in order to initialize
1018 // their reconstructor objects. In the future one has to think
1019 // of propagating AliRunInfo to the reconstructors.
1020 // AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
1024 AliError("No GRP entry found in OCDB!");
1028 TString lhcState = fGRPData->GetLHCState();
1029 if (lhcState==AliGRPObject::GetInvalidString()) {
1030 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
1031 lhcState = "UNKNOWN";
1034 TString beamType = fGRPData->GetBeamType();
1035 if (beamType==AliGRPObject::GetInvalidString()) {
1036 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
1037 beamType = "UNKNOWN";
1040 Float_t beamEnergy = fGRPData->GetBeamEnergy();
1041 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
1042 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
1045 // LHC: "multiply by 120 to get the energy in MeV"
1046 beamEnergy *= 0.120;
1048 TString runType = fGRPData->GetRunType();
1049 if (runType==AliGRPObject::GetInvalidString()) {
1050 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
1051 runType = "UNKNOWN";
1054 Int_t activeDetectors = fGRPData->GetDetectorMask();
1055 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
1056 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
1057 activeDetectors = 1074790399;
1060 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1064 // Process the list of active detectors
1065 if (activeDetectors) {
1066 UInt_t detMask = activeDetectors;
1067 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1068 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1069 fFillESD = MatchDetectorList(fFillESD,detMask);
1070 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1071 fLoadCDB.Form("%s %s %s %s",
1072 fRunLocalReconstruction.Data(),
1073 fRunTracking.Data(),
1075 fQADetectors.Data());
1076 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1077 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1)) {
1078 // switch off the vertexer
1079 AliInfo("SPD is not in the list of active detectors. Vertexer switched off.");
1080 fRunVertexFinder = kFALSE;
1082 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1083 // switch off the reading of CTP raw-data payload
1084 if (fFillTriggerESD) {
1085 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1086 fFillTriggerESD = kFALSE;
1091 AliInfo("===================================================================================");
1092 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1093 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1094 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1095 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1096 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1097 AliInfo("===================================================================================");
1099 //*** Dealing with the magnetic field map
1100 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {
1101 if (TGeoGlobalMagField::Instance()->GetField()->TestBit(AliMagF::kOverrideGRP)) {
1102 AliInfo("ExpertMode!!! GRP information will be ignored !");
1103 AliInfo("ExpertMode!!! Running with the externally locked B field !");
1106 AliInfo("Destroying existing B field instance!");
1107 delete TGeoGlobalMagField::Instance();
1110 if ( !TGeoGlobalMagField::Instance()->IsLocked() ) {
1111 // Construct the field map out of the information retrieved from GRP.
1114 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1115 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1116 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1120 Char_t l3Polarity = fGRPData->GetL3Polarity();
1121 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1122 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1127 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1128 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1129 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1133 Char_t diPolarity = fGRPData->GetDipolePolarity();
1134 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1135 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1139 // read special bits for the polarity convention and map type
1140 Int_t polConvention = fGRPData->IsPolarityConventionLHC() ? AliMagF::kConvLHC : AliMagF::kConvDCS2008;
1141 Bool_t uniformB = fGRPData->IsUniformBMap();
1144 if ( !SetFieldMap(l3Current, diCurrent, l3Polarity ? -1:1, diPolarity ? -1:1,
1145 polConvention,uniformB,beamEnergy, beamType.Data()))
1146 AliFatal("Failed to creat a B field map ! Exiting...");
1147 AliInfo("Running with the B field constructed out of GRP !");
1149 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1152 //*** Get the diamond profiles from OCDB
1153 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1155 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1157 AliError("No SPD diamond profile found in OCDB!");
1160 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1162 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1164 AliError("No diamond profile found in OCDB!");
1167 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1169 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1171 AliError("No TPC diamond profile found in OCDB!");
1174 entry = AliCDBManager::Instance()->Get("GRP/Calib/CosmicTriggers");
1176 fListOfCosmicTriggers = dynamic_cast<THashTable*>(entry->GetObject());
1178 AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
1181 if (!fListOfCosmicTriggers) {
1182 AliWarning("Can not get list of cosmic triggers from OCDB! Cosmic event specie will be effectively disabled!");
1188 //_____________________________________________________________________________
1189 Bool_t AliReconstruction::LoadCDB()
1191 AliCodeTimerAuto("");
1193 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1195 TString detStr = fLoadCDB;
1196 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1197 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1198 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1202 //_____________________________________________________________________________
1203 Bool_t AliReconstruction::LoadTriggerScalersCDB()
1205 AliCodeTimerAuto("");
1207 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
1211 AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
1212 fRunScalers = dynamic_cast<AliTriggerRunScalers*> (entry->GetObject());
1214 if (fRunScalers->CorrectScalersOverflow() == 0) AliInfo("32bit Trigger counters corrected for overflow");
1219 //_____________________________________________________________________________
1220 Bool_t AliReconstruction::Run(const char* input)
1223 AliCodeTimerAuto("");
1226 if (GetAbort() != TSelector::kContinue) return kFALSE;
1228 TChain *chain = NULL;
1229 if (fRawReader && (chain = fRawReader->GetChain())) {
1232 gProof->AddInput(this);
1233 if (!fESDOutput.IsValid()) {
1234 fESDOutput.SetProtocol("root",kTRUE);
1235 fESDOutput.SetHost(gSystem->HostName());
1236 fESDOutput.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1238 AliInfo(Form("Output file with ESDs is %s",fESDOutput.GetUrl()));
1239 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",fESDOutput.GetUrl()));
1240 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1242 chain->Process("AliReconstruction");
1245 chain->Process(this);
1250 if (GetAbort() != TSelector::kContinue) return kFALSE;
1252 if (GetAbort() != TSelector::kContinue) return kFALSE;
1253 //******* The loop over events
1254 AliInfo("Starting looping over events");
1256 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1257 (fRawReader && fRawReader->NextEvent())) {
1258 if (!ProcessEvent(iEvent)) {
1259 Abort("ProcessEvent",TSelector::kAbortFile);
1265 if (GetAbort() != TSelector::kContinue) return kFALSE;
1267 if (GetAbort() != TSelector::kContinue) return kFALSE;
1273 //_____________________________________________________________________________
1274 void AliReconstruction::InitRawReader(const char* input)
1276 AliCodeTimerAuto("");
1278 // Init raw-reader and
1279 // set the input in case of raw data
1280 if (input) fRawInput = input;
1281 fRawReader = AliRawReader::Create(fRawInput.Data());
1283 AliInfo("Reconstruction will run over digits");
1285 if (!fEquipIdMap.IsNull() && fRawReader)
1286 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1288 if (!fUseHLTData.IsNull()) {
1289 // create the RawReaderHLT which performs redirection of HLT input data for
1290 // the specified detectors
1291 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1293 fParentRawReader=fRawReader;
1294 fRawReader=pRawReader;
1296 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1299 AliSysInfo::AddStamp("CreateRawReader");
1302 //_____________________________________________________________________________
1303 void AliReconstruction::InitRun(const char* input)
1305 // Initialization of raw-reader,
1306 // run number, CDB etc.
1307 AliCodeTimerAuto("");
1308 AliSysInfo::AddStamp("Start");
1310 // Initialize raw-reader if any
1311 InitRawReader(input);
1313 // Initialize the CDB storage
1316 // Set run number in CDBManager (if it is not already set by the user)
1317 if (!SetRunNumberFromData()) {
1318 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1322 // Set CDB lock: from now on it is forbidden to reset the run number
1323 // or the default storage or to activate any further storage!
1328 //_____________________________________________________________________________
1329 void AliReconstruction::Begin(TTree *)
1331 // Initialize AlReconstruction before
1332 // going into the event loop
1333 // Should follow the TSelector convention
1334 // i.e. initialize only the object on the client side
1335 AliCodeTimerAuto("");
1337 AliReconstruction *reco = NULL;
1339 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1342 AliSysInfo::AddStamp("ReadInputInBegin");
1345 // Import ideal TGeo geometry and apply misalignment
1347 TString geom(gSystem->DirName(fGAliceFileName));
1348 geom += "/geometry.root";
1349 AliGeomManager::LoadGeometry(geom.Data());
1351 Abort("LoadGeometry", TSelector::kAbortProcess);
1354 AliSysInfo::AddStamp("LoadGeom");
1355 TString detsToCheck=fRunLocalReconstruction;
1356 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1357 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1360 AliSysInfo::AddStamp("CheckGeom");
1363 if (!MisalignGeometry(fLoadAlignData)) {
1364 Abort("MisalignGeometry", TSelector::kAbortProcess);
1367 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1368 AliSysInfo::AddStamp("MisalignGeom");
1371 Abort("InitGRP", TSelector::kAbortProcess);
1374 AliSysInfo::AddStamp("InitGRP");
1377 Abort("LoadCDB", TSelector::kAbortProcess);
1380 AliSysInfo::AddStamp("LoadCDB");
1382 if (!LoadTriggerScalersCDB()) {
1383 Abort("LoadTriggerScalersCDB", TSelector::kAbortProcess);
1386 AliSysInfo::AddStamp("LoadTriggerScalersCDB");
1389 // Read the reconstruction parameters from OCDB
1390 if (!InitRecoParams()) {
1391 AliWarning("Not all detectors have correct RecoParam objects initialized");
1393 AliSysInfo::AddStamp("InitRecoParams");
1395 if (fInput && gProof) {
1396 if (reco) *reco = *this;
1398 gProof->AddInputData(gGeoManager,kTRUE);
1400 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1401 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1402 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1403 magFieldMap->SetName("MagneticFieldMap");
1404 gProof->AddInputData(magFieldMap,kTRUE);
1409 //_____________________________________________________________________________
1410 void AliReconstruction::SlaveBegin(TTree*)
1412 // Initialization related to run-loader,
1413 // vertexer, trackers, recontructors
1414 // In proof mode it is executed on the slave
1415 AliCodeTimerAuto("");
1417 TProofOutputFile *outProofFile = NULL;
1419 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1422 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1424 AliGeomManager::SetGeometry(tgeo);
1426 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1427 Int_t runNumber = -1;
1428 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1429 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1430 man->SetCacheFlag(kTRUE);
1431 man->SetLock(kTRUE);
1435 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1436 TGeoGlobalMagField::Instance()->SetField(map);
1438 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1439 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1440 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1441 fOutput->Add(outProofFile);
1443 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1446 // get the run loader
1447 if (!InitRunLoader()) {
1448 Abort("InitRunLoader", TSelector::kAbortProcess);
1451 AliSysInfo::AddStamp("LoadLoader");
1453 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1456 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1457 Abort("CreateTrackers", TSelector::kAbortProcess);
1460 AliSysInfo::AddStamp("CreateTrackers");
1462 // create the ESD output file and tree
1463 if (!outProofFile) {
1464 ffile = TFile::Open("AliESDs.root", "RECREATE");
1465 ffile->SetCompressionLevel(2);
1466 if (!ffile->IsOpen()) {
1467 Abort("OpenESDFile", TSelector::kAbortProcess);
1472 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1473 Abort(Form("Problems opening output PROOF file: %s/%s",
1474 outProofFile->GetDir(), outProofFile->GetFileName()),
1475 TSelector::kAbortProcess);
1480 ftree = new TTree("esdTree", "Tree with ESD objects");
1481 fesd = new AliESDEvent();
1482 fesd->CreateStdContent();
1484 fesd->WriteToTree(ftree);
1485 if (fWriteESDfriend) {
1487 // Since we add the branch manually we must
1488 // book and add it after WriteToTree
1489 // otherwise it is created twice,
1490 // once via writetotree and once here.
1491 // The case for AliESDfriend is now
1492 // caught also in AlIESDEvent::WriteToTree but
1493 // be careful when changing the name (AliESDfriend is not
1494 // a TNamed so we had to hardwire it)
1495 fesdf = new AliESDfriend();
1496 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1497 br->SetFile("AliESDfriends.root");
1498 fesd->AddObject(fesdf);
1500 ftree->GetUserInfo()->Add(fesd);
1502 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1503 fhltesd = new AliESDEvent();
1504 fhltesd->CreateStdContent();
1506 // read the ESD template from CDB
1507 // HLT is allowed to put non-std content to its ESD, the non-std
1508 // objects need to be created before invocation of WriteToTree in
1509 // order to create all branches. Initialization is done from an
1510 // ESD layout template in CDB
1511 AliCDBManager* man = AliCDBManager::Instance();
1512 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1513 AliCDBEntry* hltESDConfig=NULL;
1514 if (man->GetId(hltESDConfigPath)!=NULL &&
1515 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1516 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1518 // init all internal variables from the list of objects
1519 pESDLayout->GetStdContent();
1521 // copy content and create non-std objects
1522 *fhltesd=*pESDLayout;
1525 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1526 hltESDConfigPath.GetPath().Data()));
1530 fhltesd->WriteToTree(fhlttree);
1531 fhlttree->GetUserInfo()->Add(fhltesd);
1533 ProcInfo_t procInfo;
1534 gSystem->GetProcInfo(&procInfo);
1535 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1538 //Initialize the QA and start of cycle
1539 if (fRunQA || fRunGlobalQA)
1542 //Initialize the Plane Efficiency framework
1543 if (fRunPlaneEff && !InitPlaneEff()) {
1544 Abort("InitPlaneEff", TSelector::kAbortProcess);
1548 if (strcmp(gProgName,"alieve") == 0)
1549 fRunAliEVE = InitAliEVE();
1554 //_____________________________________________________________________________
1555 Bool_t AliReconstruction::Process(Long64_t entry)
1557 // run the reconstruction over a single entry
1558 // from the chain with raw data
1559 AliCodeTimerAuto("");
1561 TTree *currTree = fChain->GetTree();
1562 AliRawVEvent *event = NULL;
1563 currTree->SetBranchAddress("rawevent",&event);
1564 currTree->GetEntry(entry);
1565 fRawReader = new AliRawReaderRoot(event);
1566 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1574 //_____________________________________________________________________________
1575 void AliReconstruction::Init(TTree *tree)
1578 AliError("The input tree is not found!");
1584 //_____________________________________________________________________________
1585 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1587 // run the reconstruction over a single event
1588 // The event loop is steered in Run method
1590 AliCodeTimerAuto("");
1592 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1593 fRunLoader->SetEventNumber(iEvent);
1594 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1596 fRunLoader->TreeE()->Fill();
1597 if (fRawReader && fRawReader->UseAutoSaveESD())
1598 fRunLoader->TreeE()->AutoSave("SaveSelf");
1601 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1605 AliInfo(Form("processing event %d", iEvent));
1607 fRunLoader->GetEvent(iEvent);
1609 // Fill Event-info object
1611 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo,fListOfCosmicTriggers);
1612 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1614 // Set the reco-params
1616 TString detStr = fLoadCDB;
1617 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1618 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1619 AliReconstructor *reconstructor = GetReconstructor(iDet);
1620 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1621 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1622 reconstructor->SetRecoParam(par);
1623 reconstructor->SetEventInfo(&fEventInfo);
1625 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1626 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(par->GetEventSpecie())) ;
1634 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1635 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1637 // local single event reconstruction
1638 if (!fRunLocalReconstruction.IsNull()) {
1639 TString detectors=fRunLocalReconstruction;
1640 // run HLT event reconstruction first
1641 // ;-( IsSelected changes the string
1642 if (IsSelected("HLT", detectors) &&
1643 !RunLocalEventReconstruction("HLT")) {
1644 if (fStopOnError) {CleanUp(); return kFALSE;}
1646 detectors=fRunLocalReconstruction;
1647 detectors.ReplaceAll("HLT", "");
1648 if (!RunLocalEventReconstruction(detectors)) {
1649 if (fStopOnError) {CleanUp(); return kFALSE;}
1653 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1654 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1655 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1656 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1658 // Set magnetic field from the tracker
1659 fesd->SetMagneticField(AliTracker::GetBz());
1660 fhltesd->SetMagneticField(AliTracker::GetBz());
1662 // Set most probable pt, for B=0 tracking
1663 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1664 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1665 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1667 // Fill raw-data error log into the ESD
1668 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1671 if (fRunVertexFinder) {
1672 if (!RunVertexFinder(fesd)) {
1673 if (fStopOnError) {CleanUp(); return kFALSE;}
1677 // For Plane Efficiency: run the SPD trackleter
1678 if (fRunPlaneEff && fSPDTrackleter) {
1679 if (!RunSPDTrackleting(fesd)) {
1680 if (fStopOnError) {CleanUp(); return kFALSE;}
1685 if (!fRunTracking.IsNull()) {
1686 if (fRunMuonTracking) {
1687 if (!RunMuonTracking(fesd)) {
1688 if (fStopOnError) {CleanUp(); return kFALSE;}
1694 if (!fRunTracking.IsNull()) {
1695 if (!RunTracking(fesd)) {
1696 if (fStopOnError) {CleanUp(); return kFALSE;}
1701 if (!fFillESD.IsNull()) {
1702 TString detectors=fFillESD;
1703 // run HLT first and on hltesd
1704 // ;-( IsSelected changes the string
1705 if (IsSelected("HLT", detectors) &&
1706 !FillESD(fhltesd, "HLT")) {
1707 if (fStopOnError) {CleanUp(); return kFALSE;}
1710 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1711 if (detectors.Contains("ALL")) {
1713 for (Int_t idet=0; idet<kNDetectors; ++idet){
1714 detectors += fgkDetectorName[idet];
1718 detectors.ReplaceAll("HLT", "");
1719 if (!FillESD(fesd, detectors)) {
1720 if (fStopOnError) {CleanUp(); return kFALSE;}
1724 // fill Event header information from the RawEventHeader
1725 if (fRawReader){FillRawEventHeaderESD(fesd);}
1726 if (fRawReader){FillRawEventHeaderESD(fhltesd);}
1729 AliESDpid::MakePID(fesd);
1731 if (fFillTriggerESD) {
1732 if (!FillTriggerESD(fesd)) {
1733 if (fStopOnError) {CleanUp(); return kFALSE;}
1736 // Always fill scalers
1737 if (!FillTriggerScalers(fesd)) {
1738 if (fStopOnError) {CleanUp(); return kFALSE;}
1745 // Propagate track to the beam pipe (if not already done by ITS)
1747 const Int_t ntracks = fesd->GetNumberOfTracks();
1748 const Double_t kBz = fesd->GetMagneticField();
1749 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1752 UShort_t *selectedIdx=new UShort_t[ntracks];
1754 for (Int_t itrack=0; itrack<ntracks; itrack++){
1755 const Double_t kMaxStep = 1; //max step over the material
1758 AliESDtrack *track = fesd->GetTrack(itrack);
1759 if (!track) continue;
1761 AliExternalTrackParam *tpcTrack =
1762 (AliExternalTrackParam *)track->GetTPCInnerParam();
1766 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1769 Int_t n=trkArray.GetEntriesFast();
1770 selectedIdx[n]=track->GetID();
1771 trkArray.AddLast(tpcTrack);
1774 //Tracks refitted by ITS should already be at the SPD vertex
1775 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1778 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1779 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1784 // Improve the reconstructed primary vertex position using the tracks
1786 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1787 if(fesd->GetPrimaryVertexSPD()) {
1788 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1789 if(vtitle.Contains("cosmics")) {
1790 runVertexFinderTracks=kFALSE;
1794 if (runVertexFinderTracks) {
1795 // TPC + ITS primary vertex
1796 ftVertexer->SetITSMode();
1797 ftVertexer->SetConstraintOff();
1798 // get cuts for vertexer from AliGRPRecoParam
1800 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1801 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1802 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1803 ftVertexer->SetCuts(cutsVertexer);
1804 delete [] cutsVertexer; cutsVertexer = NULL;
1805 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1806 ftVertexer->SetVtxStart(fDiamondProfile);
1808 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1810 if (pvtx->GetStatus()) {
1811 fesd->SetPrimaryVertexTracks(pvtx);
1812 for (Int_t i=0; i<ntracks; i++) {
1813 AliESDtrack *t = fesd->GetTrack(i);
1814 t->RelateToVertex(pvtx, kBz, kVeryBig);
1819 // TPC-only primary vertex
1820 ftVertexer->SetTPCMode();
1821 ftVertexer->SetConstraintOff();
1822 // get cuts for vertexer from AliGRPRecoParam
1824 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1825 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1826 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1827 ftVertexer->SetCuts(cutsVertexer);
1828 delete [] cutsVertexer; cutsVertexer = NULL;
1829 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1830 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1832 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1834 if (pvtx->GetStatus()) {
1835 fesd->SetPrimaryVertexTPC(pvtx);
1836 for (Int_t i=0; i<ntracks; i++) {
1837 AliESDtrack *t = fesd->GetTrack(i);
1838 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1844 delete[] selectedIdx;
1846 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1851 AliV0vertexer vtxer;
1852 vtxer.Tracks2V0vertices(fesd);
1854 if (fRunCascadeFinder) {
1856 AliCascadeVertexer cvtxer;
1857 cvtxer.V0sTracks2CascadeVertices(fesd);
1862 if (fCleanESD) CleanESD(fesd);
1865 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1866 AliQAManager::QAManager()->RunOneEvent(fesd) ;
1869 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1870 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1871 if (qadm && fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1872 qadm->Exec(AliQAv1::kESDS, fesd);
1875 if (fWriteESDfriend) {
1876 // fesdf->~AliESDfriend();
1877 // new (fesdf) AliESDfriend(); // Reset...
1878 fesd->GetESDfriend(fesdf);
1882 // Auto-save the ESD tree in case of prompt reco @P2
1883 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1884 ftree->AutoSave("SaveSelf");
1885 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1886 if (friendfile) friendfile->Save();
1893 if (fRunAliEVE) RunAliEVE();
1897 if (fWriteESDfriend) {
1898 fesdf->~AliESDfriend();
1899 new (fesdf) AliESDfriend(); // Reset...
1902 ProcInfo_t procInfo;
1903 gSystem->GetProcInfo(&procInfo);
1904 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1907 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1908 if (fReconstructor[iDet]) {
1909 fReconstructor[iDet]->SetRecoParam(NULL);
1910 fReconstructor[iDet]->SetEventInfo(NULL);
1912 if (fTracker[iDet]) fTracker[iDet]->SetEventInfo(NULL);
1915 if (fRunQA || fRunGlobalQA)
1916 AliQAManager::QAManager()->Increment() ;
1921 //_____________________________________________________________________________
1922 void AliReconstruction::SlaveTerminate()
1924 // Finalize the run on the slave side
1925 // Called after the exit
1926 // from the event loop
1927 AliCodeTimerAuto("");
1929 if (fIsNewRunLoader) { // galice.root didn't exist
1930 fRunLoader->WriteHeader("OVERWRITE");
1931 fRunLoader->CdGAFile();
1932 fRunLoader->Write(0, TObject::kOverwrite);
1935 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1936 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1938 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1939 cdbMapCopy->SetOwner(1);
1940 cdbMapCopy->SetName("cdbMap");
1941 TIter iter(cdbMap->GetTable());
1944 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1945 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1946 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1947 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1950 TList *cdbListCopy = new TList();
1951 cdbListCopy->SetOwner(1);
1952 cdbListCopy->SetName("cdbList");
1954 TIter iter2(cdbList);
1957 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1958 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1961 ftree->GetUserInfo()->Add(cdbMapCopy);
1962 ftree->GetUserInfo()->Add(cdbListCopy);
1967 if (fWriteESDfriend)
1968 ftree->SetBranchStatus("ESDfriend*",0);
1969 // we want to have only one tree version number
1970 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1971 fhlttree->Write(fhlttree->GetName(),TObject::kOverwrite);
1973 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1974 if (fRunPlaneEff && !FinishPlaneEff()) {
1975 AliWarning("Finish PlaneEff evaluation failed");
1978 // End of cycle for the in-loop
1980 AliQAManager::QAManager()->EndOfCycle() ;
1983 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1985 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS)))
1986 qadm->EndOfCycle(AliQAv1::kRECPOINTS);
1987 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1988 qadm->EndOfCycle(AliQAv1::kESDS);
1993 if (fRunQA || fRunGlobalQA) {
1995 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1996 TString qaOutputFile = outputFileName->GetTitle();
1997 qaOutputFile.ReplaceAll(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()),
1998 Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
1999 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
2000 qaProofFile->SetOutputFileName(qaOutputFile.Data());
2001 fOutput->Add(qaProofFile);
2002 MergeQA(qaProofFile->GetFileName());
2014 //_____________________________________________________________________________
2015 void AliReconstruction::Terminate()
2017 // Create tags for the events in the ESD tree (the ESD tree is always present)
2018 // In case of empty events the tags will contain dummy values
2019 AliCodeTimerAuto("");
2021 // Do not call the ESD tag creator in case of PROOF-based reconstruction
2023 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
2024 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
2027 // Cleanup of CDB manager: cache and active storages!
2028 AliCDBManager::Instance()->ClearCache();
2031 //_____________________________________________________________________________
2032 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
2034 // run the local reconstruction
2036 static Int_t eventNr=0;
2037 AliCodeTimerAuto("")
2039 TString detStr = detectors;
2040 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2041 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2042 AliReconstructor* reconstructor = GetReconstructor(iDet);
2043 if (!reconstructor) continue;
2044 AliLoader* loader = fLoader[iDet];
2045 // Matthias April 2008: temporary fix to run HLT reconstruction
2046 // although the HLT loader is missing
2047 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2049 reconstructor->Reconstruct(fRawReader, NULL);
2052 reconstructor->Reconstruct(dummy, NULL);
2057 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2060 // conversion of digits
2061 if (fRawReader && reconstructor->HasDigitConversion()) {
2062 AliInfo(Form("converting raw data digits into root objects for %s",
2063 fgkDetectorName[iDet]));
2064 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2065 // fgkDetectorName[iDet]));
2066 loader->LoadDigits("update");
2067 loader->CleanDigits();
2068 loader->MakeDigitsContainer();
2069 TTree* digitsTree = loader->TreeD();
2070 reconstructor->ConvertDigits(fRawReader, digitsTree);
2071 loader->WriteDigits("OVERWRITE");
2072 loader->UnloadDigits();
2074 // local reconstruction
2075 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2076 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2077 loader->LoadRecPoints("update");
2078 loader->CleanRecPoints();
2079 loader->MakeRecPointsContainer();
2080 TTree* clustersTree = loader->TreeR();
2081 if (fRawReader && !reconstructor->HasDigitConversion()) {
2082 reconstructor->Reconstruct(fRawReader, clustersTree);
2084 loader->LoadDigits("read");
2085 TTree* digitsTree = loader->TreeD();
2087 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2088 if (fStopOnError) return kFALSE;
2090 reconstructor->Reconstruct(digitsTree, clustersTree);
2092 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2093 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2096 loader->UnloadDigits();
2099 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2100 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2102 loader->WriteRecPoints("OVERWRITE");
2103 loader->UnloadRecPoints();
2104 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2106 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2107 AliError(Form("the following detectors were not found: %s",
2109 if (fStopOnError) return kFALSE;
2114 //_____________________________________________________________________________
2115 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2117 // run the SPD trackleting (for SPD efficiency purpouses)
2119 AliCodeTimerAuto("")
2121 Double_t vtxPos[3] = {0, 0, 0};
2122 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2124 TArrayF mcVertex(3);
2126 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2127 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2128 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2131 const AliESDVertex *vertex = esd->GetVertex();
2133 AliWarning("Vertex not found");
2136 vertex->GetXYZ(vtxPos);
2137 vertex->GetSigmaXYZ(vtxErr);
2138 if (fSPDTrackleter) {
2139 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2142 fLoader[0]->LoadRecPoints("read");
2143 TTree* tree = fLoader[0]->TreeR();
2145 AliError("Can't get the ITS cluster tree");
2148 fSPDTrackleter->LoadClusters(tree);
2149 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2151 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2152 AliError("AliITSTrackleterSPDEff Clusters2Tracks failed");
2153 // fLoader[0]->UnloadRecPoints();
2156 //fSPDTrackleter->UnloadRecPoints();
2158 AliWarning("SPDTrackleter not available");
2164 //_____________________________________________________________________________
2165 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2167 // run the barrel tracking
2169 AliCodeTimerAuto("")
2171 AliVertexer *vertexer = CreateVertexer();
2172 if (!vertexer) return kFALSE;
2174 AliInfo("running the ITS vertex finder");
2175 AliESDVertex* vertex = NULL;
2177 fLoader[0]->LoadRecPoints();
2178 TTree* cltree = fLoader[0]->TreeR();
2180 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2181 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2184 AliError("Can't get the ITS cluster tree");
2186 fLoader[0]->UnloadRecPoints();
2189 AliError("Can't get the ITS loader");
2192 AliWarning("Vertex not found");
2193 vertex = new AliESDVertex();
2194 vertex->SetName("default");
2197 vertex->SetName("reconstructed");
2202 vertex->GetXYZ(vtxPos);
2203 vertex->GetSigmaXYZ(vtxErr);
2205 esd->SetPrimaryVertexSPD(vertex);
2206 AliESDVertex *vpileup = NULL;
2207 Int_t novertices = 0;
2208 vpileup = vertexer->GetAllVertices(novertices);
2210 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2212 // if SPD multiplicity has been determined, it is stored in the ESD
2213 AliMultiplicity *mult = vertexer->GetMultiplicity();
2214 if(mult)esd->SetMultiplicity(mult);
2216 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2217 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2226 //_____________________________________________________________________________
2227 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2229 // run the HLT barrel tracking
2231 AliCodeTimerAuto("")
2234 AliError("Missing runLoader!");
2238 AliInfo("running HLT tracking");
2240 // Get a pointer to the HLT reconstructor
2241 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2242 if (!reconstructor) return kFALSE;
2245 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2246 TString detName = fgkDetectorName[iDet];
2247 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2248 reconstructor->SetOption(detName.Data());
2249 AliTracker *tracker = reconstructor->CreateTracker();
2251 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2252 if (fStopOnError) return kFALSE;
2256 Double_t vtxErr[3]={0.005,0.005,0.010};
2257 const AliESDVertex *vertex = esd->GetVertex();
2258 vertex->GetXYZ(vtxPos);
2259 tracker->SetVertex(vtxPos,vtxErr);
2261 fLoader[iDet]->LoadRecPoints("read");
2262 TTree* tree = fLoader[iDet]->TreeR();
2264 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2267 tracker->LoadClusters(tree);
2269 if (tracker->Clusters2Tracks(esd) != 0) {
2270 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2274 tracker->UnloadClusters();
2282 //_____________________________________________________________________________
2283 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2285 // run the muon spectrometer tracking
2287 AliCodeTimerAuto("")
2290 AliError("Missing runLoader!");
2293 Int_t iDet = 7; // for MUON
2295 AliInfo("is running...");
2297 // Get a pointer to the MUON reconstructor
2298 AliReconstructor *reconstructor = GetReconstructor(iDet);
2299 if (!reconstructor) return kFALSE;
2302 TString detName = fgkDetectorName[iDet];
2303 AliDebug(1, Form("%s tracking", detName.Data()));
2304 AliTracker *tracker = reconstructor->CreateTracker();
2306 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2311 fLoader[iDet]->LoadRecPoints("read");
2313 tracker->LoadClusters(fLoader[iDet]->TreeR());
2315 Int_t rv = tracker->Clusters2Tracks(esd);
2319 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2323 fLoader[iDet]->UnloadRecPoints();
2325 tracker->UnloadClusters();
2333 //_____________________________________________________________________________
2334 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2336 // run the barrel tracking
2337 static Int_t eventNr=0;
2338 AliCodeTimerAuto("")
2340 AliInfo("running tracking");
2342 // Set the event info which is used
2343 // by the trackers in order to obtain
2344 // information about read-out detectors,
2346 AliDebug(1, "Setting event info");
2347 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2348 if (!fTracker[iDet]) continue;
2349 fTracker[iDet]->SetEventInfo(&fEventInfo);
2352 //Fill the ESD with the T0 info (will be used by the TOF)
2353 if (fReconstructor[11] && fLoader[11]) {
2354 fLoader[11]->LoadRecPoints("READ");
2355 TTree *treeR = fLoader[11]->TreeR();
2357 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2361 // pass 1: TPC + ITS inwards
2362 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2363 if (!fTracker[iDet]) continue;
2364 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2367 fLoader[iDet]->LoadRecPoints("read");
2368 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2369 TTree* tree = fLoader[iDet]->TreeR();
2371 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2374 fTracker[iDet]->LoadClusters(tree);
2375 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2377 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2378 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2381 // preliminary PID in TPC needed by the ITS tracker
2383 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2384 AliESDpid::MakePID(esd);
2386 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2389 // pass 2: ALL backwards
2391 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2392 if (!fTracker[iDet]) continue;
2393 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2396 if (iDet > 1) { // all except ITS, TPC
2398 fLoader[iDet]->LoadRecPoints("read");
2399 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2400 tree = fLoader[iDet]->TreeR();
2402 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2405 fTracker[iDet]->LoadClusters(tree);
2406 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2410 if (iDet>1) // start filling residuals for the "outer" detectors
2412 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2413 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2415 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2416 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2417 if ( elem && (! elem->At(0)) ) {
2418 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2419 if (qadm) qadm->InitRecPointsForTracker() ;
2423 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2424 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2429 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2430 fTracker[iDet]->UnloadClusters();
2431 fLoader[iDet]->UnloadRecPoints();
2433 // updated PID in TPC needed by the ITS tracker -MI
2435 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2436 AliESDpid::MakePID(esd);
2438 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2440 //stop filling residuals for the "outer" detectors
2441 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2443 // pass 3: TRD + TPC + ITS refit inwards
2445 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2446 if (!fTracker[iDet]) continue;
2447 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2450 if (iDet<2) // start filling residuals for TPC and ITS
2452 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2453 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2455 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2456 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2457 if ( elem && (! elem->At(0)) ) {
2458 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2459 if (qadm) qadm->InitRecPointsForTracker() ;
2464 if (fTracker[iDet]->RefitInward(esd) != 0) {
2465 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2468 // run postprocessing
2469 if (fTracker[iDet]->PostProcess(esd) != 0) {
2470 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2473 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2476 // write space-points to the ESD in case alignment data output
2478 if (fWriteAlignmentData)
2479 WriteAlignmentData(esd);
2481 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2482 if (!fTracker[iDet]) continue;
2484 fTracker[iDet]->UnloadClusters();
2485 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2486 fLoader[iDet]->UnloadRecPoints();
2487 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2489 // stop filling residuals for TPC and ITS
2490 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2496 //_____________________________________________________________________________
2497 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2499 // Remove the data which are not needed for the physics analysis.
2502 Int_t nTracks=esd->GetNumberOfTracks();
2503 Int_t nV0s=esd->GetNumberOfV0s();
2505 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2507 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2508 Bool_t rc=esd->Clean(cleanPars);
2510 nTracks=esd->GetNumberOfTracks();
2511 nV0s=esd->GetNumberOfV0s();
2513 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2518 //_____________________________________________________________________________
2519 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2521 // fill the event summary data
2523 AliCodeTimerAuto("")
2524 static Int_t eventNr=0;
2525 TString detStr = detectors;
2527 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2528 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2529 AliReconstructor* reconstructor = GetReconstructor(iDet);
2530 if (!reconstructor) continue;
2531 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2532 TTree* clustersTree = NULL;
2533 if (fLoader[iDet]) {
2534 fLoader[iDet]->LoadRecPoints("read");
2535 clustersTree = fLoader[iDet]->TreeR();
2536 if (!clustersTree) {
2537 AliError(Form("Can't get the %s clusters tree",
2538 fgkDetectorName[iDet]));
2539 if (fStopOnError) return kFALSE;
2542 if (fRawReader && !reconstructor->HasDigitConversion()) {
2543 reconstructor->FillESD(fRawReader, clustersTree, esd);
2545 TTree* digitsTree = NULL;
2546 if (fLoader[iDet]) {
2547 fLoader[iDet]->LoadDigits("read");
2548 digitsTree = fLoader[iDet]->TreeD();
2550 AliError(Form("Can't get the %s digits tree",
2551 fgkDetectorName[iDet]));
2552 if (fStopOnError) return kFALSE;
2555 reconstructor->FillESD(digitsTree, clustersTree, esd);
2556 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2558 if (fLoader[iDet]) {
2559 fLoader[iDet]->UnloadRecPoints();
2563 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2564 AliError(Form("the following detectors were not found: %s",
2566 if (fStopOnError) return kFALSE;
2568 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2573 //_____________________________________________________________________________
2574 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2576 // Reads the trigger decision which is
2577 // stored in Trigger.root file and fills
2578 // the corresponding esd entries
2580 AliCodeTimerAuto("")
2582 AliInfo("Filling trigger information into the ESD");
2585 AliCTPRawStream input(fRawReader);
2586 if (!input.Next()) {
2587 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2590 if (esd->GetTriggerMask() != input.GetClassMask())
2591 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2592 input.GetClassMask(),esd->GetTriggerMask()));
2593 if (esd->GetOrbitNumber() != input.GetOrbitID())
2594 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2595 input.GetOrbitID(),esd->GetOrbitNumber()));
2596 if (esd->GetBunchCrossNumber() != input.GetBCID())
2597 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2598 input.GetBCID(),esd->GetBunchCrossNumber()));
2599 AliESDHeader* esdheader = esd->GetHeader();
2600 esdheader->SetL0TriggerInputs(input.GetL0Inputs());
2601 esdheader->SetL1TriggerInputs(input.GetL1Inputs());
2602 esdheader->SetL2TriggerInputs(input.GetL2Inputs());
2604 UInt_t orbit=input.GetOrbitID();
2605 for(Int_t i=0 ; i<input.GetNIRs() ; i++ )
2606 if(TMath::Abs(Int_t(orbit-(input.GetIR(i))->GetOrbit()))<=1){
2607 esdheader->AddTriggerIR(input.GetIR(i));
2613 //_____________________________________________________________________________
2614 Bool_t AliReconstruction::FillTriggerScalers(AliESDEvent*& esd)
2617 //fRunScalers->Print();
2618 if(fRunScalers && fRunScalers->CheckRunScalers()){
2619 AliTimeStamp* timestamp = new AliTimeStamp(esd->GetOrbitNumber(), esd->GetPeriodNumber(), esd->GetBunchCrossNumber());
2620 //AliTimeStamp* timestamp = new AliTimeStamp(10308000, 0, (ULong64_t)486238);
2621 AliESDHeader* esdheader = fesd->GetHeader();
2622 for(Int_t i=0;i<50;i++){
2623 if((1<<i) & esd->GetTriggerMask()){
2624 AliTriggerScalersESD* scalesd = fRunScalers->GetScalersForEventClass( timestamp, i+1);
2625 if(scalesd)esdheader->SetTriggerScalersRecord(scalesd);
2631 //_____________________________________________________________________________
2632 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2635 // Filling information from RawReader Header
2638 if (!fRawReader) return kFALSE;
2640 AliInfo("Filling information from RawReader Header");
2642 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2643 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2644 esd->SetPeriodNumber(fRawReader->GetPeriod());
2646 esd->SetTimeStamp(fRawReader->GetTimestamp());
2647 esd->SetEventType(fRawReader->GetType());
2653 //_____________________________________________________________________________
2654 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2656 // check whether detName is contained in detectors
2657 // if yes, it is removed from detectors
2659 // check if all detectors are selected
2660 if ((detectors.CompareTo("ALL") == 0) ||
2661 detectors.BeginsWith("ALL ") ||
2662 detectors.EndsWith(" ALL") ||
2663 detectors.Contains(" ALL ")) {
2668 // search for the given detector
2669 Bool_t result = kFALSE;
2670 if ((detectors.CompareTo(detName) == 0) ||
2671 detectors.BeginsWith(detName+" ") ||
2672 detectors.EndsWith(" "+detName) ||
2673 detectors.Contains(" "+detName+" ")) {
2674 detectors.ReplaceAll(detName, "");
2678 // clean up the detectors string
2679 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2680 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2681 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2686 //_____________________________________________________________________________
2687 Bool_t AliReconstruction::InitRunLoader()
2689 // get or create the run loader
2691 if (gAlice) delete gAlice;
2694 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2695 // load all base libraries to get the loader classes
2696 TString libs = gSystem->GetLibraries();
2697 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2698 TString detName = fgkDetectorName[iDet];
2699 if (detName == "HLT") continue;
2700 if (libs.Contains("lib" + detName + "base.so")) continue;
2701 gSystem->Load("lib" + detName + "base.so");
2703 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2705 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2710 fRunLoader->CdGAFile();
2711 fRunLoader->LoadgAlice();
2713 //PH This is a temporary fix to give access to the kinematics
2714 //PH that is needed for the labels of ITS clusters
2715 fRunLoader->LoadHeader();
2716 fRunLoader->LoadKinematics();
2718 } else { // galice.root does not exist
2720 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2722 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2723 AliConfig::GetDefaultEventFolderName(),
2726 AliError(Form("could not create run loader in file %s",
2727 fGAliceFileName.Data()));
2731 fIsNewRunLoader = kTRUE;
2732 fRunLoader->MakeTree("E");
2734 if (fNumberOfEventsPerFile > 0)
2735 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2737 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2743 //_____________________________________________________________________________
2744 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2746 // get the reconstructor object and the loader for a detector
2748 if (fReconstructor[iDet]) {
2749 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2750 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2751 fReconstructor[iDet]->SetRecoParam(par);
2752 fReconstructor[iDet]->SetRunInfo(fRunInfo);
2754 return fReconstructor[iDet];
2757 // load the reconstructor object
2758 TPluginManager* pluginManager = gROOT->GetPluginManager();
2759 TString detName = fgkDetectorName[iDet];
2760 TString recName = "Ali" + detName + "Reconstructor";
2762 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2764 AliReconstructor* reconstructor = NULL;
2765 // first check if a plugin is defined for the reconstructor
2766 TPluginHandler* pluginHandler =
2767 pluginManager->FindHandler("AliReconstructor", detName);
2768 // if not, add a plugin for it
2769 if (!pluginHandler) {
2770 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2771 TString libs = gSystem->GetLibraries();
2772 if (libs.Contains("lib" + detName + "base.so") ||
2773 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2774 pluginManager->AddHandler("AliReconstructor", detName,
2775 recName, detName + "rec", recName + "()");
2777 pluginManager->AddHandler("AliReconstructor", detName,
2778 recName, detName, recName + "()");
2780 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2782 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2783 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2785 if (reconstructor) {
2786 TObject* obj = fOptions.FindObject(detName.Data());
2787 if (obj) reconstructor->SetOption(obj->GetTitle());
2788 reconstructor->SetRunInfo(fRunInfo);
2789 reconstructor->Init();
2790 fReconstructor[iDet] = reconstructor;
2793 // get or create the loader
2794 if (detName != "HLT") {
2795 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2796 if (!fLoader[iDet]) {
2797 AliConfig::Instance()
2798 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2800 // first check if a plugin is defined for the loader
2802 pluginManager->FindHandler("AliLoader", detName);
2803 // if not, add a plugin for it
2804 if (!pluginHandler) {
2805 TString loaderName = "Ali" + detName + "Loader";
2806 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2807 pluginManager->AddHandler("AliLoader", detName,
2808 loaderName, detName + "base",
2809 loaderName + "(const char*, TFolder*)");
2810 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2812 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2814 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2815 fRunLoader->GetEventFolder());
2817 if (!fLoader[iDet]) { // use default loader
2818 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2820 if (!fLoader[iDet]) {
2821 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2822 if (fStopOnError) return NULL;
2824 fRunLoader->AddLoader(fLoader[iDet]);
2825 fRunLoader->CdGAFile();
2826 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2827 fRunLoader->Write(0, TObject::kOverwrite);
2832 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2833 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2834 reconstructor->SetRecoParam(par);
2835 reconstructor->SetRunInfo(fRunInfo);
2837 return reconstructor;
2840 //_____________________________________________________________________________
2841 AliVertexer* AliReconstruction::CreateVertexer()
2843 // create the vertexer
2844 // Please note that the caller is the owner of the
2847 AliVertexer* vertexer = NULL;
2848 AliReconstructor* itsReconstructor = GetReconstructor(0);
2849 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) || fRunTracking.Contains("ITS"))) {
2850 vertexer = itsReconstructor->CreateVertexer();
2853 AliWarning("couldn't create a vertexer for ITS");
2859 //_____________________________________________________________________________
2860 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2862 // create the trackers
2863 AliInfo("Creating trackers");
2865 TString detStr = detectors;
2866 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2867 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2868 AliReconstructor* reconstructor = GetReconstructor(iDet);
2869 if (!reconstructor) continue;
2870 TString detName = fgkDetectorName[iDet];
2871 if (detName == "HLT") {
2872 fRunHLTTracking = kTRUE;
2875 if (detName == "MUON") {
2876 fRunMuonTracking = kTRUE;
2881 fTracker[iDet] = reconstructor->CreateTracker();
2882 if (!fTracker[iDet] && (iDet < 7)) {
2883 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2884 if (fStopOnError) return kFALSE;
2886 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2892 //_____________________________________________________________________________
2893 void AliReconstruction::CleanUp()
2895 // delete trackers and the run loader and close and delete the file
2897 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2898 delete fReconstructor[iDet];
2899 fReconstructor[iDet] = NULL;
2900 fLoader[iDet] = NULL;
2901 delete fTracker[iDet];
2902 fTracker[iDet] = NULL;
2907 delete fSPDTrackleter;
2908 fSPDTrackleter = NULL;
2917 delete fParentRawReader;
2918 fParentRawReader=NULL;
2926 if (AliQAManager::QAManager())
2927 AliQAManager::QAManager()->ShowQA() ;
2928 AliQAManager::Destroy() ;
2930 TGeoGlobalMagField::Instance()->SetField(NULL);
2933 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2935 // Write space-points which are then used in the alignment procedures
2936 // For the moment only ITS, TPC, TRD and TOF
2938 Int_t ntracks = esd->GetNumberOfTracks();
2939 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2941 AliESDtrack *track = esd->GetTrack(itrack);
2944 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
2945 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2946 nsp += track->GetNcls(iDet);
2948 if (iDet==0) { // ITS "extra" clusters
2949 track->GetClusters(iDet,idx);
2950 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2955 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2956 track->SetTrackPointArray(sp);
2958 for (Int_t iDet = 5; iDet >= 0; iDet--) {
2959 AliTracker *tracker = fTracker[iDet];
2960 if (!tracker) continue;
2961 Int_t nspdet = track->GetClusters(iDet,idx);
2963 if (iDet==0) // ITS "extra" clusters
2964 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2966 if (nspdet <= 0) continue;
2970 while (isp2 < nspdet) {
2971 Bool_t isvalid=kTRUE;
2973 Int_t index=idx[isp++];
2974 if (index < 0) continue;
2976 TString dets = fgkDetectorName[iDet];
2977 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2978 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2979 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2980 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2981 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2983 isvalid = tracker->GetTrackPoint(index,p);
2986 if (!isvalid) continue;
2987 if (iDet==0 && (isp-1)>=6) p.SetExtra();
2988 sp->AddPoint(isptrack,&p); isptrack++;
2995 //_____________________________________________________________________________
2996 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2998 // The method reads the raw-data error log
2999 // accumulated within the rawReader.
3000 // It extracts the raw-data errors related to
3001 // the current event and stores them into
3002 // a TClonesArray inside the esd object.
3004 if (!fRawReader) return;
3006 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
3008 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
3010 if (iEvent != log->GetEventNumber()) continue;
3012 esd->AddRawDataErrorLog(log);
3017 //_____________________________________________________________________________
3018 void AliReconstruction::CheckQA()
3020 // check the QA of SIM for this run and remove the detectors
3021 // with status Fatal
3023 // TString newRunLocalReconstruction ;
3024 // TString newRunTracking ;
3025 // TString newFillESD ;
3027 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
3028 // TString detName(AliQAv1::GetDetName(iDet)) ;
3029 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
3030 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
3031 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
3032 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
3034 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
3035 // fRunLocalReconstruction.Contains("ALL") ) {
3036 // newRunLocalReconstruction += detName ;
3037 // newRunLocalReconstruction += " " ;
3039 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
3040 // fRunTracking.Contains("ALL") ) {
3041 // newRunTracking += detName ;
3042 // newRunTracking += " " ;
3044 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
3045 // fFillESD.Contains("ALL") ) {
3046 // newFillESD += detName ;
3047 // newFillESD += " " ;
3051 // fRunLocalReconstruction = newRunLocalReconstruction ;
3052 // fRunTracking = newRunTracking ;
3053 // fFillESD = newFillESD ;
3056 //_____________________________________________________________________________
3057 Int_t AliReconstruction::GetDetIndex(const char* detector)
3059 // return the detector index corresponding to detector
3061 for (index = 0; index < kNDetectors ; index++) {
3062 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3067 //_____________________________________________________________________________
3068 Bool_t AliReconstruction::FinishPlaneEff() {
3070 // Here execute all the necessary operationis, at the end of the tracking phase,
3071 // in case that evaluation of PlaneEfficiencies was required for some detector.
3072 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3074 // This Preliminary version works only FOR ITS !!!!!
3075 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3078 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3081 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3082 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3083 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3084 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
3085 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3086 TString name=planeeff->GetName();
3088 TFile* pefile = TFile::Open(name, "RECREATE");
3089 ret=(Bool_t)planeeff->Write();
3091 if(planeeff->GetCreateHistos()) {
3092 TString hname=planeeff->GetName();
3093 hname+="Histo.root";
3094 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3097 if(fSPDTrackleter) {
3098 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3099 TString name="AliITSPlaneEffSPDtracklet.root";
3100 TFile* pefile = TFile::Open(name, "RECREATE");
3101 ret=(Bool_t)planeeff->Write();
3103 AliESDEvent *dummy=NULL;
3104 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3109 //_____________________________________________________________________________
3110 Bool_t AliReconstruction::InitPlaneEff() {
3112 // Here execute all the necessary operations, before of the tracking phase,
3113 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3114 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3115 // which should be updated/recalculated.
3117 // This Preliminary version will work only FOR ITS !!!!!
3118 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3121 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3123 AliWarning(Form("Implementation of this method not yet completed !! Method return kTRUE"));
3125 fSPDTrackleter = NULL;
3126 AliReconstructor* itsReconstructor = GetReconstructor(0);
3127 if (itsReconstructor) {
3128 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3130 if (fSPDTrackleter) {
3131 AliInfo("Trackleter for SPD has been created");
3137 //_____________________________________________________________________________
3138 Bool_t AliReconstruction::InitAliEVE()
3140 // This method should be called only in case
3141 // AliReconstruction is run
3142 // within the alieve environment.
3143 // It will initialize AliEVE in a way
3144 // so that it can visualize event processed
3145 // by AliReconstruction.
3146 // The return flag shows whenever the
3147 // AliEVE initialization was successful or not.
3150 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3151 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3152 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3154 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3155 gROOT->ProcessLine("alieve_online_init()");
3160 //_____________________________________________________________________________
3161 void AliReconstruction::RunAliEVE()
3163 // Runs AliEVE visualisation of
3164 // the current event.
3165 // Should be executed only after
3166 // successful initialization of AliEVE.
3168 AliInfo("Running AliEVE...");
3169 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3173 //_____________________________________________________________________________
3174 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3176 // Allows to run QA for a selected set of detectors
3177 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3178 // all selected detectors run the same selected tasks
3180 if (!detAndAction.Contains(":")) {
3181 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3185 Int_t colon = detAndAction.Index(":") ;
3186 fQADetectors = detAndAction(0, colon) ;
3187 if (fQADetectors.Contains("ALL") )
3188 fQADetectors = fFillESD ;
3189 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3190 if (fQATasks.Contains("ALL") ) {
3191 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3193 fQATasks.ToUpper() ;
3195 if ( fQATasks.Contains("RAW") )
3196 tempo = Form("%d ", AliQAv1::kRAWS) ;
3197 if ( fQATasks.Contains("DIGIT") )
3198 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3199 if ( fQATasks.Contains("RECPOINT") )
3200 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3201 if ( fQATasks.Contains("ESD") )
3202 tempo += Form("%d ", AliQAv1::kESDS) ;
3204 if (fQATasks.IsNull()) {
3205 AliInfo("No QA requested\n") ;
3210 TString tempo(fQATasks) ;
3211 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3212 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3213 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3214 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3215 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3220 //_____________________________________________________________________________
3221 Bool_t AliReconstruction::InitRecoParams()
3223 // The method accesses OCDB and retrieves all
3224 // the available reco-param objects from there.
3226 Bool_t isOK = kTRUE;
3228 if (fRecoParam.GetDetRecoParamArray(kNDetectors)) {
3229 AliInfo("Using custom GRP reconstruction parameters");
3232 AliInfo("Loading GRP reconstruction parameter objects");
3234 AliCDBPath path("GRP","Calib","RecoParam");
3235 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3237 AliWarning("Couldn't find GRP RecoParam entry in OCDB");
3241 TObject *recoParamObj = entry->GetObject();
3242 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3243 // GRP has a normal TobjArray of AliDetectorRecoParam objects
3244 // Registering them in AliRecoParam
3245 fRecoParam.AddDetRecoParamArray(kNDetectors,dynamic_cast<TObjArray*>(recoParamObj));
3247 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3248 // GRP has only onse set of reco parameters
3249 // Registering it in AliRecoParam
3250 AliInfo("Single set of GRP reconstruction parameters found");
3251 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3252 fRecoParam.AddDetRecoParam(kNDetectors,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3255 AliError("No valid GRP RecoParam object found in the OCDB");
3262 TString detStr = fLoadCDB;
3263 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3265 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3267 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3268 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3272 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3274 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3275 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3277 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3281 TObject *recoParamObj = entry->GetObject();
3282 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3283 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3284 // Registering them in AliRecoParam
3285 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3287 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3288 // The detector has only onse set of reco parameters
3289 // Registering it in AliRecoParam
3290 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3291 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3292 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3295 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3299 // FIX ME: We have to disable the unloading of reco-param CDB
3300 // entries because QA framework is using them. Has to be fix in
3301 // a way that the QA takes the objects already constructed in
3303 // AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3307 if (AliDebugLevel() > 0) fRecoParam.Print();
3312 //_____________________________________________________________________________
3313 Bool_t AliReconstruction::GetEventInfo()
3315 // Fill the event info object
3317 AliCodeTimerAuto("")
3319 AliCentralTrigger *aCTP = NULL;
3321 fEventInfo.SetEventType(fRawReader->GetType());
3323 ULong64_t mask = fRawReader->GetClassMask();
3324 fEventInfo.SetTriggerMask(mask);
3325 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3326 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3328 aCTP = new AliCentralTrigger();
3329 TString configstr("");
3330 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3331 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3335 aCTP->SetClassMask(mask);
3336 aCTP->SetClusterMask(clmask);
3339 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3341 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3342 aCTP = fRunLoader->GetTrigger();
3343 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3344 // get inputs from actp - just get
3345 AliESDHeader* esdheader = fesd->GetHeader();
3346 esdheader->SetL0TriggerInputs(aCTP->GetL0TriggerInputs());
3347 esdheader->SetL1TriggerInputs(aCTP->GetL1TriggerInputs());
3348 esdheader->SetL2TriggerInputs(aCTP->GetL2TriggerInputs());
3349 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3352 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3357 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3359 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3360 if (fRawReader) delete aCTP;
3364 UChar_t clustmask = 0;
3366 ULong64_t trmask = fEventInfo.GetTriggerMask();
3367 const TObjArray& classesArray = config->GetClasses();
3368 Int_t nclasses = classesArray.GetEntriesFast();
3369 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3370 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3372 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3373 fesd->SetTriggerClass(trclass->GetName(),trindex);
3374 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3375 if (trmask & (1ull << trindex)) {
3377 trclasses += trclass->GetName();
3379 clustmask |= trclass->GetCluster()->GetClusterMask();
3383 fEventInfo.SetTriggerClasses(trclasses);
3385 // Set the information in ESD
3386 fesd->SetTriggerMask(trmask);
3387 fesd->SetTriggerCluster(clustmask);
3389 if (!aCTP->CheckTriggeredDetectors()) {
3390 if (fRawReader) delete aCTP;
3394 if (fRawReader) delete aCTP;
3396 // We have to fill also the HLT decision here!!
3402 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3404 // Match the detector list found in the rec.C or the default 'ALL'
3405 // to the list found in the GRP (stored there by the shuttle PP which
3406 // gets the information from ECS)
3407 static TString resultList;
3408 TString detList = detectorList;
3412 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3413 if ((detectorMask >> iDet) & 0x1) {
3414 TString det = AliDAQ::OfflineModuleName(iDet);
3415 if ((detList.CompareTo("ALL") == 0) ||
3416 ((detList.BeginsWith("ALL ") ||
3417 detList.EndsWith(" ALL") ||
3418 detList.Contains(" ALL ")) &&
3419 !(detList.BeginsWith("-"+det+" ") ||
3420 detList.EndsWith(" -"+det) ||
3421 detList.Contains(" -"+det+" "))) ||
3422 (detList.CompareTo(det) == 0) ||
3423 detList.BeginsWith(det+" ") ||
3424 detList.EndsWith(" "+det) ||
3425 detList.Contains( " "+det+" " )) {
3426 if (!resultList.EndsWith(det + " ")) {
3435 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3436 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3437 if ((detList.CompareTo("ALL") == 0) ||
3438 ((detList.BeginsWith("ALL ") ||
3439 detList.EndsWith(" ALL") ||
3440 detList.Contains(" ALL ")) &&
3441 !(detList.BeginsWith("-"+hltDet+" ") ||
3442 detList.EndsWith(" -"+hltDet) ||
3443 detList.Contains(" -"+hltDet+" "))) ||
3444 (detList.CompareTo(hltDet) == 0) ||
3445 detList.BeginsWith(hltDet+" ") ||
3446 detList.EndsWith(" "+hltDet) ||
3447 detList.Contains( " "+hltDet+" " )) {
3448 resultList += hltDet;
3452 return resultList.Data();
3456 //______________________________________________________________________________
3457 void AliReconstruction::Abort(const char *method, EAbort what)
3459 // Abort processing. If what = kAbortProcess, the Process() loop will be
3460 // aborted. If what = kAbortFile, the current file in a chain will be
3461 // aborted and the processing will continue with the next file, if there
3462 // is no next file then Process() will be aborted. Abort() can also be
3463 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3464 // the SlaveTerminate() and Terminate() are always called. The abort flag
3465 // can be checked in these methods using GetAbort().
3467 // The method is overwritten in AliReconstruction for better handling of
3468 // reco specific errors
3470 if (!fStopOnError) return;
3474 TString whyMess = method;
3475 whyMess += " failed! Aborting...";
3477 AliError(whyMess.Data());
3480 TString mess = "Abort";
3481 if (fAbort == kAbortProcess)
3482 mess = "AbortProcess";
3483 else if (fAbort == kAbortFile)
3486 Info(mess, whyMess.Data());