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 <TMessage.h>
133 #include "AliAlignObj.h"
134 #include "AliCDBEntry.h"
135 #include "AliCDBManager.h"
136 #include "AliCDBStorage.h"
137 #include "AliCTPRawStream.h"
138 #include "AliCascadeVertexer.h"
139 #include "AliCentralTrigger.h"
140 #include "AliCodeTimer.h"
142 #include "AliDetectorRecoParam.h"
143 #include "AliESDCaloCells.h"
144 #include "AliESDCaloCluster.h"
145 #include "AliESDEvent.h"
146 #include "AliESDMuonTrack.h"
147 #include "AliESDPmdTrack.h"
148 #include "AliESDTagCreator.h"
149 #include "AliESDVertex.h"
150 #include "AliESDcascade.h"
151 #include "AliESDfriend.h"
152 #include "AliESDkink.h"
153 #include "AliESDpid.h"
154 #include "AliESDtrack.h"
155 #include "AliESDtrack.h"
156 #include "AliEventInfo.h"
157 #include "AliGRPObject.h"
158 #include "AliGRPRecoParam.h"
159 #include "AliGenEventHeader.h"
160 #include "AliGeomManager.h"
161 #include "AliGlobalQADataMaker.h"
162 #include "AliHeader.h"
165 #include "AliMultiplicity.h"
167 #include "AliPlaneEff.h"
169 #include "AliQADataMakerRec.h"
170 #include "AliQAManager.h"
171 #include "AliRawVEvent.h"
172 #include "AliRawEventHeaderBase.h"
173 #include "AliRawHLTManager.h"
174 #include "AliRawReaderDate.h"
175 #include "AliRawReaderFile.h"
176 #include "AliRawReaderRoot.h"
177 #include "AliReconstruction.h"
178 #include "AliReconstructor.h"
180 #include "AliRunInfo.h"
181 #include "AliRunLoader.h"
182 #include "AliSysInfo.h" // memory snapshots
183 #include "AliTrackPointArray.h"
184 #include "AliTracker.h"
185 #include "AliTriggerClass.h"
186 #include "AliTriggerCluster.h"
187 #include "AliTriggerIR.h"
188 #include "AliTriggerConfiguration.h"
189 #include "AliV0vertexer.h"
190 #include "AliVertexer.h"
191 #include "AliVertexerTracks.h"
193 ClassImp(AliReconstruction)
195 //_____________________________________________________________________________
196 const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
198 //_____________________________________________________________________________
199 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
201 fRunVertexFinder(kTRUE),
202 fRunVertexFinderTracks(kTRUE),
203 fRunHLTTracking(kFALSE),
204 fRunMuonTracking(kFALSE),
206 fRunCascadeFinder(kTRUE),
207 fStopOnError(kFALSE),
208 fWriteAlignmentData(kFALSE),
209 fWriteESDfriend(kFALSE),
210 fFillTriggerESD(kTRUE),
218 fRunLocalReconstruction("ALL"),
222 fUseTrackingErrorsForAlignment(""),
223 fGAliceFileName(gAliceFilename),
229 fNumberOfEventsPerFile((UInt_t)-1),
231 fLoadAlignFromCDB(kTRUE),
232 fLoadAlignData("ALL"),
240 fParentRawReader(NULL),
244 fSPDTrackleter(NULL),
246 fDiamondProfileSPD(NULL),
247 fDiamondProfile(NULL),
248 fDiamondProfileTPC(NULL),
249 fListOfCosmicTriggers(NULL),
253 fAlignObjArray(NULL),
257 fInitCDBCalled(kFALSE),
258 fSetRunNumberFromDataCalled(kFALSE),
263 fSameQACycle(kFALSE),
264 fInitQACalled(kFALSE),
265 fWriteQAExpertData(kTRUE),
266 fRunPlaneEff(kFALSE),
275 fIsNewRunLoader(kFALSE),
279 // create reconstruction object with default parameters
282 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
283 fReconstructor[iDet] = NULL;
284 fLoader[iDet] = NULL;
285 fTracker[iDet] = NULL;
287 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
288 fQACycles[iDet] = 999999 ;
289 fQAWriteExpert[iDet] = kFALSE ;
295 //_____________________________________________________________________________
296 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
298 fRunVertexFinder(rec.fRunVertexFinder),
299 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
300 fRunHLTTracking(rec.fRunHLTTracking),
301 fRunMuonTracking(rec.fRunMuonTracking),
302 fRunV0Finder(rec.fRunV0Finder),
303 fRunCascadeFinder(rec.fRunCascadeFinder),
304 fStopOnError(rec.fStopOnError),
305 fWriteAlignmentData(rec.fWriteAlignmentData),
306 fWriteESDfriend(rec.fWriteESDfriend),
307 fFillTriggerESD(rec.fFillTriggerESD),
309 fCleanESD(rec.fCleanESD),
310 fV0DCAmax(rec.fV0DCAmax),
311 fV0CsPmin(rec.fV0CsPmin),
315 fRunLocalReconstruction(rec.fRunLocalReconstruction),
316 fRunTracking(rec.fRunTracking),
317 fFillESD(rec.fFillESD),
318 fLoadCDB(rec.fLoadCDB),
319 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
320 fGAliceFileName(rec.fGAliceFileName),
321 fRawInput(rec.fRawInput),
322 fESDOutput(rec.fESDOutput),
323 fEquipIdMap(rec.fEquipIdMap),
324 fFirstEvent(rec.fFirstEvent),
325 fLastEvent(rec.fLastEvent),
326 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
328 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
329 fLoadAlignData(rec.fLoadAlignData),
330 fUseHLTData(rec.fUseHLTData),
337 fParentRawReader(NULL),
339 fRecoParam(rec.fRecoParam),
341 fSPDTrackleter(NULL),
343 fDiamondProfileSPD(rec.fDiamondProfileSPD),
344 fDiamondProfile(rec.fDiamondProfile),
345 fDiamondProfileTPC(rec.fDiamondProfileTPC),
346 fListOfCosmicTriggers(NULL),
350 fAlignObjArray(rec.fAlignObjArray),
351 fCDBUri(rec.fCDBUri),
352 fQARefUri(rec.fQARefUri),
354 fInitCDBCalled(rec.fInitCDBCalled),
355 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
356 fQADetectors(rec.fQADetectors),
357 fQATasks(rec.fQATasks),
359 fRunGlobalQA(rec.fRunGlobalQA),
360 fSameQACycle(rec.fSameQACycle),
361 fInitQACalled(rec.fInitQACalled),
362 fWriteQAExpertData(rec.fWriteQAExpertData),
363 fRunPlaneEff(rec.fRunPlaneEff),
372 fIsNewRunLoader(rec.fIsNewRunLoader),
378 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
379 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
381 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
382 fReconstructor[iDet] = NULL;
383 fLoader[iDet] = NULL;
384 fTracker[iDet] = NULL;
387 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
388 fQACycles[iDet] = rec.fQACycles[iDet];
389 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
392 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
393 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
397 //_____________________________________________________________________________
398 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
400 // assignment operator
401 // Used in PROOF mode
402 // Be very careful while modifing it!
403 // Simple rules to follow:
404 // for persistent data members - use their assignment operators
405 // for non-persistent ones - do nothing or take the default values from constructor
406 // TSelector members should not be touched
407 if(&rec == this) return *this;
409 fRunVertexFinder = rec.fRunVertexFinder;
410 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
411 fRunHLTTracking = rec.fRunHLTTracking;
412 fRunMuonTracking = rec.fRunMuonTracking;
413 fRunV0Finder = rec.fRunV0Finder;
414 fRunCascadeFinder = rec.fRunCascadeFinder;
415 fStopOnError = rec.fStopOnError;
416 fWriteAlignmentData = rec.fWriteAlignmentData;
417 fWriteESDfriend = rec.fWriteESDfriend;
418 fFillTriggerESD = rec.fFillTriggerESD;
420 fCleanESD = rec.fCleanESD;
421 fV0DCAmax = rec.fV0DCAmax;
422 fV0CsPmin = rec.fV0CsPmin;
426 fRunLocalReconstruction = rec.fRunLocalReconstruction;
427 fRunTracking = rec.fRunTracking;
428 fFillESD = rec.fFillESD;
429 fLoadCDB = rec.fLoadCDB;
430 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
431 fGAliceFileName = rec.fGAliceFileName;
432 fRawInput = rec.fRawInput;
433 fESDOutput = rec.fESDOutput;
434 fEquipIdMap = rec.fEquipIdMap;
435 fFirstEvent = rec.fFirstEvent;
436 fLastEvent = rec.fLastEvent;
437 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
439 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
440 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
443 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
444 fLoadAlignData = rec.fLoadAlignData;
445 fUseHLTData = rec.fUseHLTData;
447 delete fRunInfo; fRunInfo = NULL;
448 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
450 fEventInfo = rec.fEventInfo;
452 delete fRunScalers; fRunScalers = NULL;
453 if (rec.fRunScalers) fRunScalers = new AliTriggerRunScalers(*rec.fRunScalers);
458 fParentRawReader = NULL;
460 fRecoParam = rec.fRecoParam;
462 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
463 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
464 delete fLoader[iDet]; fLoader[iDet] = NULL;
465 delete fTracker[iDet]; fTracker[iDet] = NULL;
468 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
469 fQACycles[iDet] = rec.fQACycles[iDet];
470 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
473 delete fSPDTrackleter; fSPDTrackleter = NULL;
475 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
476 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
477 delete fDiamondProfile; fDiamondProfile = NULL;
478 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
479 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
480 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
482 delete fListOfCosmicTriggers; fListOfCosmicTriggers = NULL;
483 if (rec.fListOfCosmicTriggers) fListOfCosmicTriggers = (THashTable*)((rec.fListOfCosmicTriggers)->Clone());
485 delete fGRPData; fGRPData = NULL;
486 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
487 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
489 delete fAlignObjArray; fAlignObjArray = NULL;
492 fQARefUri = rec.fQARefUri;
493 fSpecCDBUri.Delete();
494 fInitCDBCalled = rec.fInitCDBCalled;
495 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
496 fQADetectors = rec.fQADetectors;
497 fQATasks = rec.fQATasks;
499 fRunGlobalQA = rec.fRunGlobalQA;
500 fSameQACycle = rec.fSameQACycle;
501 fInitQACalled = rec.fInitQACalled;
502 fWriteQAExpertData = rec.fWriteQAExpertData;
503 fRunPlaneEff = rec.fRunPlaneEff;
512 fIsNewRunLoader = rec.fIsNewRunLoader;
519 //_____________________________________________________________________________
520 AliReconstruction::~AliReconstruction()
525 if (fListOfCosmicTriggers) {
526 fListOfCosmicTriggers->Delete();
527 delete fListOfCosmicTriggers;
532 if (fAlignObjArray) {
533 fAlignObjArray->Delete();
534 delete fAlignObjArray;
536 fSpecCDBUri.Delete();
538 AliCodeTimer::Instance()->Print();
541 //_____________________________________________________________________________
542 void AliReconstruction::InitQA()
544 //Initialize the QA and start of cycle
545 AliCodeTimerAuto("");
547 if (fInitQACalled) return;
548 fInitQACalled = kTRUE;
550 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kRECMODE) ;
551 if (fWriteQAExpertData)
552 qam->SetWriteExpert() ;
554 if (qam->IsDefaultStorageSet()) {
555 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
556 AliWarning("Default QA reference storage has been already set !");
557 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fQARefUri.Data()));
558 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
559 fQARefUri = qam->GetDefaultStorage()->GetURI();
561 if (fQARefUri.Length() > 0) {
562 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
563 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
564 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
566 fQARefUri="local://$ALICE_ROOT/QAref";
567 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
568 AliWarning("Default QA refeference storage not yet set !!!!");
569 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
570 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
573 qam->SetDefaultStorage(fQARefUri);
577 qam->SetActiveDetectors(fQADetectors) ;
578 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) {
579 qam->SetCycleLength(AliQAv1::DETECTORINDEX_t(det), fQACycles[det]) ;
580 qam->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
582 if (!fRawReader && !fInput && fQATasks.Contains(AliQAv1::kRAWS))
583 fQATasks.ReplaceAll(Form("%d",AliQAv1::kRAWS), "") ;
584 qam->SetTasks(fQATasks) ;
585 qam->InitQADataMaker(AliCDBManager::Instance()->GetRun()) ;
588 Bool_t sameCycle = kFALSE ;
589 AliQADataMaker *qadm = qam->GetQADataMaker(AliQAv1::kGLOBAL);
590 AliInfo(Form("Initializing the global QA data maker"));
591 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS))) {
592 qadm->StartOfCycle(AliQAv1::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
593 TObjArray **arr=qadm->Init(AliQAv1::kRECPOINTS);
594 AliTracker::SetResidualsArray(arr);
597 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS))) {
598 qadm->StartOfCycle(AliQAv1::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
599 qadm->Init(AliQAv1::kESDS);
602 AliSysInfo::AddStamp("InitQA") ;
605 //_____________________________________________________________________________
606 void AliReconstruction::MergeQA(const char *fileName)
608 //Initialize the QA and start of cycle
609 AliCodeTimerAuto("") ;
610 AliQAManager::QAManager()->Merge(AliCDBManager::Instance()->GetRun(),fileName) ;
611 AliSysInfo::AddStamp("MergeQA") ;
614 //_____________________________________________________________________________
615 void AliReconstruction::InitCDB()
617 // activate a default CDB storage
618 // First check if we have any CDB storage set, because it is used
619 // to retrieve the calibration and alignment constants
620 AliCodeTimerAuto("");
622 if (fInitCDBCalled) return;
623 fInitCDBCalled = kTRUE;
625 AliCDBManager* man = AliCDBManager::Instance();
626 if (man->IsDefaultStorageSet())
628 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
629 AliWarning("Default CDB storage has been already set !");
630 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
631 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
632 fCDBUri = man->GetDefaultStorage()->GetURI();
635 if (fCDBUri.Length() > 0)
637 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
638 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
639 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
641 fCDBUri="local://$ALICE_ROOT/OCDB";
642 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
643 AliWarning("Default CDB storage not yet set !!!!");
644 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
645 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
648 man->SetDefaultStorage(fCDBUri);
651 // Now activate the detector specific CDB storage locations
652 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
653 TObject* obj = fSpecCDBUri[i];
655 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
656 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
657 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
658 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
660 AliSysInfo::AddStamp("InitCDB");
663 //_____________________________________________________________________________
664 void AliReconstruction::SetDefaultStorage(const char* uri) {
665 // Store the desired default CDB storage location
666 // Activate it later within the Run() method
672 //_____________________________________________________________________________
673 void AliReconstruction::SetQARefDefaultStorage(const char* uri) {
674 // Store the desired default CDB storage location
675 // Activate it later within the Run() method
678 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
681 //_____________________________________________________________________________
682 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
683 // Store a detector-specific CDB storage location
684 // Activate it later within the Run() method
686 AliCDBPath aPath(calibType);
687 if(!aPath.IsValid()){
688 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
689 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
690 if(!strcmp(calibType, fgkDetectorName[iDet])) {
691 aPath.SetPath(Form("%s/*", calibType));
692 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
696 if(!aPath.IsValid()){
697 AliError(Form("Not a valid path or detector: %s", calibType));
702 // // check that calibType refers to a "valid" detector name
703 // Bool_t isDetector = kFALSE;
704 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
705 // TString detName = fgkDetectorName[iDet];
706 // if(aPath.GetLevel0() == detName) {
707 // isDetector = kTRUE;
713 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
717 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
718 if (obj) fSpecCDBUri.Remove(obj);
719 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
723 //_____________________________________________________________________________
724 Bool_t AliReconstruction::SetRunNumberFromData()
726 // The method is called in Run() in order
727 // to set a correct run number.
728 // In case of raw data reconstruction the
729 // run number is taken from the raw data header
731 if (fSetRunNumberFromDataCalled) return kTRUE;
732 fSetRunNumberFromDataCalled = kTRUE;
734 AliCDBManager* man = AliCDBManager::Instance();
737 if(fRawReader->NextEvent()) {
738 if(man->GetRun() > 0) {
739 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
741 man->SetRun(fRawReader->GetRunNumber());
742 fRawReader->RewindEvents();
745 if(man->GetRun() > 0) {
746 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
749 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
755 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
757 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
762 // read run number from gAlice
763 if(rl->GetHeader()) {
764 man->SetRun(rl->GetHeader()->GetRun());
769 AliError("Neither run-loader header nor RawReader objects are found !");
781 //_____________________________________________________________________________
782 void AliReconstruction::SetCDBLock() {
783 // Set CDB lock: from now on it is forbidden to reset the run number
784 // or the default storage or to activate any further storage!
786 AliCDBManager::Instance()->SetLock(1);
789 //_____________________________________________________________________________
790 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
792 // Read the alignment objects from CDB.
793 // Each detector is supposed to have the
794 // alignment objects in DET/Align/Data CDB path.
795 // All the detector objects are then collected,
796 // sorted by geometry level (starting from ALIC) and
797 // then applied to the TGeo geometry.
798 // Finally an overlaps check is performed.
800 // Load alignment data from CDB and fill fAlignObjArray
801 if(fLoadAlignFromCDB){
803 TString detStr = detectors;
804 TString loadAlObjsListOfDets = "";
806 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
807 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
808 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
810 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
812 loadAlObjsListOfDets += fgkDetectorName[iDet];
813 loadAlObjsListOfDets += " ";
815 } // end loop over detectors
817 if(AliGeomManager::GetNalignable("GRP") != 0)
818 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
819 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
820 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
822 // Check if the array with alignment objects was
823 // provided by the user. If yes, apply the objects
824 // to the present TGeo geometry
825 if (fAlignObjArray) {
826 if (gGeoManager && gGeoManager->IsClosed()) {
827 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
828 AliError("The misalignment of one or more volumes failed!"
829 "Compare the list of simulated detectors and the list of detector alignment data!");
834 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
840 if (fAlignObjArray) {
841 fAlignObjArray->Delete();
842 delete fAlignObjArray; fAlignObjArray=NULL;
848 //_____________________________________________________________________________
849 void AliReconstruction::SetGAliceFile(const char* fileName)
851 // set the name of the galice file
853 fGAliceFileName = fileName;
856 //_____________________________________________________________________________
857 void AliReconstruction::SetInput(const char* input)
859 // In case the input string starts with 'mem://', we run in an online mode
860 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
861 // file is assumed. One can give as an input:
862 // mem://: - events taken from DAQ monitoring libs online
864 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
865 if (input) fRawInput = input;
868 //_____________________________________________________________________________
869 void AliReconstruction::SetOutput(const char* output)
871 // Set the output ESD filename
872 // 'output' is a normalt ROOT url
873 // The method is used in case of raw-data reco with PROOF
874 if (output) fESDOutput.SetUrl(output);
877 //_____________________________________________________________________________
878 void AliReconstruction::SetOption(const char* detector, const char* option)
880 // set options for the reconstruction of a detector
882 TObject* obj = fOptions.FindObject(detector);
883 if (obj) fOptions.Remove(obj);
884 fOptions.Add(new TNamed(detector, option));
887 //_____________________________________________________________________________
888 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
890 // Set custom reconstruction parameters for a given detector
891 // Single set of parameters for all the events
893 // First check if the reco-params are global
894 if(!strcmp(detector, "GRP")) {
896 fRecoParam.AddDetRecoParam(kNDetectors,par);
900 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
901 if(!strcmp(detector, fgkDetectorName[iDet])) {
903 fRecoParam.AddDetRecoParam(iDet,par);
910 //_____________________________________________________________________________
911 Bool_t AliReconstruction::SetFieldMap(Float_t l3Cur, Float_t diCur, Float_t l3Pol,
912 Float_t diPol, Int_t convention, Bool_t uniform,
913 Float_t beamenergy, const Char_t *beamtype, const Char_t *path)
915 //------------------------------------------------
916 // The magnetic field map, defined externally...
917 // L3 current 30000 A -> 0.5 T
918 // L3 current 12000 A -> 0.2 T
919 // dipole current 6000 A
920 // The polarities must match the convention (LHC or DCS2008)
921 // unless the special uniform map was used for MC
922 //------------------------------------------------
923 const Float_t l3NominalCurrent1=30000.; // (A)
924 const Float_t l3NominalCurrent2=12000.; // (A)
925 const Float_t diNominalCurrent =6000. ; // (A)
927 const Float_t tolerance=0.03; // relative current tolerance
928 const Float_t zero=77.; // "zero" current (A)
933 l3Cur = TMath::Abs(l3Cur);
934 diCur = TMath::Abs(diCur);
936 if (TMath::Abs((sclDip=diCur/diNominalCurrent)-1.) > tolerance && !uniform) {
937 if (diCur <= zero) sclDip = 0.; // some small current.. -> Dipole OFF
939 AliError(Form("Wrong dipole current (%f A)!",diCur));
945 // special treatment of special MC with uniform mag field (normalized to 0.5 T)
946 // no check for scaling/polarities are done
947 map = AliMagF::k5kGUniform;
948 sclL3 = l3Cur/l3NominalCurrent1;
951 if (TMath::Abs((sclL3=l3Cur/l3NominalCurrent1)-1.) < tolerance) map = AliMagF::k5kG;
952 else if (TMath::Abs((sclL3=l3Cur/l3NominalCurrent2)-1.) < tolerance) map = AliMagF::k2kG;
953 else if (l3Cur <= zero) { sclL3 = 0; map = AliMagF::k5kGUniform;}
955 AliError(Form("Wrong L3 current (%f A)!",l3Cur));
960 if (sclDip!=0 && (map==AliMagF::k5kG || map==AliMagF::k2kG) &&
961 ((convention==AliMagF::kConvLHC && l3Pol!=diPol) ||
962 (convention==AliMagF::kConvDCS2008 && l3Pol==diPol)) ) {
963 AliError(Form("Wrong combination for L3/Dipole polarities (%c/%c) in %s convention",
964 l3Pol>0?'+':'-',diPol>0?'+':'-',convention==AliMagF::kConvDCS2008?"DCS2008":"LHC"));
968 if (l3Pol<0) sclL3 = -sclL3;
969 if (diPol<0) sclDip = -sclDip;
971 AliMagF::BeamType_t btype = AliMagF::kNoBeamField;
972 TString btypestr = beamtype;
974 TPRegexp protonBeam("(proton|p)\\s*-?\\s*\\1");
975 TPRegexp ionBeam("(lead|pb|ion|a)\\s*-?\\s*\\1");
976 if (btypestr.Contains(ionBeam)) btype = AliMagF::kBeamTypeAA;
977 else if (btypestr.Contains(protonBeam)) btype = AliMagF::kBeamTypepp;
978 else AliInfo(Form("Assume no LHC magnet field for the beam type %s, ",beamtype));
981 sprintf(ttl,"L3: %+5d Dip: %+4d kA; %s | Polarities in %s convention",(int)TMath::Sign(l3Cur,float(sclL3)),
982 (int)TMath::Sign(diCur,float(sclDip)),uniform ? " Constant":"",
983 convention==AliMagF::kConvLHC ? "LHC":"DCS2008");
984 // LHC and DCS08 conventions have opposite dipole polarities
985 if ( AliMagF::GetPolarityConvention() != convention) sclDip = -sclDip;
987 AliMagF* fld = new AliMagF("MagneticFieldMap", ttl, 2, sclL3, sclDip, 10., map, path,
989 TGeoGlobalMagField::Instance()->SetField( fld );
990 TGeoGlobalMagField::Instance()->Lock();
995 Bool_t AliReconstruction::InitGRP() {
996 //------------------------------------
997 // Initialization of the GRP entry
998 //------------------------------------
999 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
1003 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
1006 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
1008 fGRPData = new AliGRPObject();
1009 fGRPData->ReadValuesFromMap(m);
1013 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
1014 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
1018 // FIX ME: The unloading of GRP entry is temporarily disabled
1019 // because ZDC and VZERO are using it in order to initialize
1020 // their reconstructor objects. In the future one has to think
1021 // of propagating AliRunInfo to the reconstructors.
1022 // AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
1026 AliError("No GRP entry found in OCDB!");
1030 TString lhcState = fGRPData->GetLHCState();
1031 if (lhcState==AliGRPObject::GetInvalidString()) {
1032 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
1033 lhcState = "UNKNOWN";
1036 TString beamType = fGRPData->GetBeamType();
1037 if (beamType==AliGRPObject::GetInvalidString()) {
1038 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
1039 beamType = "UNKNOWN";
1042 Float_t beamEnergy = fGRPData->GetBeamEnergy();
1043 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
1044 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
1047 // LHC: "multiply by 120 to get the energy in MeV"
1048 beamEnergy *= 0.120;
1050 TString runType = fGRPData->GetRunType();
1051 if (runType==AliGRPObject::GetInvalidString()) {
1052 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
1053 runType = "UNKNOWN";
1056 Int_t activeDetectors = fGRPData->GetDetectorMask();
1057 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
1058 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
1059 activeDetectors = 1074790399;
1062 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1066 // Process the list of active detectors
1067 if (activeDetectors) {
1068 UInt_t detMask = activeDetectors;
1069 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1070 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1071 fFillESD = MatchDetectorList(fFillESD,detMask);
1072 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1073 fLoadCDB.Form("%s %s %s %s",
1074 fRunLocalReconstruction.Data(),
1075 fRunTracking.Data(),
1077 fQADetectors.Data());
1078 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1079 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1)) {
1080 // switch off the vertexer
1081 AliInfo("SPD is not in the list of active detectors. Vertexer switched off.");
1082 fRunVertexFinder = kFALSE;
1084 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1085 // switch off the reading of CTP raw-data payload
1086 if (fFillTriggerESD) {
1087 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1088 fFillTriggerESD = kFALSE;
1093 AliInfo("===================================================================================");
1094 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1095 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1096 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1097 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1098 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1099 AliInfo("===================================================================================");
1101 //*** Dealing with the magnetic field map
1102 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {
1103 if (TGeoGlobalMagField::Instance()->GetField()->TestBit(AliMagF::kOverrideGRP)) {
1104 AliInfo("ExpertMode!!! GRP information will be ignored !");
1105 AliInfo("ExpertMode!!! Running with the externally locked B field !");
1108 AliInfo("Destroying existing B field instance!");
1109 delete TGeoGlobalMagField::Instance();
1112 if ( !TGeoGlobalMagField::Instance()->IsLocked() ) {
1113 // Construct the field map out of the information retrieved from GRP.
1116 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1117 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1118 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1122 Char_t l3Polarity = fGRPData->GetL3Polarity();
1123 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1124 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1129 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1130 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1131 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1135 Char_t diPolarity = fGRPData->GetDipolePolarity();
1136 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1137 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1141 // read special bits for the polarity convention and map type
1142 Int_t polConvention = fGRPData->IsPolarityConventionLHC() ? AliMagF::kConvLHC : AliMagF::kConvDCS2008;
1143 Bool_t uniformB = fGRPData->IsUniformBMap();
1146 if ( !SetFieldMap(l3Current, diCurrent, l3Polarity ? -1:1, diPolarity ? -1:1,
1147 polConvention,uniformB,beamEnergy, beamType.Data()))
1148 AliFatal("Failed to creat a B field map ! Exiting...");
1149 AliInfo("Running with the B field constructed out of GRP !");
1151 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1154 //*** Get the diamond profiles from OCDB
1155 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1157 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1159 AliError("No SPD diamond profile found in OCDB!");
1162 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1164 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1166 AliError("No diamond profile found in OCDB!");
1169 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1171 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1173 AliError("No TPC diamond profile found in OCDB!");
1176 entry = AliCDBManager::Instance()->Get("GRP/Calib/CosmicTriggers");
1178 fListOfCosmicTriggers = dynamic_cast<THashTable*>(entry->GetObject());
1180 AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
1183 if (!fListOfCosmicTriggers) {
1184 AliWarning("Can not get list of cosmic triggers from OCDB! Cosmic event specie will be effectively disabled!");
1190 //_____________________________________________________________________________
1191 Bool_t AliReconstruction::LoadCDB()
1193 AliCodeTimerAuto("");
1195 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1197 TString detStr = fLoadCDB;
1198 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1199 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1200 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1204 //_____________________________________________________________________________
1205 Bool_t AliReconstruction::LoadTriggerScalersCDB()
1207 AliCodeTimerAuto("");
1209 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
1213 AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
1214 fRunScalers = dynamic_cast<AliTriggerRunScalers*> (entry->GetObject());
1216 if (fRunScalers->CorrectScalersOverflow() == 0) AliInfo("32bit Trigger counters corrected for overflow");
1221 //_____________________________________________________________________________
1222 Bool_t AliReconstruction::Run(const char* input)
1225 AliCodeTimerAuto("");
1228 if (GetAbort() != TSelector::kContinue) return kFALSE;
1230 TChain *chain = NULL;
1231 if (fRawReader && (chain = fRawReader->GetChain())) {
1236 gProof->Exec("TGrid::Connect(\"alien://\")",kTRUE);
1238 TMessage::EnableSchemaEvolutionForAll(kTRUE);
1239 gProof->Exec("TMessage::EnableSchemaEvolutionForAll(kTRUE)",kTRUE);
1241 gProof->AddInput(this);
1242 if (!fESDOutput.IsValid()) {
1243 fESDOutput.SetProtocol("root",kTRUE);
1244 fESDOutput.SetHost(gSystem->HostName());
1245 fESDOutput.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1247 AliInfo(Form("Output file with ESDs is %s",fESDOutput.GetUrl()));
1248 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",fESDOutput.GetUrl()));
1249 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1251 chain->Process("AliReconstruction");
1254 chain->Process(this);
1259 if (GetAbort() != TSelector::kContinue) return kFALSE;
1261 if (GetAbort() != TSelector::kContinue) return kFALSE;
1262 //******* The loop over events
1263 AliInfo("Starting looping over events");
1265 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1266 (fRawReader && fRawReader->NextEvent())) {
1267 if (!ProcessEvent(iEvent)) {
1268 Abort("ProcessEvent",TSelector::kAbortFile);
1274 if (GetAbort() != TSelector::kContinue) return kFALSE;
1276 if (GetAbort() != TSelector::kContinue) return kFALSE;
1282 //_____________________________________________________________________________
1283 void AliReconstruction::InitRawReader(const char* input)
1285 AliCodeTimerAuto("");
1287 // Init raw-reader and
1288 // set the input in case of raw data
1289 if (input) fRawInput = input;
1290 fRawReader = AliRawReader::Create(fRawInput.Data());
1292 AliInfo("Reconstruction will run over digits");
1294 if (!fEquipIdMap.IsNull() && fRawReader)
1295 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1297 if (!fUseHLTData.IsNull()) {
1298 // create the RawReaderHLT which performs redirection of HLT input data for
1299 // the specified detectors
1300 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1302 fParentRawReader=fRawReader;
1303 fRawReader=pRawReader;
1305 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1308 AliSysInfo::AddStamp("CreateRawReader");
1311 //_____________________________________________________________________________
1312 void AliReconstruction::InitRun(const char* input)
1314 // Initialization of raw-reader,
1315 // run number, CDB etc.
1316 AliCodeTimerAuto("");
1317 AliSysInfo::AddStamp("Start");
1319 // Initialize raw-reader if any
1320 InitRawReader(input);
1322 // Initialize the CDB storage
1325 // Set run number in CDBManager (if it is not already set by the user)
1326 if (!SetRunNumberFromData()) {
1327 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1331 // Set CDB lock: from now on it is forbidden to reset the run number
1332 // or the default storage or to activate any further storage!
1337 //_____________________________________________________________________________
1338 void AliReconstruction::Begin(TTree *)
1340 // Initialize AlReconstruction before
1341 // going into the event loop
1342 // Should follow the TSelector convention
1343 // i.e. initialize only the object on the client side
1344 AliCodeTimerAuto("");
1346 AliReconstruction *reco = NULL;
1348 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1351 AliSysInfo::AddStamp("ReadInputInBegin");
1354 // Import ideal TGeo geometry and apply misalignment
1356 TString geom(gSystem->DirName(fGAliceFileName));
1357 geom += "/geometry.root";
1358 AliGeomManager::LoadGeometry(geom.Data());
1360 Abort("LoadGeometry", TSelector::kAbortProcess);
1363 AliSysInfo::AddStamp("LoadGeom");
1364 TString detsToCheck=fRunLocalReconstruction;
1365 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1366 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1369 AliSysInfo::AddStamp("CheckGeom");
1372 if (!MisalignGeometry(fLoadAlignData)) {
1373 Abort("MisalignGeometry", TSelector::kAbortProcess);
1376 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1377 AliSysInfo::AddStamp("MisalignGeom");
1380 Abort("InitGRP", TSelector::kAbortProcess);
1383 AliSysInfo::AddStamp("InitGRP");
1386 Abort("LoadCDB", TSelector::kAbortProcess);
1389 AliSysInfo::AddStamp("LoadCDB");
1391 if (!LoadTriggerScalersCDB()) {
1392 Abort("LoadTriggerScalersCDB", TSelector::kAbortProcess);
1395 AliSysInfo::AddStamp("LoadTriggerScalersCDB");
1398 // Read the reconstruction parameters from OCDB
1399 if (!InitRecoParams()) {
1400 AliWarning("Not all detectors have correct RecoParam objects initialized");
1402 AliSysInfo::AddStamp("InitRecoParams");
1404 if (fInput && gProof) {
1405 if (reco) *reco = *this;
1407 gProof->AddInputData(gGeoManager,kTRUE);
1409 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1410 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1411 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1412 magFieldMap->SetName("MagneticFieldMap");
1413 gProof->AddInputData(magFieldMap,kTRUE);
1418 //_____________________________________________________________________________
1419 void AliReconstruction::SlaveBegin(TTree*)
1421 // Initialization related to run-loader,
1422 // vertexer, trackers, recontructors
1423 // In proof mode it is executed on the slave
1424 AliCodeTimerAuto("");
1426 TProofOutputFile *outProofFile = NULL;
1428 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1431 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1433 AliGeomManager::SetGeometry(tgeo);
1435 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1436 Int_t runNumber = -1;
1437 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1438 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1439 man->SetCacheFlag(kTRUE);
1440 man->SetLock(kTRUE);
1444 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1445 TGeoGlobalMagField::Instance()->SetField(map);
1447 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1448 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1449 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1450 fOutput->Add(outProofFile);
1452 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1455 // get the run loader
1456 if (!InitRunLoader()) {
1457 Abort("InitRunLoader", TSelector::kAbortProcess);
1460 AliSysInfo::AddStamp("LoadLoader");
1462 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1465 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1466 Abort("CreateTrackers", TSelector::kAbortProcess);
1469 AliSysInfo::AddStamp("CreateTrackers");
1471 // create the ESD output file and tree
1472 if (!outProofFile) {
1473 ffile = TFile::Open("AliESDs.root", "RECREATE");
1474 ffile->SetCompressionLevel(2);
1475 if (!ffile->IsOpen()) {
1476 Abort("OpenESDFile", TSelector::kAbortProcess);
1481 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1482 Abort(Form("Problems opening output PROOF file: %s/%s",
1483 outProofFile->GetDir(), outProofFile->GetFileName()),
1484 TSelector::kAbortProcess);
1489 ftree = new TTree("esdTree", "Tree with ESD objects");
1490 fesd = new AliESDEvent();
1491 fesd->CreateStdContent();
1493 fesd->WriteToTree(ftree);
1494 if (fWriteESDfriend) {
1496 // Since we add the branch manually we must
1497 // book and add it after WriteToTree
1498 // otherwise it is created twice,
1499 // once via writetotree and once here.
1500 // The case for AliESDfriend is now
1501 // caught also in AlIESDEvent::WriteToTree but
1502 // be careful when changing the name (AliESDfriend is not
1503 // a TNamed so we had to hardwire it)
1504 fesdf = new AliESDfriend();
1505 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1506 br->SetFile("AliESDfriends.root");
1507 fesd->AddObject(fesdf);
1509 ftree->GetUserInfo()->Add(fesd);
1511 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1512 fhltesd = new AliESDEvent();
1513 fhltesd->CreateStdContent();
1515 // read the ESD template from CDB
1516 // HLT is allowed to put non-std content to its ESD, the non-std
1517 // objects need to be created before invocation of WriteToTree in
1518 // order to create all branches. Initialization is done from an
1519 // ESD layout template in CDB
1520 AliCDBManager* man = AliCDBManager::Instance();
1521 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1522 AliCDBEntry* hltESDConfig=NULL;
1523 if (man->GetId(hltESDConfigPath)!=NULL &&
1524 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1525 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1527 // init all internal variables from the list of objects
1528 pESDLayout->GetStdContent();
1530 // copy content and create non-std objects
1531 *fhltesd=*pESDLayout;
1534 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1535 hltESDConfigPath.GetPath().Data()));
1539 fhltesd->WriteToTree(fhlttree);
1540 fhlttree->GetUserInfo()->Add(fhltesd);
1542 ProcInfo_t procInfo;
1543 gSystem->GetProcInfo(&procInfo);
1544 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1547 //Initialize the QA and start of cycle
1548 if (fRunQA || fRunGlobalQA)
1551 //Initialize the Plane Efficiency framework
1552 if (fRunPlaneEff && !InitPlaneEff()) {
1553 Abort("InitPlaneEff", TSelector::kAbortProcess);
1557 if (strcmp(gProgName,"alieve") == 0)
1558 fRunAliEVE = InitAliEVE();
1563 //_____________________________________________________________________________
1564 Bool_t AliReconstruction::Process(Long64_t entry)
1566 // run the reconstruction over a single entry
1567 // from the chain with raw data
1568 AliCodeTimerAuto("");
1570 TTree *currTree = fChain->GetTree();
1571 AliRawVEvent *event = NULL;
1572 currTree->SetBranchAddress("rawevent",&event);
1573 currTree->GetEntry(entry);
1574 fRawReader = new AliRawReaderRoot(event);
1575 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1583 //_____________________________________________________________________________
1584 void AliReconstruction::Init(TTree *tree)
1587 AliError("The input tree is not found!");
1593 //_____________________________________________________________________________
1594 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1596 // run the reconstruction over a single event
1597 // The event loop is steered in Run method
1599 AliCodeTimerAuto("");
1601 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1602 fRunLoader->SetEventNumber(iEvent);
1603 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1605 fRunLoader->TreeE()->Fill();
1606 if (fRawReader && fRawReader->UseAutoSaveESD())
1607 fRunLoader->TreeE()->AutoSave("SaveSelf");
1610 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1614 AliInfo(Form("processing event %d", iEvent));
1616 fRunLoader->GetEvent(iEvent);
1618 // Fill Event-info object
1620 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo,fListOfCosmicTriggers);
1621 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1623 // Set the reco-params
1625 TString detStr = fLoadCDB;
1626 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1627 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1628 AliReconstructor *reconstructor = GetReconstructor(iDet);
1629 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1630 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1631 reconstructor->SetRecoParam(par);
1632 reconstructor->SetEventInfo(&fEventInfo);
1634 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1635 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(par->GetEventSpecie())) ;
1643 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1644 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1646 // local single event reconstruction
1647 if (!fRunLocalReconstruction.IsNull()) {
1648 TString detectors=fRunLocalReconstruction;
1649 // run HLT event reconstruction first
1650 // ;-( IsSelected changes the string
1651 if (IsSelected("HLT", detectors) &&
1652 !RunLocalEventReconstruction("HLT")) {
1653 if (fStopOnError) {CleanUp(); return kFALSE;}
1655 detectors=fRunLocalReconstruction;
1656 detectors.ReplaceAll("HLT", "");
1657 if (!RunLocalEventReconstruction(detectors)) {
1658 if (fStopOnError) {CleanUp(); return kFALSE;}
1662 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1663 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1664 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1665 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1667 // Set magnetic field from the tracker
1668 fesd->SetMagneticField(AliTracker::GetBz());
1669 fhltesd->SetMagneticField(AliTracker::GetBz());
1671 // Set most probable pt, for B=0 tracking
1672 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1673 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1674 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1676 // Fill raw-data error log into the ESD
1677 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1680 if (fRunVertexFinder) {
1681 if (!RunVertexFinder(fesd)) {
1682 if (fStopOnError) {CleanUp(); return kFALSE;}
1686 // For Plane Efficiency: run the SPD trackleter
1687 if (fRunPlaneEff && fSPDTrackleter) {
1688 if (!RunSPDTrackleting(fesd)) {
1689 if (fStopOnError) {CleanUp(); return kFALSE;}
1694 if (!fRunTracking.IsNull()) {
1695 if (fRunMuonTracking) {
1696 if (!RunMuonTracking(fesd)) {
1697 if (fStopOnError) {CleanUp(); return kFALSE;}
1703 if (!fRunTracking.IsNull()) {
1704 if (!RunTracking(fesd)) {
1705 if (fStopOnError) {CleanUp(); return kFALSE;}
1710 if (!fFillESD.IsNull()) {
1711 TString detectors=fFillESD;
1712 // run HLT first and on hltesd
1713 // ;-( IsSelected changes the string
1714 if (IsSelected("HLT", detectors) &&
1715 !FillESD(fhltesd, "HLT")) {
1716 if (fStopOnError) {CleanUp(); return kFALSE;}
1719 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1720 if (detectors.Contains("ALL")) {
1722 for (Int_t idet=0; idet<kNDetectors; ++idet){
1723 detectors += fgkDetectorName[idet];
1727 detectors.ReplaceAll("HLT", "");
1728 if (!FillESD(fesd, detectors)) {
1729 if (fStopOnError) {CleanUp(); return kFALSE;}
1733 // fill Event header information from the RawEventHeader
1734 if (fRawReader){FillRawEventHeaderESD(fesd);}
1735 if (fRawReader){FillRawEventHeaderESD(fhltesd);}
1738 AliESDpid::MakePID(fesd);
1740 if (fFillTriggerESD) {
1741 if (!FillTriggerESD(fesd)) {
1742 if (fStopOnError) {CleanUp(); return kFALSE;}
1745 // Always fill scalers
1746 if (!FillTriggerScalers(fesd)) {
1747 if (fStopOnError) {CleanUp(); return kFALSE;}
1754 // Propagate track to the beam pipe (if not already done by ITS)
1756 const Int_t ntracks = fesd->GetNumberOfTracks();
1757 const Double_t kBz = fesd->GetMagneticField();
1758 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1761 UShort_t *selectedIdx=new UShort_t[ntracks];
1763 for (Int_t itrack=0; itrack<ntracks; itrack++){
1764 const Double_t kMaxStep = 1; //max step over the material
1767 AliESDtrack *track = fesd->GetTrack(itrack);
1768 if (!track) continue;
1770 AliExternalTrackParam *tpcTrack =
1771 (AliExternalTrackParam *)track->GetTPCInnerParam();
1775 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1778 Int_t n=trkArray.GetEntriesFast();
1779 selectedIdx[n]=track->GetID();
1780 trkArray.AddLast(tpcTrack);
1783 //Tracks refitted by ITS should already be at the SPD vertex
1784 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1787 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1788 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1793 // Improve the reconstructed primary vertex position using the tracks
1795 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1796 if(fesd->GetPrimaryVertexSPD()) {
1797 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1798 if(vtitle.Contains("cosmics")) {
1799 runVertexFinderTracks=kFALSE;
1803 if (runVertexFinderTracks) {
1804 // TPC + ITS primary vertex
1805 ftVertexer->SetITSMode();
1806 ftVertexer->SetConstraintOff();
1807 // get cuts for vertexer from AliGRPRecoParam
1809 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1810 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1811 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1812 ftVertexer->SetCuts(cutsVertexer);
1813 delete [] cutsVertexer; cutsVertexer = NULL;
1814 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1815 ftVertexer->SetVtxStart(fDiamondProfile);
1817 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1819 if (pvtx->GetStatus()) {
1820 fesd->SetPrimaryVertexTracks(pvtx);
1821 for (Int_t i=0; i<ntracks; i++) {
1822 AliESDtrack *t = fesd->GetTrack(i);
1823 t->RelateToVertex(pvtx, kBz, kVeryBig);
1828 // TPC-only primary vertex
1829 ftVertexer->SetTPCMode();
1830 ftVertexer->SetConstraintOff();
1831 // get cuts for vertexer from AliGRPRecoParam
1833 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1834 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1835 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1836 ftVertexer->SetCuts(cutsVertexer);
1837 delete [] cutsVertexer; cutsVertexer = NULL;
1838 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1839 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1841 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1843 if (pvtx->GetStatus()) {
1844 fesd->SetPrimaryVertexTPC(pvtx);
1845 for (Int_t i=0; i<ntracks; i++) {
1846 AliESDtrack *t = fesd->GetTrack(i);
1847 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1853 delete[] selectedIdx;
1855 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1860 AliV0vertexer vtxer;
1861 vtxer.Tracks2V0vertices(fesd);
1863 if (fRunCascadeFinder) {
1865 AliCascadeVertexer cvtxer;
1866 cvtxer.V0sTracks2CascadeVertices(fesd);
1871 if (fCleanESD) CleanESD(fesd);
1874 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1875 AliQAManager::QAManager()->RunOneEvent(fesd) ;
1878 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1879 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1880 if (qadm && fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1881 qadm->Exec(AliQAv1::kESDS, fesd);
1884 if (fWriteESDfriend) {
1885 // fesdf->~AliESDfriend();
1886 // new (fesdf) AliESDfriend(); // Reset...
1887 fesd->GetESDfriend(fesdf);
1891 // Auto-save the ESD tree in case of prompt reco @P2
1892 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1893 ftree->AutoSave("SaveSelf");
1894 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1895 if (friendfile) friendfile->Save();
1902 if (fRunAliEVE) RunAliEVE();
1906 if (fWriteESDfriend) {
1907 fesdf->~AliESDfriend();
1908 new (fesdf) AliESDfriend(); // Reset...
1911 ProcInfo_t procInfo;
1912 gSystem->GetProcInfo(&procInfo);
1913 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1916 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1917 if (fReconstructor[iDet]) {
1918 fReconstructor[iDet]->SetRecoParam(NULL);
1919 fReconstructor[iDet]->SetEventInfo(NULL);
1921 if (fTracker[iDet]) fTracker[iDet]->SetEventInfo(NULL);
1924 if (fRunQA || fRunGlobalQA)
1925 AliQAManager::QAManager()->Increment() ;
1930 //_____________________________________________________________________________
1931 void AliReconstruction::SlaveTerminate()
1933 // Finalize the run on the slave side
1934 // Called after the exit
1935 // from the event loop
1936 AliCodeTimerAuto("");
1938 if (fIsNewRunLoader) { // galice.root didn't exist
1939 fRunLoader->WriteHeader("OVERWRITE");
1940 fRunLoader->CdGAFile();
1941 fRunLoader->Write(0, TObject::kOverwrite);
1944 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1945 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1947 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1948 cdbMapCopy->SetOwner(1);
1949 cdbMapCopy->SetName("cdbMap");
1950 TIter iter(cdbMap->GetTable());
1953 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1954 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1955 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1956 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1959 TList *cdbListCopy = new TList();
1960 cdbListCopy->SetOwner(1);
1961 cdbListCopy->SetName("cdbList");
1963 TIter iter2(cdbList);
1966 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1967 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1970 ftree->GetUserInfo()->Add(cdbMapCopy);
1971 ftree->GetUserInfo()->Add(cdbListCopy);
1976 if (fWriteESDfriend)
1977 ftree->SetBranchStatus("ESDfriend*",0);
1978 // we want to have only one tree version number
1979 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1980 fhlttree->Write(fhlttree->GetName(),TObject::kOverwrite);
1982 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1983 if (fRunPlaneEff && !FinishPlaneEff()) {
1984 AliWarning("Finish PlaneEff evaluation failed");
1987 // End of cycle for the in-loop
1989 AliQAManager::QAManager()->EndOfCycle() ;
1992 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1994 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS)))
1995 qadm->EndOfCycle(AliQAv1::kRECPOINTS);
1996 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1997 qadm->EndOfCycle(AliQAv1::kESDS);
2002 if (fRunQA || fRunGlobalQA) {
2004 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
2005 TString qaOutputFile = outputFileName->GetTitle();
2006 qaOutputFile.ReplaceAll(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()),
2007 Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
2008 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
2009 qaProofFile->SetOutputFileName(qaOutputFile.Data());
2010 fOutput->Add(qaProofFile);
2011 MergeQA(qaProofFile->GetFileName());
2023 //_____________________________________________________________________________
2024 void AliReconstruction::Terminate()
2026 // Create tags for the events in the ESD tree (the ESD tree is always present)
2027 // In case of empty events the tags will contain dummy values
2028 AliCodeTimerAuto("");
2030 // Do not call the ESD tag creator in case of PROOF-based reconstruction
2032 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
2033 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
2036 // Cleanup of CDB manager: cache and active storages!
2037 AliCDBManager::Instance()->ClearCache();
2040 //_____________________________________________________________________________
2041 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
2043 // run the local reconstruction
2045 static Int_t eventNr=0;
2046 AliCodeTimerAuto("")
2048 TString detStr = detectors;
2049 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2050 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2051 AliReconstructor* reconstructor = GetReconstructor(iDet);
2052 if (!reconstructor) continue;
2053 AliLoader* loader = fLoader[iDet];
2054 // Matthias April 2008: temporary fix to run HLT reconstruction
2055 // although the HLT loader is missing
2056 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2058 reconstructor->Reconstruct(fRawReader, NULL);
2061 reconstructor->Reconstruct(dummy, NULL);
2066 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2069 // conversion of digits
2070 if (fRawReader && reconstructor->HasDigitConversion()) {
2071 AliInfo(Form("converting raw data digits into root objects for %s",
2072 fgkDetectorName[iDet]));
2073 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2074 // fgkDetectorName[iDet]));
2075 loader->LoadDigits("update");
2076 loader->CleanDigits();
2077 loader->MakeDigitsContainer();
2078 TTree* digitsTree = loader->TreeD();
2079 reconstructor->ConvertDigits(fRawReader, digitsTree);
2080 loader->WriteDigits("OVERWRITE");
2081 loader->UnloadDigits();
2083 // local reconstruction
2084 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2085 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2086 loader->LoadRecPoints("update");
2087 loader->CleanRecPoints();
2088 loader->MakeRecPointsContainer();
2089 TTree* clustersTree = loader->TreeR();
2090 if (fRawReader && !reconstructor->HasDigitConversion()) {
2091 reconstructor->Reconstruct(fRawReader, clustersTree);
2093 loader->LoadDigits("read");
2094 TTree* digitsTree = loader->TreeD();
2096 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2097 if (fStopOnError) return kFALSE;
2099 reconstructor->Reconstruct(digitsTree, clustersTree);
2101 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2102 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2105 loader->UnloadDigits();
2108 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2109 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2111 loader->WriteRecPoints("OVERWRITE");
2112 loader->UnloadRecPoints();
2113 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2115 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2116 AliError(Form("the following detectors were not found: %s",
2118 if (fStopOnError) return kFALSE;
2123 //_____________________________________________________________________________
2124 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2126 // run the SPD trackleting (for SPD efficiency purpouses)
2128 AliCodeTimerAuto("")
2130 Double_t vtxPos[3] = {0, 0, 0};
2131 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2133 TArrayF mcVertex(3);
2135 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2136 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2137 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2140 const AliESDVertex *vertex = esd->GetVertex();
2142 AliWarning("Vertex not found");
2145 vertex->GetXYZ(vtxPos);
2146 vertex->GetSigmaXYZ(vtxErr);
2147 if (fSPDTrackleter) {
2148 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2151 fLoader[0]->LoadRecPoints("read");
2152 TTree* tree = fLoader[0]->TreeR();
2154 AliError("Can't get the ITS cluster tree");
2157 fSPDTrackleter->LoadClusters(tree);
2158 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2160 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2161 AliError("AliITSTrackleterSPDEff Clusters2Tracks failed");
2162 // fLoader[0]->UnloadRecPoints();
2165 //fSPDTrackleter->UnloadRecPoints();
2167 AliWarning("SPDTrackleter not available");
2173 //_____________________________________________________________________________
2174 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2176 // run the barrel tracking
2178 AliCodeTimerAuto("")
2180 AliVertexer *vertexer = CreateVertexer();
2181 if (!vertexer) return kFALSE;
2183 AliInfo("running the ITS vertex finder");
2184 AliESDVertex* vertex = NULL;
2186 fLoader[0]->LoadRecPoints();
2187 TTree* cltree = fLoader[0]->TreeR();
2189 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2190 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2193 AliError("Can't get the ITS cluster tree");
2195 fLoader[0]->UnloadRecPoints();
2198 AliError("Can't get the ITS loader");
2201 AliWarning("Vertex not found");
2202 vertex = new AliESDVertex();
2203 vertex->SetName("default");
2206 vertex->SetName("reconstructed");
2211 vertex->GetXYZ(vtxPos);
2212 vertex->GetSigmaXYZ(vtxErr);
2214 esd->SetPrimaryVertexSPD(vertex);
2215 AliESDVertex *vpileup = NULL;
2216 Int_t novertices = 0;
2217 vpileup = vertexer->GetAllVertices(novertices);
2219 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2221 // if SPD multiplicity has been determined, it is stored in the ESD
2222 AliMultiplicity *mult = vertexer->GetMultiplicity();
2223 if(mult)esd->SetMultiplicity(mult);
2225 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2226 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2235 //_____________________________________________________________________________
2236 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2238 // run the HLT barrel tracking
2240 AliCodeTimerAuto("")
2243 AliError("Missing runLoader!");
2247 AliInfo("running HLT tracking");
2249 // Get a pointer to the HLT reconstructor
2250 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2251 if (!reconstructor) return kFALSE;
2254 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2255 TString detName = fgkDetectorName[iDet];
2256 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2257 reconstructor->SetOption(detName.Data());
2258 AliTracker *tracker = reconstructor->CreateTracker();
2260 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2261 if (fStopOnError) return kFALSE;
2265 Double_t vtxErr[3]={0.005,0.005,0.010};
2266 const AliESDVertex *vertex = esd->GetVertex();
2267 vertex->GetXYZ(vtxPos);
2268 tracker->SetVertex(vtxPos,vtxErr);
2270 fLoader[iDet]->LoadRecPoints("read");
2271 TTree* tree = fLoader[iDet]->TreeR();
2273 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2276 tracker->LoadClusters(tree);
2278 if (tracker->Clusters2Tracks(esd) != 0) {
2279 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2283 tracker->UnloadClusters();
2291 //_____________________________________________________________________________
2292 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2294 // run the muon spectrometer tracking
2296 AliCodeTimerAuto("")
2299 AliError("Missing runLoader!");
2302 Int_t iDet = 7; // for MUON
2304 AliInfo("is running...");
2306 // Get a pointer to the MUON reconstructor
2307 AliReconstructor *reconstructor = GetReconstructor(iDet);
2308 if (!reconstructor) return kFALSE;
2311 TString detName = fgkDetectorName[iDet];
2312 AliDebug(1, Form("%s tracking", detName.Data()));
2313 AliTracker *tracker = reconstructor->CreateTracker();
2315 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2320 fLoader[iDet]->LoadRecPoints("read");
2322 tracker->LoadClusters(fLoader[iDet]->TreeR());
2324 Int_t rv = tracker->Clusters2Tracks(esd);
2328 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2332 fLoader[iDet]->UnloadRecPoints();
2334 tracker->UnloadClusters();
2342 //_____________________________________________________________________________
2343 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2345 // run the barrel tracking
2346 static Int_t eventNr=0;
2347 AliCodeTimerAuto("")
2349 AliInfo("running tracking");
2351 // Set the event info which is used
2352 // by the trackers in order to obtain
2353 // information about read-out detectors,
2355 AliDebug(1, "Setting event info");
2356 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2357 if (!fTracker[iDet]) continue;
2358 fTracker[iDet]->SetEventInfo(&fEventInfo);
2361 //Fill the ESD with the T0 info (will be used by the TOF)
2362 if (fReconstructor[11] && fLoader[11]) {
2363 fLoader[11]->LoadRecPoints("READ");
2364 TTree *treeR = fLoader[11]->TreeR();
2366 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2370 // pass 1: TPC + ITS inwards
2371 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2372 if (!fTracker[iDet]) continue;
2373 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2376 fLoader[iDet]->LoadRecPoints("read");
2377 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2378 TTree* tree = fLoader[iDet]->TreeR();
2380 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2383 fTracker[iDet]->LoadClusters(tree);
2384 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2386 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2387 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2390 // preliminary PID in TPC needed by the ITS tracker
2392 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2393 AliESDpid::MakePID(esd);
2395 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2398 // pass 2: ALL backwards
2400 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2401 if (!fTracker[iDet]) continue;
2402 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2405 if (iDet > 1) { // all except ITS, TPC
2407 fLoader[iDet]->LoadRecPoints("read");
2408 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2409 tree = fLoader[iDet]->TreeR();
2411 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2414 fTracker[iDet]->LoadClusters(tree);
2415 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2419 if (iDet>1) // start filling residuals for the "outer" detectors
2421 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2422 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2424 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2425 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2426 if ( elem && (! elem->At(0)) ) {
2427 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2428 if (qadm) qadm->InitRecPointsForTracker() ;
2432 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2433 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2438 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2439 fTracker[iDet]->UnloadClusters();
2440 fLoader[iDet]->UnloadRecPoints();
2442 // updated PID in TPC needed by the ITS tracker -MI
2444 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2445 AliESDpid::MakePID(esd);
2447 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2449 //stop filling residuals for the "outer" detectors
2450 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2452 // pass 3: TRD + TPC + ITS refit inwards
2454 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2455 if (!fTracker[iDet]) continue;
2456 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2459 if (iDet<2) // start filling residuals for TPC and ITS
2461 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2462 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2464 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2465 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2466 if ( elem && (! elem->At(0)) ) {
2467 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2468 if (qadm) qadm->InitRecPointsForTracker() ;
2473 if (fTracker[iDet]->RefitInward(esd) != 0) {
2474 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2477 // run postprocessing
2478 if (fTracker[iDet]->PostProcess(esd) != 0) {
2479 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2482 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2485 // write space-points to the ESD in case alignment data output
2487 if (fWriteAlignmentData)
2488 WriteAlignmentData(esd);
2490 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2491 if (!fTracker[iDet]) continue;
2493 fTracker[iDet]->UnloadClusters();
2494 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2495 fLoader[iDet]->UnloadRecPoints();
2496 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2498 // stop filling residuals for TPC and ITS
2499 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2505 //_____________________________________________________________________________
2506 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2508 // Remove the data which are not needed for the physics analysis.
2511 Int_t nTracks=esd->GetNumberOfTracks();
2512 Int_t nV0s=esd->GetNumberOfV0s();
2514 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2516 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2517 Bool_t rc=esd->Clean(cleanPars);
2519 nTracks=esd->GetNumberOfTracks();
2520 nV0s=esd->GetNumberOfV0s();
2522 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2527 //_____________________________________________________________________________
2528 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2530 // fill the event summary data
2532 AliCodeTimerAuto("")
2533 static Int_t eventNr=0;
2534 TString detStr = detectors;
2536 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2537 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2538 AliReconstructor* reconstructor = GetReconstructor(iDet);
2539 if (!reconstructor) continue;
2540 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2541 TTree* clustersTree = NULL;
2542 if (fLoader[iDet]) {
2543 fLoader[iDet]->LoadRecPoints("read");
2544 clustersTree = fLoader[iDet]->TreeR();
2545 if (!clustersTree) {
2546 AliError(Form("Can't get the %s clusters tree",
2547 fgkDetectorName[iDet]));
2548 if (fStopOnError) return kFALSE;
2551 if (fRawReader && !reconstructor->HasDigitConversion()) {
2552 reconstructor->FillESD(fRawReader, clustersTree, esd);
2554 TTree* digitsTree = NULL;
2555 if (fLoader[iDet]) {
2556 fLoader[iDet]->LoadDigits("read");
2557 digitsTree = fLoader[iDet]->TreeD();
2559 AliError(Form("Can't get the %s digits tree",
2560 fgkDetectorName[iDet]));
2561 if (fStopOnError) return kFALSE;
2564 reconstructor->FillESD(digitsTree, clustersTree, esd);
2565 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2567 if (fLoader[iDet]) {
2568 fLoader[iDet]->UnloadRecPoints();
2572 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2573 AliError(Form("the following detectors were not found: %s",
2575 if (fStopOnError) return kFALSE;
2577 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2582 //_____________________________________________________________________________
2583 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2585 // Reads the trigger decision which is
2586 // stored in Trigger.root file and fills
2587 // the corresponding esd entries
2589 AliCodeTimerAuto("")
2591 AliInfo("Filling trigger information into the ESD");
2594 AliCTPRawStream input(fRawReader);
2595 if (!input.Next()) {
2596 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2599 if (esd->GetTriggerMask() != input.GetClassMask())
2600 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2601 input.GetClassMask(),esd->GetTriggerMask()));
2602 if (esd->GetOrbitNumber() != input.GetOrbitID())
2603 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2604 input.GetOrbitID(),esd->GetOrbitNumber()));
2605 if (esd->GetBunchCrossNumber() != input.GetBCID())
2606 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2607 input.GetBCID(),esd->GetBunchCrossNumber()));
2608 AliESDHeader* esdheader = esd->GetHeader();
2609 esdheader->SetL0TriggerInputs(input.GetL0Inputs());
2610 esdheader->SetL1TriggerInputs(input.GetL1Inputs());
2611 esdheader->SetL2TriggerInputs(input.GetL2Inputs());
2613 UInt_t orbit=input.GetOrbitID();
2614 for(Int_t i=0 ; i<input.GetNIRs() ; i++ )
2615 if(TMath::Abs(Int_t(orbit-(input.GetIR(i))->GetOrbit()))<=1){
2616 esdheader->AddTriggerIR(input.GetIR(i));
2622 //_____________________________________________________________________________
2623 Bool_t AliReconstruction::FillTriggerScalers(AliESDEvent*& esd)
2626 //fRunScalers->Print();
2627 if(fRunScalers && fRunScalers->CheckRunScalers()){
2628 AliTimeStamp* timestamp = new AliTimeStamp(esd->GetOrbitNumber(), esd->GetPeriodNumber(), esd->GetBunchCrossNumber());
2629 //AliTimeStamp* timestamp = new AliTimeStamp(10308000, 0, (ULong64_t)486238);
2630 AliESDHeader* esdheader = fesd->GetHeader();
2631 for(Int_t i=0;i<50;i++){
2632 if((1<<i) & esd->GetTriggerMask()){
2633 AliTriggerScalersESD* scalesd = fRunScalers->GetScalersForEventClass( timestamp, i+1);
2634 if(scalesd)esdheader->SetTriggerScalersRecord(scalesd);
2640 //_____________________________________________________________________________
2641 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2644 // Filling information from RawReader Header
2647 if (!fRawReader) return kFALSE;
2649 AliInfo("Filling information from RawReader Header");
2651 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2652 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2653 esd->SetPeriodNumber(fRawReader->GetPeriod());
2655 esd->SetTimeStamp(fRawReader->GetTimestamp());
2656 esd->SetEventType(fRawReader->GetType());
2662 //_____________________________________________________________________________
2663 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2665 // check whether detName is contained in detectors
2666 // if yes, it is removed from detectors
2668 // check if all detectors are selected
2669 if ((detectors.CompareTo("ALL") == 0) ||
2670 detectors.BeginsWith("ALL ") ||
2671 detectors.EndsWith(" ALL") ||
2672 detectors.Contains(" ALL ")) {
2677 // search for the given detector
2678 Bool_t result = kFALSE;
2679 if ((detectors.CompareTo(detName) == 0) ||
2680 detectors.BeginsWith(detName+" ") ||
2681 detectors.EndsWith(" "+detName) ||
2682 detectors.Contains(" "+detName+" ")) {
2683 detectors.ReplaceAll(detName, "");
2687 // clean up the detectors string
2688 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2689 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2690 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2695 //_____________________________________________________________________________
2696 Bool_t AliReconstruction::InitRunLoader()
2698 // get or create the run loader
2700 if (gAlice) delete gAlice;
2703 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2704 // load all base libraries to get the loader classes
2705 TString libs = gSystem->GetLibraries();
2706 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2707 TString detName = fgkDetectorName[iDet];
2708 if (detName == "HLT") continue;
2709 if (libs.Contains("lib" + detName + "base.so")) continue;
2710 gSystem->Load("lib" + detName + "base.so");
2712 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2714 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2719 fRunLoader->CdGAFile();
2720 fRunLoader->LoadgAlice();
2722 //PH This is a temporary fix to give access to the kinematics
2723 //PH that is needed for the labels of ITS clusters
2724 fRunLoader->LoadHeader();
2725 fRunLoader->LoadKinematics();
2727 } else { // galice.root does not exist
2729 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2731 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2732 AliConfig::GetDefaultEventFolderName(),
2735 AliError(Form("could not create run loader in file %s",
2736 fGAliceFileName.Data()));
2740 fIsNewRunLoader = kTRUE;
2741 fRunLoader->MakeTree("E");
2743 if (fNumberOfEventsPerFile > 0)
2744 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2746 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2752 //_____________________________________________________________________________
2753 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2755 // get the reconstructor object and the loader for a detector
2757 if (fReconstructor[iDet]) {
2758 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2759 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2760 fReconstructor[iDet]->SetRecoParam(par);
2761 fReconstructor[iDet]->SetRunInfo(fRunInfo);
2763 return fReconstructor[iDet];
2766 // load the reconstructor object
2767 TPluginManager* pluginManager = gROOT->GetPluginManager();
2768 TString detName = fgkDetectorName[iDet];
2769 TString recName = "Ali" + detName + "Reconstructor";
2771 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2773 AliReconstructor* reconstructor = NULL;
2774 // first check if a plugin is defined for the reconstructor
2775 TPluginHandler* pluginHandler =
2776 pluginManager->FindHandler("AliReconstructor", detName);
2777 // if not, add a plugin for it
2778 if (!pluginHandler) {
2779 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2780 TString libs = gSystem->GetLibraries();
2781 if (libs.Contains("lib" + detName + "base.so") ||
2782 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2783 pluginManager->AddHandler("AliReconstructor", detName,
2784 recName, detName + "rec", recName + "()");
2786 pluginManager->AddHandler("AliReconstructor", detName,
2787 recName, detName, recName + "()");
2789 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2791 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2792 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2794 if (reconstructor) {
2795 TObject* obj = fOptions.FindObject(detName.Data());
2796 if (obj) reconstructor->SetOption(obj->GetTitle());
2797 reconstructor->SetRunInfo(fRunInfo);
2798 reconstructor->Init();
2799 fReconstructor[iDet] = reconstructor;
2802 // get or create the loader
2803 if (detName != "HLT") {
2804 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2805 if (!fLoader[iDet]) {
2806 AliConfig::Instance()
2807 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2809 // first check if a plugin is defined for the loader
2811 pluginManager->FindHandler("AliLoader", detName);
2812 // if not, add a plugin for it
2813 if (!pluginHandler) {
2814 TString loaderName = "Ali" + detName + "Loader";
2815 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2816 pluginManager->AddHandler("AliLoader", detName,
2817 loaderName, detName + "base",
2818 loaderName + "(const char*, TFolder*)");
2819 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2821 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2823 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2824 fRunLoader->GetEventFolder());
2826 if (!fLoader[iDet]) { // use default loader
2827 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2829 if (!fLoader[iDet]) {
2830 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2831 if (fStopOnError) return NULL;
2833 fRunLoader->AddLoader(fLoader[iDet]);
2834 fRunLoader->CdGAFile();
2835 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2836 fRunLoader->Write(0, TObject::kOverwrite);
2841 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2842 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2843 reconstructor->SetRecoParam(par);
2844 reconstructor->SetRunInfo(fRunInfo);
2846 return reconstructor;
2849 //_____________________________________________________________________________
2850 AliVertexer* AliReconstruction::CreateVertexer()
2852 // create the vertexer
2853 // Please note that the caller is the owner of the
2856 AliVertexer* vertexer = NULL;
2857 AliReconstructor* itsReconstructor = GetReconstructor(0);
2858 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) || fRunTracking.Contains("ITS"))) {
2859 vertexer = itsReconstructor->CreateVertexer();
2862 AliWarning("couldn't create a vertexer for ITS");
2868 //_____________________________________________________________________________
2869 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2871 // create the trackers
2872 AliInfo("Creating trackers");
2874 TString detStr = detectors;
2875 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2876 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2877 AliReconstructor* reconstructor = GetReconstructor(iDet);
2878 if (!reconstructor) continue;
2879 TString detName = fgkDetectorName[iDet];
2880 if (detName == "HLT") {
2881 fRunHLTTracking = kTRUE;
2884 if (detName == "MUON") {
2885 fRunMuonTracking = kTRUE;
2890 fTracker[iDet] = reconstructor->CreateTracker();
2891 if (!fTracker[iDet] && (iDet < 7)) {
2892 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2893 if (fStopOnError) return kFALSE;
2895 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2901 //_____________________________________________________________________________
2902 void AliReconstruction::CleanUp()
2904 // delete trackers and the run loader and close and delete the file
2906 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2907 delete fReconstructor[iDet];
2908 fReconstructor[iDet] = NULL;
2909 fLoader[iDet] = NULL;
2910 delete fTracker[iDet];
2911 fTracker[iDet] = NULL;
2916 delete fSPDTrackleter;
2917 fSPDTrackleter = NULL;
2926 delete fParentRawReader;
2927 fParentRawReader=NULL;
2935 if (AliQAManager::QAManager())
2936 AliQAManager::QAManager()->ShowQA() ;
2937 AliQAManager::Destroy() ;
2939 TGeoGlobalMagField::Instance()->SetField(NULL);
2942 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2944 // Write space-points which are then used in the alignment procedures
2945 // For the moment only ITS, TPC, TRD and TOF
2947 Int_t ntracks = esd->GetNumberOfTracks();
2948 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2950 AliESDtrack *track = esd->GetTrack(itrack);
2953 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
2954 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2955 nsp += track->GetNcls(iDet);
2957 if (iDet==0) { // ITS "extra" clusters
2958 track->GetClusters(iDet,idx);
2959 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2964 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2965 track->SetTrackPointArray(sp);
2967 for (Int_t iDet = 5; iDet >= 0; iDet--) {
2968 AliTracker *tracker = fTracker[iDet];
2969 if (!tracker) continue;
2970 Int_t nspdet = track->GetClusters(iDet,idx);
2972 if (iDet==0) // ITS "extra" clusters
2973 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2975 if (nspdet <= 0) continue;
2979 while (isp2 < nspdet) {
2980 Bool_t isvalid=kTRUE;
2982 Int_t index=idx[isp++];
2983 if (index < 0) continue;
2985 TString dets = fgkDetectorName[iDet];
2986 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2987 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2988 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2989 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2990 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2992 isvalid = tracker->GetTrackPoint(index,p);
2995 if (!isvalid) continue;
2996 if (iDet==0 && (isp-1)>=6) p.SetExtra();
2997 sp->AddPoint(isptrack,&p); isptrack++;
3004 //_____________________________________________________________________________
3005 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
3007 // The method reads the raw-data error log
3008 // accumulated within the rawReader.
3009 // It extracts the raw-data errors related to
3010 // the current event and stores them into
3011 // a TClonesArray inside the esd object.
3013 if (!fRawReader) return;
3015 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
3017 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
3019 if (iEvent != log->GetEventNumber()) continue;
3021 esd->AddRawDataErrorLog(log);
3026 //_____________________________________________________________________________
3027 void AliReconstruction::CheckQA()
3029 // check the QA of SIM for this run and remove the detectors
3030 // with status Fatal
3032 // TString newRunLocalReconstruction ;
3033 // TString newRunTracking ;
3034 // TString newFillESD ;
3036 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
3037 // TString detName(AliQAv1::GetDetName(iDet)) ;
3038 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
3039 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
3040 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
3041 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
3043 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
3044 // fRunLocalReconstruction.Contains("ALL") ) {
3045 // newRunLocalReconstruction += detName ;
3046 // newRunLocalReconstruction += " " ;
3048 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
3049 // fRunTracking.Contains("ALL") ) {
3050 // newRunTracking += detName ;
3051 // newRunTracking += " " ;
3053 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
3054 // fFillESD.Contains("ALL") ) {
3055 // newFillESD += detName ;
3056 // newFillESD += " " ;
3060 // fRunLocalReconstruction = newRunLocalReconstruction ;
3061 // fRunTracking = newRunTracking ;
3062 // fFillESD = newFillESD ;
3065 //_____________________________________________________________________________
3066 Int_t AliReconstruction::GetDetIndex(const char* detector)
3068 // return the detector index corresponding to detector
3070 for (index = 0; index < kNDetectors ; index++) {
3071 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3076 //_____________________________________________________________________________
3077 Bool_t AliReconstruction::FinishPlaneEff() {
3079 // Here execute all the necessary operationis, at the end of the tracking phase,
3080 // in case that evaluation of PlaneEfficiencies was required for some detector.
3081 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3083 // This Preliminary version works only FOR ITS !!!!!
3084 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3087 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3090 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3091 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3092 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3093 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
3094 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3095 TString name=planeeff->GetName();
3097 TFile* pefile = TFile::Open(name, "RECREATE");
3098 ret=(Bool_t)planeeff->Write();
3100 if(planeeff->GetCreateHistos()) {
3101 TString hname=planeeff->GetName();
3102 hname+="Histo.root";
3103 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3106 if(fSPDTrackleter) {
3107 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3108 TString name="AliITSPlaneEffSPDtracklet.root";
3109 TFile* pefile = TFile::Open(name, "RECREATE");
3110 ret=(Bool_t)planeeff->Write();
3112 AliESDEvent *dummy=NULL;
3113 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3118 //_____________________________________________________________________________
3119 Bool_t AliReconstruction::InitPlaneEff() {
3121 // Here execute all the necessary operations, before of the tracking phase,
3122 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3123 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3124 // which should be updated/recalculated.
3126 // This Preliminary version will work only FOR ITS !!!!!
3127 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3130 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3132 AliWarning(Form("Implementation of this method not yet completed !! Method return kTRUE"));
3134 fSPDTrackleter = NULL;
3135 AliReconstructor* itsReconstructor = GetReconstructor(0);
3136 if (itsReconstructor) {
3137 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3139 if (fSPDTrackleter) {
3140 AliInfo("Trackleter for SPD has been created");
3146 //_____________________________________________________________________________
3147 Bool_t AliReconstruction::InitAliEVE()
3149 // This method should be called only in case
3150 // AliReconstruction is run
3151 // within the alieve environment.
3152 // It will initialize AliEVE in a way
3153 // so that it can visualize event processed
3154 // by AliReconstruction.
3155 // The return flag shows whenever the
3156 // AliEVE initialization was successful or not.
3159 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3160 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3161 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3163 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3164 gROOT->ProcessLine("alieve_online_init()");
3169 //_____________________________________________________________________________
3170 void AliReconstruction::RunAliEVE()
3172 // Runs AliEVE visualisation of
3173 // the current event.
3174 // Should be executed only after
3175 // successful initialization of AliEVE.
3177 AliInfo("Running AliEVE...");
3178 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3182 //_____________________________________________________________________________
3183 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3185 // Allows to run QA for a selected set of detectors
3186 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3187 // all selected detectors run the same selected tasks
3189 if (!detAndAction.Contains(":")) {
3190 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3194 Int_t colon = detAndAction.Index(":") ;
3195 fQADetectors = detAndAction(0, colon) ;
3196 if (fQADetectors.Contains("ALL") )
3197 fQADetectors = fFillESD ;
3198 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3199 if (fQATasks.Contains("ALL") ) {
3200 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3202 fQATasks.ToUpper() ;
3204 if ( fQATasks.Contains("RAW") )
3205 tempo = Form("%d ", AliQAv1::kRAWS) ;
3206 if ( fQATasks.Contains("DIGIT") )
3207 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3208 if ( fQATasks.Contains("RECPOINT") )
3209 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3210 if ( fQATasks.Contains("ESD") )
3211 tempo += Form("%d ", AliQAv1::kESDS) ;
3213 if (fQATasks.IsNull()) {
3214 AliInfo("No QA requested\n") ;
3219 TString tempo(fQATasks) ;
3220 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3221 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3222 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3223 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3224 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3229 //_____________________________________________________________________________
3230 Bool_t AliReconstruction::InitRecoParams()
3232 // The method accesses OCDB and retrieves all
3233 // the available reco-param objects from there.
3235 Bool_t isOK = kTRUE;
3237 if (fRecoParam.GetDetRecoParamArray(kNDetectors)) {
3238 AliInfo("Using custom GRP reconstruction parameters");
3241 AliInfo("Loading GRP reconstruction parameter objects");
3243 AliCDBPath path("GRP","Calib","RecoParam");
3244 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3246 AliWarning("Couldn't find GRP RecoParam entry in OCDB");
3250 TObject *recoParamObj = entry->GetObject();
3251 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3252 // GRP has a normal TobjArray of AliDetectorRecoParam objects
3253 // Registering them in AliRecoParam
3254 fRecoParam.AddDetRecoParamArray(kNDetectors,dynamic_cast<TObjArray*>(recoParamObj));
3256 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3257 // GRP has only onse set of reco parameters
3258 // Registering it in AliRecoParam
3259 AliInfo("Single set of GRP reconstruction parameters found");
3260 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3261 fRecoParam.AddDetRecoParam(kNDetectors,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3264 AliError("No valid GRP RecoParam object found in the OCDB");
3271 TString detStr = fLoadCDB;
3272 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3274 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3276 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3277 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3281 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3283 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3284 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3286 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3290 TObject *recoParamObj = entry->GetObject();
3291 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3292 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3293 // Registering them in AliRecoParam
3294 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3296 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3297 // The detector has only onse set of reco parameters
3298 // Registering it in AliRecoParam
3299 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3300 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3301 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3304 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3308 // FIX ME: We have to disable the unloading of reco-param CDB
3309 // entries because QA framework is using them. Has to be fix in
3310 // a way that the QA takes the objects already constructed in
3312 // AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3316 if (AliDebugLevel() > 0) fRecoParam.Print();
3321 //_____________________________________________________________________________
3322 Bool_t AliReconstruction::GetEventInfo()
3324 // Fill the event info object
3326 AliCodeTimerAuto("")
3328 AliCentralTrigger *aCTP = NULL;
3330 fEventInfo.SetEventType(fRawReader->GetType());
3332 ULong64_t mask = fRawReader->GetClassMask();
3333 fEventInfo.SetTriggerMask(mask);
3334 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3335 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3337 aCTP = new AliCentralTrigger();
3338 TString configstr("");
3339 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3340 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3344 aCTP->SetClassMask(mask);
3345 aCTP->SetClusterMask(clmask);
3348 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3350 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3351 aCTP = fRunLoader->GetTrigger();
3352 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3353 // get inputs from actp - just get
3354 AliESDHeader* esdheader = fesd->GetHeader();
3355 esdheader->SetL0TriggerInputs(aCTP->GetL0TriggerInputs());
3356 esdheader->SetL1TriggerInputs(aCTP->GetL1TriggerInputs());
3357 esdheader->SetL2TriggerInputs(aCTP->GetL2TriggerInputs());
3358 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3361 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3366 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3368 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3369 if (fRawReader) delete aCTP;
3373 UChar_t clustmask = 0;
3375 ULong64_t trmask = fEventInfo.GetTriggerMask();
3376 const TObjArray& classesArray = config->GetClasses();
3377 Int_t nclasses = classesArray.GetEntriesFast();
3378 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3379 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3381 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3382 fesd->SetTriggerClass(trclass->GetName(),trindex);
3383 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3384 if (trmask & (1ull << trindex)) {
3386 trclasses += trclass->GetName();
3388 clustmask |= trclass->GetCluster()->GetClusterMask();
3392 fEventInfo.SetTriggerClasses(trclasses);
3394 // Set the information in ESD
3395 fesd->SetTriggerMask(trmask);
3396 fesd->SetTriggerCluster(clustmask);
3398 if (!aCTP->CheckTriggeredDetectors()) {
3399 if (fRawReader) delete aCTP;
3403 if (fRawReader) delete aCTP;
3405 // We have to fill also the HLT decision here!!
3411 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3413 // Match the detector list found in the rec.C or the default 'ALL'
3414 // to the list found in the GRP (stored there by the shuttle PP which
3415 // gets the information from ECS)
3416 static TString resultList;
3417 TString detList = detectorList;
3421 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3422 if ((detectorMask >> iDet) & 0x1) {
3423 TString det = AliDAQ::OfflineModuleName(iDet);
3424 if ((detList.CompareTo("ALL") == 0) ||
3425 ((detList.BeginsWith("ALL ") ||
3426 detList.EndsWith(" ALL") ||
3427 detList.Contains(" ALL ")) &&
3428 !(detList.BeginsWith("-"+det+" ") ||
3429 detList.EndsWith(" -"+det) ||
3430 detList.Contains(" -"+det+" "))) ||
3431 (detList.CompareTo(det) == 0) ||
3432 detList.BeginsWith(det+" ") ||
3433 detList.EndsWith(" "+det) ||
3434 detList.Contains( " "+det+" " )) {
3435 if (!resultList.EndsWith(det + " ")) {
3444 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3445 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3446 if ((detList.CompareTo("ALL") == 0) ||
3447 ((detList.BeginsWith("ALL ") ||
3448 detList.EndsWith(" ALL") ||
3449 detList.Contains(" ALL ")) &&
3450 !(detList.BeginsWith("-"+hltDet+" ") ||
3451 detList.EndsWith(" -"+hltDet) ||
3452 detList.Contains(" -"+hltDet+" "))) ||
3453 (detList.CompareTo(hltDet) == 0) ||
3454 detList.BeginsWith(hltDet+" ") ||
3455 detList.EndsWith(" "+hltDet) ||
3456 detList.Contains( " "+hltDet+" " )) {
3457 resultList += hltDet;
3461 return resultList.Data();
3465 //______________________________________________________________________________
3466 void AliReconstruction::Abort(const char *method, EAbort what)
3468 // Abort processing. If what = kAbortProcess, the Process() loop will be
3469 // aborted. If what = kAbortFile, the current file in a chain will be
3470 // aborted and the processing will continue with the next file, if there
3471 // is no next file then Process() will be aborted. Abort() can also be
3472 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3473 // the SlaveTerminate() and Terminate() are always called. The abort flag
3474 // can be checked in these methods using GetAbort().
3476 // The method is overwritten in AliReconstruction for better handling of
3477 // reco specific errors
3479 if (!fStopOnError) return;
3483 TString whyMess = method;
3484 whyMess += " failed! Aborting...";
3486 AliError(whyMess.Data());
3489 TString mess = "Abort";
3490 if (fAbort == kAbortProcess)
3491 mess = "AbortProcess";
3492 else if (fAbort == kAbortFile)
3495 Info(mess, whyMess.Data());
3498 //______________________________________________________________________________
3499 Bool_t AliReconstruction::ProcessEvent(void* event)
3501 // Method that is used in case the event loop
3502 // is steered from outside, for example by AMORE
3503 // 'event' is a pointer to the DATE event in the memory
3505 if (fRawReader) delete fRawReader;
3506 fRawReader = new AliRawReaderDate(event);
3507 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());