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>
130 #include "AliAlignObj.h"
131 #include "AliCDBEntry.h"
132 #include "AliCDBManager.h"
133 #include "AliCDBStorage.h"
134 #include "AliCTPRawStream.h"
135 #include "AliCascadeVertexer.h"
136 #include "AliCentralTrigger.h"
137 #include "AliCodeTimer.h"
139 #include "AliDetectorRecoParam.h"
140 #include "AliESDCaloCells.h"
141 #include "AliESDCaloCluster.h"
142 #include "AliESDEvent.h"
143 #include "AliESDMuonTrack.h"
144 #include "AliESDPmdTrack.h"
145 #include "AliESDTagCreator.h"
146 #include "AliESDVertex.h"
147 #include "AliESDcascade.h"
148 #include "AliESDfriend.h"
149 #include "AliESDkink.h"
150 #include "AliESDpid.h"
151 #include "AliESDtrack.h"
152 #include "AliESDtrack.h"
153 #include "AliEventInfo.h"
154 #include "AliGRPObject.h"
155 #include "AliGRPRecoParam.h"
156 #include "AliGenEventHeader.h"
157 #include "AliGeomManager.h"
158 #include "AliGlobalQADataMaker.h"
159 #include "AliHeader.h"
162 #include "AliMultiplicity.h"
164 #include "AliPlaneEff.h"
166 #include "AliQADataMakerRec.h"
167 #include "AliQAManager.h"
168 #include "AliRawVEvent.h"
169 #include "AliRawEventHeaderBase.h"
170 #include "AliRawHLTManager.h"
171 #include "AliRawReaderDate.h"
172 #include "AliRawReaderFile.h"
173 #include "AliRawReaderRoot.h"
174 #include "AliReconstruction.h"
175 #include "AliReconstructor.h"
177 #include "AliRunInfo.h"
178 #include "AliRunLoader.h"
179 #include "AliSysInfo.h" // memory snapshots
180 #include "AliTrackPointArray.h"
181 #include "AliTracker.h"
182 #include "AliTriggerClass.h"
183 #include "AliTriggerCluster.h"
184 #include "AliTriggerConfiguration.h"
185 #include "AliV0vertexer.h"
186 #include "AliVertexer.h"
187 #include "AliVertexerTracks.h"
189 ClassImp(AliReconstruction)
191 //_____________________________________________________________________________
192 const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
194 //_____________________________________________________________________________
195 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
197 fUniformField(kFALSE),
198 fRunVertexFinder(kTRUE),
199 fRunVertexFinderTracks(kTRUE),
200 fRunHLTTracking(kFALSE),
201 fRunMuonTracking(kFALSE),
203 fRunCascadeFinder(kTRUE),
204 fStopOnError(kFALSE),
205 fWriteAlignmentData(kFALSE),
206 fWriteESDfriend(kFALSE),
207 fFillTriggerESD(kTRUE),
215 fRunLocalReconstruction("ALL"),
219 fUseTrackingErrorsForAlignment(""),
220 fGAliceFileName(gAliceFilename),
226 fNumberOfEventsPerFile((UInt_t)-1),
228 fLoadAlignFromCDB(kTRUE),
229 fLoadAlignData("ALL"),
236 fParentRawReader(NULL),
240 fSPDTrackleter(NULL),
242 fDiamondProfileSPD(NULL),
243 fDiamondProfile(NULL),
244 fDiamondProfileTPC(NULL),
248 fAlignObjArray(NULL),
252 fInitCDBCalled(kFALSE),
253 fSetRunNumberFromDataCalled(kFALSE),
259 fSameQACycle(kFALSE),
260 fInitQACalled(kFALSE),
261 fWriteQAExpertData(kTRUE),
262 fRunPlaneEff(kFALSE),
271 fIsNewRunLoader(kFALSE),
275 // create reconstruction object with default parameters
278 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
279 fReconstructor[iDet] = NULL;
280 fLoader[iDet] = NULL;
281 fTracker[iDet] = NULL;
283 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
284 fQACycles[iDet] = 999999 ;
285 fQAWriteExpert[iDet] = kFALSE ;
291 //_____________________________________________________________________________
292 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
294 fUniformField(rec.fUniformField),
295 fRunVertexFinder(rec.fRunVertexFinder),
296 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
297 fRunHLTTracking(rec.fRunHLTTracking),
298 fRunMuonTracking(rec.fRunMuonTracking),
299 fRunV0Finder(rec.fRunV0Finder),
300 fRunCascadeFinder(rec.fRunCascadeFinder),
301 fStopOnError(rec.fStopOnError),
302 fWriteAlignmentData(rec.fWriteAlignmentData),
303 fWriteESDfriend(rec.fWriteESDfriend),
304 fFillTriggerESD(rec.fFillTriggerESD),
306 fCleanESD(rec.fCleanESD),
307 fV0DCAmax(rec.fV0DCAmax),
308 fV0CsPmin(rec.fV0CsPmin),
312 fRunLocalReconstruction(rec.fRunLocalReconstruction),
313 fRunTracking(rec.fRunTracking),
314 fFillESD(rec.fFillESD),
315 fLoadCDB(rec.fLoadCDB),
316 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
317 fGAliceFileName(rec.fGAliceFileName),
318 fRawInput(rec.fRawInput),
319 fESDOutput(rec.fESDOutput),
320 fEquipIdMap(rec.fEquipIdMap),
321 fFirstEvent(rec.fFirstEvent),
322 fLastEvent(rec.fLastEvent),
323 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
325 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
326 fLoadAlignData(rec.fLoadAlignData),
327 fUseHLTData(rec.fUseHLTData),
333 fParentRawReader(NULL),
335 fRecoParam(rec.fRecoParam),
337 fSPDTrackleter(NULL),
339 fDiamondProfileSPD(rec.fDiamondProfileSPD),
340 fDiamondProfile(rec.fDiamondProfile),
341 fDiamondProfileTPC(rec.fDiamondProfileTPC),
345 fAlignObjArray(rec.fAlignObjArray),
346 fCDBUri(rec.fCDBUri),
347 fQARefUri(rec.fQARefUri),
349 fInitCDBCalled(rec.fInitCDBCalled),
350 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
351 fQADetectors(rec.fQADetectors),
353 fQATasks(rec.fQATasks),
355 fRunGlobalQA(rec.fRunGlobalQA),
356 fSameQACycle(rec.fSameQACycle),
357 fInitQACalled(rec.fInitQACalled),
358 fWriteQAExpertData(rec.fWriteQAExpertData),
359 fRunPlaneEff(rec.fRunPlaneEff),
368 fIsNewRunLoader(rec.fIsNewRunLoader),
374 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
375 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
377 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
378 fReconstructor[iDet] = NULL;
379 fLoader[iDet] = NULL;
380 fTracker[iDet] = NULL;
383 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
384 fQACycles[iDet] = rec.fQACycles[iDet];
385 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
388 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
389 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
393 //_____________________________________________________________________________
394 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
396 // assignment operator
397 // Used in PROOF mode
398 // Be very careful while modifing it!
399 // Simple rules to follow:
400 // for persistent data members - use their assignment operators
401 // for non-persistent ones - do nothing or take the default values from constructor
402 // TSelector members should not be touched
403 if(&rec == this) return *this;
405 fUniformField = rec.fUniformField;
406 fRunVertexFinder = rec.fRunVertexFinder;
407 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
408 fRunHLTTracking = rec.fRunHLTTracking;
409 fRunMuonTracking = rec.fRunMuonTracking;
410 fRunV0Finder = rec.fRunV0Finder;
411 fRunCascadeFinder = rec.fRunCascadeFinder;
412 fStopOnError = rec.fStopOnError;
413 fWriteAlignmentData = rec.fWriteAlignmentData;
414 fWriteESDfriend = rec.fWriteESDfriend;
415 fFillTriggerESD = rec.fFillTriggerESD;
417 fCleanESD = rec.fCleanESD;
418 fV0DCAmax = rec.fV0DCAmax;
419 fV0CsPmin = rec.fV0CsPmin;
423 fRunLocalReconstruction = rec.fRunLocalReconstruction;
424 fRunTracking = rec.fRunTracking;
425 fFillESD = rec.fFillESD;
426 fLoadCDB = rec.fLoadCDB;
427 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
428 fGAliceFileName = rec.fGAliceFileName;
429 fRawInput = rec.fRawInput;
430 fESDOutput = rec.fESDOutput;
431 fEquipIdMap = rec.fEquipIdMap;
432 fFirstEvent = rec.fFirstEvent;
433 fLastEvent = rec.fLastEvent;
434 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
436 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
437 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
440 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
441 fLoadAlignData = rec.fLoadAlignData;
442 fUseHLTData = rec.fUseHLTData;
444 delete fRunInfo; fRunInfo = NULL;
445 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
447 fEventInfo = rec.fEventInfo;
451 fParentRawReader = NULL;
453 fRecoParam = rec.fRecoParam;
455 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
456 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
457 delete fLoader[iDet]; fLoader[iDet] = NULL;
458 delete fTracker[iDet]; fTracker[iDet] = NULL;
461 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
462 fQACycles[iDet] = rec.fQACycles[iDet];
463 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
466 delete fSPDTrackleter; fSPDTrackleter = NULL;
468 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
469 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
470 delete fDiamondProfile; fDiamondProfile = NULL;
471 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
472 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
473 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
475 delete fGRPData; fGRPData = NULL;
476 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
477 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
479 delete fAlignObjArray; fAlignObjArray = NULL;
482 fQARefUri = rec.fQARefUri;
483 fSpecCDBUri.Delete();
484 fInitCDBCalled = rec.fInitCDBCalled;
485 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
486 fQADetectors = rec.fQADetectors;
488 fQATasks = rec.fQATasks;
490 fRunGlobalQA = rec.fRunGlobalQA;
491 fSameQACycle = rec.fSameQACycle;
492 fInitQACalled = rec.fInitQACalled;
493 fWriteQAExpertData = rec.fWriteQAExpertData;
494 fRunPlaneEff = rec.fRunPlaneEff;
503 fIsNewRunLoader = rec.fIsNewRunLoader;
510 //_____________________________________________________________________________
511 AliReconstruction::~AliReconstruction()
518 if (fAlignObjArray) {
519 fAlignObjArray->Delete();
520 delete fAlignObjArray;
522 fSpecCDBUri.Delete();
524 AliCodeTimer::Instance()->Print();
527 //_____________________________________________________________________________
528 void AliReconstruction::InitQA()
530 //Initialize the QA and start of cycle
531 AliCodeTimerAuto("");
533 if (fInitQACalled) return;
534 fInitQACalled = kTRUE;
536 fQAManager = AliQAManager::QAManager("rec") ;
537 if (fWriteQAExpertData)
538 fQAManager->SetWriteExpert() ;
540 if (fQAManager->IsDefaultStorageSet()) {
541 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
542 AliWarning("Default QA reference storage has been already set !");
543 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fQARefUri.Data()));
544 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
545 fQARefUri = fQAManager->GetDefaultStorage()->GetURI();
547 if (fQARefUri.Length() > 0) {
548 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
549 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
550 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
552 fQARefUri="local://$ALICE_ROOT/QAref";
553 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
554 AliWarning("Default QA refeference storage not yet set !!!!");
555 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
556 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
559 fQAManager->SetDefaultStorage(fQARefUri);
563 fQAManager->SetActiveDetectors(fQADetectors) ;
564 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) {
565 fQAManager->SetCycleLength(AliQAv1::DETECTORINDEX_t(det), fQACycles[det]) ;
566 fQAManager->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
568 if (!fRawReader && !fInput && fQATasks.Contains(AliQAv1::kRAWS))
569 fQATasks.ReplaceAll(Form("%d",AliQAv1::kRAWS), "") ;
570 fQAManager->SetTasks(fQATasks) ;
571 fQAManager->InitQADataMaker(AliCDBManager::Instance()->GetRun()) ;
574 Bool_t sameCycle = kFALSE ;
575 AliQADataMaker *qadm = fQAManager->GetQADataMaker(AliQAv1::kGLOBAL);
576 AliInfo(Form("Initializing the global QA data maker"));
577 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS))) {
578 qadm->StartOfCycle(AliQAv1::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
579 TObjArray **arr=qadm->Init(AliQAv1::kRECPOINTS);
580 AliTracker::SetResidualsArray(arr);
583 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS))) {
584 qadm->StartOfCycle(AliQAv1::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
585 qadm->Init(AliQAv1::kESDS);
588 AliSysInfo::AddStamp("InitQA") ;
591 //_____________________________________________________________________________
592 void AliReconstruction::MergeQA(const char *fileName)
594 //Initialize the QA and start of cycle
595 AliCodeTimerAuto("") ;
596 if ( ! fQAManager ) {
597 AliFatal("Hum... this should not happen") ;
599 fQAManager->Merge(AliCDBManager::Instance()->GetRun(),fileName) ;
601 AliSysInfo::AddStamp("MergeQA") ;
604 //_____________________________________________________________________________
605 void AliReconstruction::InitCDB()
607 // activate a default CDB storage
608 // First check if we have any CDB storage set, because it is used
609 // to retrieve the calibration and alignment constants
610 AliCodeTimerAuto("");
612 if (fInitCDBCalled) return;
613 fInitCDBCalled = kTRUE;
615 AliCDBManager* man = AliCDBManager::Instance();
616 if (man->IsDefaultStorageSet())
618 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
619 AliWarning("Default CDB storage has been already set !");
620 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
621 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
622 fCDBUri = man->GetDefaultStorage()->GetURI();
625 if (fCDBUri.Length() > 0)
627 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
628 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
629 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
631 fCDBUri="local://$ALICE_ROOT/OCDB";
632 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
633 AliWarning("Default CDB storage not yet set !!!!");
634 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
635 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
638 man->SetDefaultStorage(fCDBUri);
641 // Now activate the detector specific CDB storage locations
642 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
643 TObject* obj = fSpecCDBUri[i];
645 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
646 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
647 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
648 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
650 AliSysInfo::AddStamp("InitCDB");
653 //_____________________________________________________________________________
654 void AliReconstruction::SetDefaultStorage(const char* uri) {
655 // Store the desired default CDB storage location
656 // Activate it later within the Run() method
662 //_____________________________________________________________________________
663 void AliReconstruction::SetQARefDefaultStorage(const char* uri) {
664 // Store the desired default CDB storage location
665 // Activate it later within the Run() method
668 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
671 //_____________________________________________________________________________
672 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
673 // Store a detector-specific CDB storage location
674 // Activate it later within the Run() method
676 AliCDBPath aPath(calibType);
677 if(!aPath.IsValid()){
678 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
679 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
680 if(!strcmp(calibType, fgkDetectorName[iDet])) {
681 aPath.SetPath(Form("%s/*", calibType));
682 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
686 if(!aPath.IsValid()){
687 AliError(Form("Not a valid path or detector: %s", calibType));
692 // // check that calibType refers to a "valid" detector name
693 // Bool_t isDetector = kFALSE;
694 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
695 // TString detName = fgkDetectorName[iDet];
696 // if(aPath.GetLevel0() == detName) {
697 // isDetector = kTRUE;
703 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
707 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
708 if (obj) fSpecCDBUri.Remove(obj);
709 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
713 //_____________________________________________________________________________
714 Bool_t AliReconstruction::SetRunNumberFromData()
716 // The method is called in Run() in order
717 // to set a correct run number.
718 // In case of raw data reconstruction the
719 // run number is taken from the raw data header
721 if (fSetRunNumberFromDataCalled) return kTRUE;
722 fSetRunNumberFromDataCalled = kTRUE;
724 AliCDBManager* man = AliCDBManager::Instance();
727 if(fRawReader->NextEvent()) {
728 if(man->GetRun() > 0) {
729 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
731 man->SetRun(fRawReader->GetRunNumber());
732 fRawReader->RewindEvents();
735 if(man->GetRun() > 0) {
736 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
739 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
745 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
747 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
752 // read run number from gAlice
753 if(rl->GetHeader()) {
754 man->SetRun(rl->GetHeader()->GetRun());
759 AliError("Neither run-loader header nor RawReader objects are found !");
771 //_____________________________________________________________________________
772 void AliReconstruction::SetCDBLock() {
773 // Set CDB lock: from now on it is forbidden to reset the run number
774 // or the default storage or to activate any further storage!
776 AliCDBManager::Instance()->SetLock(1);
779 //_____________________________________________________________________________
780 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
782 // Read the alignment objects from CDB.
783 // Each detector is supposed to have the
784 // alignment objects in DET/Align/Data CDB path.
785 // All the detector objects are then collected,
786 // sorted by geometry level (starting from ALIC) and
787 // then applied to the TGeo geometry.
788 // Finally an overlaps check is performed.
790 // Load alignment data from CDB and fill fAlignObjArray
791 if(fLoadAlignFromCDB){
793 TString detStr = detectors;
794 TString loadAlObjsListOfDets = "";
796 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
797 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
798 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
800 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
802 loadAlObjsListOfDets += fgkDetectorName[iDet];
803 loadAlObjsListOfDets += " ";
805 } // end loop over detectors
807 if(AliGeomManager::GetNalignable("GRP") != 0)
808 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
809 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
810 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
812 // Check if the array with alignment objects was
813 // provided by the user. If yes, apply the objects
814 // to the present TGeo geometry
815 if (fAlignObjArray) {
816 if (gGeoManager && gGeoManager->IsClosed()) {
817 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
818 AliError("The misalignment of one or more volumes failed!"
819 "Compare the list of simulated detectors and the list of detector alignment data!");
824 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
830 if (fAlignObjArray) {
831 fAlignObjArray->Delete();
832 delete fAlignObjArray; fAlignObjArray=NULL;
838 //_____________________________________________________________________________
839 void AliReconstruction::SetGAliceFile(const char* fileName)
841 // set the name of the galice file
843 fGAliceFileName = fileName;
846 //_____________________________________________________________________________
847 void AliReconstruction::SetInput(const char* input)
849 // In case the input string starts with 'mem://', we run in an online mode
850 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
851 // file is assumed. One can give as an input:
852 // mem://: - events taken from DAQ monitoring libs online
854 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
855 if (input) fRawInput = input;
858 //_____________________________________________________________________________
859 void AliReconstruction::SetOutput(const char* output)
861 // Set the output ESD filename
862 // 'output' is a normalt ROOT url
863 // The method is used in case of raw-data reco with PROOF
864 if (output) fESDOutput.SetUrl(output);
867 //_____________________________________________________________________________
868 void AliReconstruction::SetOption(const char* detector, const char* option)
870 // set options for the reconstruction of a detector
872 TObject* obj = fOptions.FindObject(detector);
873 if (obj) fOptions.Remove(obj);
874 fOptions.Add(new TNamed(detector, option));
877 //_____________________________________________________________________________
878 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
880 // Set custom reconstruction parameters for a given detector
881 // Single set of parameters for all the events
883 // First check if the reco-params are global
884 if(!strcmp(detector, "GRP")) {
886 fRecoParam.AddDetRecoParam(kNDetectors,par);
890 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
891 if(!strcmp(detector, fgkDetectorName[iDet])) {
893 fRecoParam.AddDetRecoParam(iDet,par);
900 //_____________________________________________________________________________
901 Bool_t AliReconstruction::SetFieldMap(Float_t l3Cur, Float_t diCur, Float_t l3Pol,
902 Float_t diPol, Float_t beamenergy,
903 const Char_t *beamtype, const Char_t *path)
905 //------------------------------------------------
906 // The magnetic field map, defined externally...
907 // L3 current 30000 A -> 0.5 T
908 // L3 current 12000 A -> 0.2 T
909 // dipole current 6000 A
910 // The polarities must be the same
911 //------------------------------------------------
912 const Float_t l3NominalCurrent1=30000.; // (A)
913 const Float_t l3NominalCurrent2=12000.; // (A)
914 const Float_t diNominalCurrent =6000. ; // (A)
916 const Float_t tolerance=0.03; // relative current tolerance
917 const Float_t zero=77.; // "zero" current (A)
919 TString s=(l3Pol < 0) ? "L3: -" : "L3: +";
921 AliMagF::BMap_t map = AliMagF::k5kG;
925 l3Cur = TMath::Abs(l3Cur);
926 if (TMath::Abs(l3Cur-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
927 fcL3 = l3Cur/l3NominalCurrent1;
930 } else if (TMath::Abs(l3Cur-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
931 fcL3 = l3Cur/l3NominalCurrent2;
934 } else if (l3Cur <= zero) {
936 map = AliMagF::k5kGUniform;
938 fUniformField=kTRUE; // track with the uniform (zero) B field
940 AliError(Form("Wrong L3 current (%f A)!",l3Cur));
944 diCur = TMath::Abs(diCur);
945 if (TMath::Abs(diCur-diNominalCurrent)/diNominalCurrent < tolerance) {
946 // 3% current tolerance...
947 fcDip = diCur/diNominalCurrent;
949 } else if (diCur <= zero) { // some small current..
953 AliError(Form("Wrong dipole current (%f A)!",diCur));
957 if (l3Pol!=diPol && (map==AliMagF::k5kG || map==AliMagF::k2kG) && fcDip!=0) {
958 AliError("L3 and Dipole polarities must be the same");
962 if (l3Pol<0) fcL3 = -fcL3;
963 if (diPol<0) fcDip = -fcDip;
965 AliMagF::BeamType_t btype = AliMagF::kNoBeamField;
966 TString btypestr = beamtype;
968 TPRegexp protonBeam("(proton|p)\\s*-?\\s*\\1");
969 TPRegexp ionBeam("(lead|pb|ion|a)\\s*-?\\s*\\1");
970 if (btypestr.Contains(ionBeam)) btype = AliMagF::kBeamTypeAA;
971 else if (btypestr.Contains(protonBeam)) btype = AliMagF::kBeamTypepp;
973 AliInfo(Form("Cannot determine the beam type from %s, assume no LHC magnet field",beamtype));
976 AliMagF* fld = new AliMagF("MagneticFieldMap", s.Data(), 2, fcL3, fcDip, 10., map, path,
978 TGeoGlobalMagField::Instance()->SetField( fld );
979 TGeoGlobalMagField::Instance()->Lock();
985 Bool_t AliReconstruction::InitGRP() {
986 //------------------------------------
987 // Initialization of the GRP entry
988 //------------------------------------
989 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
993 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
996 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
998 fGRPData = new AliGRPObject();
999 fGRPData->ReadValuesFromMap(m);
1003 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
1004 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
1008 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
1012 AliError("No GRP entry found in OCDB!");
1016 TString lhcState = fGRPData->GetLHCState();
1017 if (lhcState==AliGRPObject::GetInvalidString()) {
1018 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
1019 lhcState = "UNKNOWN";
1022 TString beamType = fGRPData->GetBeamType();
1023 if (beamType==AliGRPObject::GetInvalidString()) {
1024 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
1025 beamType = "UNKNOWN";
1028 Float_t beamEnergy = fGRPData->GetBeamEnergy();
1029 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
1030 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
1033 // energy is provided in MeV*120
1034 beamEnergy /= 120E3;
1036 TString runType = fGRPData->GetRunType();
1037 if (runType==AliGRPObject::GetInvalidString()) {
1038 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
1039 runType = "UNKNOWN";
1042 Int_t activeDetectors = fGRPData->GetDetectorMask();
1043 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
1044 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
1045 activeDetectors = 1074790399;
1048 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1052 // Process the list of active detectors
1053 if (activeDetectors) {
1054 UInt_t detMask = activeDetectors;
1055 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1056 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1057 fFillESD = MatchDetectorList(fFillESD,detMask);
1058 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1059 fLoadCDB.Form("%s %s %s %s",
1060 fRunLocalReconstruction.Data(),
1061 fRunTracking.Data(),
1063 fQADetectors.Data());
1064 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1065 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1)) {
1066 // switch off the vertexer
1067 AliInfo("SPD is not in the list of active detectors. Vertexer switched off.");
1068 fRunVertexFinder = kFALSE;
1070 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1071 // switch off the reading of CTP raw-data payload
1072 if (fFillTriggerESD) {
1073 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1074 fFillTriggerESD = kFALSE;
1079 AliInfo("===================================================================================");
1080 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1081 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1082 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1083 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1084 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1085 AliInfo("===================================================================================");
1087 //*** Dealing with the magnetic field map
1088 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {AliInfo("Running with the externally locked B field !");}
1090 // Construct the field map out of the information retrieved from GRP.
1093 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1094 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1095 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1099 Char_t l3Polarity = fGRPData->GetL3Polarity();
1100 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1101 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1106 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1107 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1108 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1112 Char_t diPolarity = fGRPData->GetDipolePolarity();
1113 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1114 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1119 TObjString *l3Current=
1120 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
1122 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1125 TObjString *l3Polarity=
1126 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
1128 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1133 TObjString *diCurrent=
1134 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
1136 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1139 TObjString *diPolarity=
1140 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
1142 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1148 if ( !SetFieldMap(l3Current, diCurrent, l3Polarity ? -1:1, diPolarity ? -1:1) )
1149 AliFatal("Failed to creat a B field map ! Exiting...");
1150 AliInfo("Running with the B field constructed out of GRP !");
1152 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1156 //*** Get the diamond profiles from OCDB
1157 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1159 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1161 AliError("No SPD diamond profile found in OCDB!");
1164 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1166 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1168 AliError("No diamond profile found in OCDB!");
1171 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1173 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1175 AliError("No TPC diamond profile found in OCDB!");
1181 //_____________________________________________________________________________
1182 Bool_t AliReconstruction::LoadCDB()
1184 AliCodeTimerAuto("");
1186 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1188 TString detStr = fLoadCDB;
1189 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1190 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1191 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1196 //_____________________________________________________________________________
1197 Bool_t AliReconstruction::Run(const char* input)
1200 AliCodeTimerAuto("");
1203 if (GetAbort() != TSelector::kContinue) return kFALSE;
1205 TChain *chain = NULL;
1206 if (fRawReader && (chain = fRawReader->GetChain())) {
1209 gProof->AddInput(this);
1210 if (!fESDOutput.IsValid()) {
1211 fESDOutput.SetProtocol("root",kTRUE);
1212 fESDOutput.SetHost(gSystem->HostName());
1213 fESDOutput.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1215 AliInfo(Form("Output file with ESDs is %s",fESDOutput.GetUrl()));
1216 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",fESDOutput.GetUrl()));
1217 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1219 chain->Process("AliReconstruction");
1222 chain->Process(this);
1227 if (GetAbort() != TSelector::kContinue) return kFALSE;
1229 if (GetAbort() != TSelector::kContinue) return kFALSE;
1230 //******* The loop over events
1231 AliInfo("Starting looping over events");
1233 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1234 (fRawReader && fRawReader->NextEvent())) {
1235 if (!ProcessEvent(iEvent)) {
1236 Abort("ProcessEvent",TSelector::kAbortFile);
1242 if (GetAbort() != TSelector::kContinue) return kFALSE;
1244 if (GetAbort() != TSelector::kContinue) return kFALSE;
1250 //_____________________________________________________________________________
1251 void AliReconstruction::InitRawReader(const char* input)
1253 AliCodeTimerAuto("");
1255 // Init raw-reader and
1256 // set the input in case of raw data
1257 if (input) fRawInput = input;
1258 fRawReader = AliRawReader::Create(fRawInput.Data());
1260 AliInfo("Reconstruction will run over digits");
1262 if (!fEquipIdMap.IsNull() && fRawReader)
1263 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1265 if (!fUseHLTData.IsNull()) {
1266 // create the RawReaderHLT which performs redirection of HLT input data for
1267 // the specified detectors
1268 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1270 fParentRawReader=fRawReader;
1271 fRawReader=pRawReader;
1273 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1276 AliSysInfo::AddStamp("CreateRawReader");
1279 //_____________________________________________________________________________
1280 void AliReconstruction::InitRun(const char* input)
1282 // Initialization of raw-reader,
1283 // run number, CDB etc.
1284 AliCodeTimerAuto("");
1285 AliSysInfo::AddStamp("Start");
1287 // Initialize raw-reader if any
1288 InitRawReader(input);
1290 // Initialize the CDB storage
1293 // Set run number in CDBManager (if it is not already set by the user)
1294 if (!SetRunNumberFromData()) {
1295 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1299 // Set CDB lock: from now on it is forbidden to reset the run number
1300 // or the default storage or to activate any further storage!
1305 //_____________________________________________________________________________
1306 void AliReconstruction::Begin(TTree *)
1308 // Initialize AlReconstruction before
1309 // going into the event loop
1310 // Should follow the TSelector convention
1311 // i.e. initialize only the object on the client side
1312 AliCodeTimerAuto("");
1314 AliReconstruction *reco = NULL;
1316 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1319 AliSysInfo::AddStamp("ReadInputInBegin");
1322 // Import ideal TGeo geometry and apply misalignment
1324 TString geom(gSystem->DirName(fGAliceFileName));
1325 geom += "/geometry.root";
1326 AliGeomManager::LoadGeometry(geom.Data());
1328 Abort("LoadGeometry", TSelector::kAbortProcess);
1331 AliSysInfo::AddStamp("LoadGeom");
1332 TString detsToCheck=fRunLocalReconstruction;
1333 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1334 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1337 AliSysInfo::AddStamp("CheckGeom");
1340 if (!MisalignGeometry(fLoadAlignData)) {
1341 Abort("MisalignGeometry", TSelector::kAbortProcess);
1344 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1345 AliSysInfo::AddStamp("MisalignGeom");
1348 Abort("InitGRP", TSelector::kAbortProcess);
1351 AliSysInfo::AddStamp("InitGRP");
1354 Abort("LoadCDB", TSelector::kAbortProcess);
1357 AliSysInfo::AddStamp("LoadCDB");
1359 // Read the reconstruction parameters from OCDB
1360 if (!InitRecoParams()) {
1361 AliWarning("Not all detectors have correct RecoParam objects initialized");
1363 AliSysInfo::AddStamp("InitRecoParams");
1365 if (fInput && gProof) {
1366 if (reco) *reco = *this;
1368 gProof->AddInputData(gGeoManager,kTRUE);
1370 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1371 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1372 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1373 magFieldMap->SetName("MagneticFieldMap");
1374 gProof->AddInputData(magFieldMap,kTRUE);
1379 //_____________________________________________________________________________
1380 void AliReconstruction::SlaveBegin(TTree*)
1382 // Initialization related to run-loader,
1383 // vertexer, trackers, recontructors
1384 // In proof mode it is executed on the slave
1385 AliCodeTimerAuto("");
1387 TProofOutputFile *outProofFile = NULL;
1389 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1392 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1394 AliGeomManager::SetGeometry(tgeo);
1396 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1397 Int_t runNumber = -1;
1398 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1399 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1400 man->SetCacheFlag(kTRUE);
1401 man->SetLock(kTRUE);
1405 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1406 TGeoGlobalMagField::Instance()->SetField(map);
1408 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1409 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1410 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1411 fOutput->Add(outProofFile);
1413 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1416 // get the run loader
1417 if (!InitRunLoader()) {
1418 Abort("InitRunLoader", TSelector::kAbortProcess);
1421 AliSysInfo::AddStamp("LoadLoader");
1423 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1426 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1427 Abort("CreateTrackers", TSelector::kAbortProcess);
1430 AliSysInfo::AddStamp("CreateTrackers");
1432 // create the ESD output file and tree
1433 if (!outProofFile) {
1434 ffile = TFile::Open("AliESDs.root", "RECREATE");
1435 ffile->SetCompressionLevel(2);
1436 if (!ffile->IsOpen()) {
1437 Abort("OpenESDFile", TSelector::kAbortProcess);
1442 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1443 Abort(Form("Problems opening output PROOF file: %s/%s",
1444 outProofFile->GetDir(), outProofFile->GetFileName()),
1445 TSelector::kAbortProcess);
1450 ftree = new TTree("esdTree", "Tree with ESD objects");
1451 fesd = new AliESDEvent();
1452 fesd->CreateStdContent();
1454 fesd->WriteToTree(ftree);
1455 if (fWriteESDfriend) {
1457 // Since we add the branch manually we must
1458 // book and add it after WriteToTree
1459 // otherwise it is created twice,
1460 // once via writetotree and once here.
1461 // The case for AliESDfriend is now
1462 // caught also in AlIESDEvent::WriteToTree but
1463 // be careful when changing the name (AliESDfriend is not
1464 // a TNamed so we had to hardwire it)
1465 fesdf = new AliESDfriend();
1466 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1467 br->SetFile("AliESDfriends.root");
1468 fesd->AddObject(fesdf);
1470 ftree->GetUserInfo()->Add(fesd);
1472 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1473 fhltesd = new AliESDEvent();
1474 fhltesd->CreateStdContent();
1476 // read the ESD template from CDB
1477 // HLT is allowed to put non-std content to its ESD, the non-std
1478 // objects need to be created before invocation of WriteToTree in
1479 // order to create all branches. Initialization is done from an
1480 // ESD layout template in CDB
1481 AliCDBManager* man = AliCDBManager::Instance();
1482 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1483 AliCDBEntry* hltESDConfig=NULL;
1484 if (man->GetId(hltESDConfigPath)!=NULL &&
1485 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1486 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1488 // init all internal variables from the list of objects
1489 pESDLayout->GetStdContent();
1491 // copy content and create non-std objects
1492 *fhltesd=*pESDLayout;
1495 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1496 hltESDConfigPath.GetPath().Data()));
1500 fhltesd->WriteToTree(fhlttree);
1501 fhlttree->GetUserInfo()->Add(fhltesd);
1503 ProcInfo_t procInfo;
1504 gSystem->GetProcInfo(&procInfo);
1505 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1508 //Initialize the QA and start of cycle
1509 if (fRunQA || fRunGlobalQA)
1512 //Initialize the Plane Efficiency framework
1513 if (fRunPlaneEff && !InitPlaneEff()) {
1514 Abort("InitPlaneEff", TSelector::kAbortProcess);
1518 if (strcmp(gProgName,"alieve") == 0)
1519 fRunAliEVE = InitAliEVE();
1524 //_____________________________________________________________________________
1525 Bool_t AliReconstruction::Process(Long64_t entry)
1527 // run the reconstruction over a single entry
1528 // from the chain with raw data
1529 AliCodeTimerAuto("");
1531 TTree *currTree = fChain->GetTree();
1532 AliRawVEvent *event = NULL;
1533 currTree->SetBranchAddress("rawevent",&event);
1534 currTree->GetEntry(entry);
1535 fRawReader = new AliRawReaderRoot(event);
1536 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1544 //_____________________________________________________________________________
1545 void AliReconstruction::Init(TTree *tree)
1548 AliError("The input tree is not found!");
1554 //_____________________________________________________________________________
1555 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1557 // run the reconstruction over a single event
1558 // The event loop is steered in Run method
1560 AliCodeTimerAuto("");
1562 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1563 fRunLoader->SetEventNumber(iEvent);
1564 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1566 fRunLoader->TreeE()->Fill();
1567 if (fRawReader && fRawReader->UseAutoSaveESD())
1568 fRunLoader->TreeE()->AutoSave("SaveSelf");
1571 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1575 AliInfo(Form("processing event %d", iEvent));
1577 fRunLoader->GetEvent(iEvent);
1579 // Fill Event-info object
1581 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo);
1582 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1584 // Set the reco-params
1586 TString detStr = fLoadCDB;
1587 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1588 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1589 AliReconstructor *reconstructor = GetReconstructor(iDet);
1590 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1591 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1592 reconstructor->SetRecoParam(par);
1594 fQAManager->SetRecoParam(iDet, par) ;
1602 fQAManager->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1603 fQAManager->RunOneEvent(fRawReader) ;
1605 // local single event reconstruction
1606 if (!fRunLocalReconstruction.IsNull()) {
1607 TString detectors=fRunLocalReconstruction;
1608 // run HLT event reconstruction first
1609 // ;-( IsSelected changes the string
1610 if (IsSelected("HLT", detectors) &&
1611 !RunLocalEventReconstruction("HLT")) {
1612 if (fStopOnError) {CleanUp(); return kFALSE;}
1614 detectors=fRunLocalReconstruction;
1615 detectors.ReplaceAll("HLT", "");
1616 if (!RunLocalEventReconstruction(detectors)) {
1617 if (fStopOnError) {CleanUp(); return kFALSE;}
1621 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1622 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1623 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1624 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1626 // Set magnetic field from the tracker
1627 fesd->SetMagneticField(AliTracker::GetBz());
1628 fhltesd->SetMagneticField(AliTracker::GetBz());
1630 // Set most probable pt, for B=0 tracking
1631 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1632 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1633 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1635 // Fill raw-data error log into the ESD
1636 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1639 if (fRunVertexFinder) {
1640 if (!RunVertexFinder(fesd)) {
1641 if (fStopOnError) {CleanUp(); return kFALSE;}
1645 // For Plane Efficiency: run the SPD trackleter
1646 if (fRunPlaneEff && fSPDTrackleter) {
1647 if (!RunSPDTrackleting(fesd)) {
1648 if (fStopOnError) {CleanUp(); return kFALSE;}
1653 if (!fRunTracking.IsNull()) {
1654 if (fRunMuonTracking) {
1655 if (!RunMuonTracking(fesd)) {
1656 if (fStopOnError) {CleanUp(); return kFALSE;}
1662 if (!fRunTracking.IsNull()) {
1663 if (!RunTracking(fesd)) {
1664 if (fStopOnError) {CleanUp(); return kFALSE;}
1669 if (!fFillESD.IsNull()) {
1670 TString detectors=fFillESD;
1671 // run HLT first and on hltesd
1672 // ;-( IsSelected changes the string
1673 if (IsSelected("HLT", detectors) &&
1674 !FillESD(fhltesd, "HLT")) {
1675 if (fStopOnError) {CleanUp(); return kFALSE;}
1678 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1679 if (detectors.Contains("ALL")) {
1681 for (Int_t idet=0; idet<kNDetectors; ++idet){
1682 detectors += fgkDetectorName[idet];
1686 detectors.ReplaceAll("HLT", "");
1687 if (!FillESD(fesd, detectors)) {
1688 if (fStopOnError) {CleanUp(); return kFALSE;}
1692 // fill Event header information from the RawEventHeader
1693 if (fRawReader){FillRawEventHeaderESD(fesd);}
1696 AliESDpid::MakePID(fesd);
1698 if (fFillTriggerESD) {
1699 if (!FillTriggerESD(fesd)) {
1700 if (fStopOnError) {CleanUp(); return kFALSE;}
1707 // Propagate track to the beam pipe (if not already done by ITS)
1709 const Int_t ntracks = fesd->GetNumberOfTracks();
1710 const Double_t kBz = fesd->GetMagneticField();
1711 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1714 UShort_t *selectedIdx=new UShort_t[ntracks];
1716 for (Int_t itrack=0; itrack<ntracks; itrack++){
1717 const Double_t kMaxStep = 1; //max step over the material
1720 AliESDtrack *track = fesd->GetTrack(itrack);
1721 if (!track) continue;
1723 AliExternalTrackParam *tpcTrack =
1724 (AliExternalTrackParam *)track->GetTPCInnerParam();
1728 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1731 Int_t n=trkArray.GetEntriesFast();
1732 selectedIdx[n]=track->GetID();
1733 trkArray.AddLast(tpcTrack);
1736 //Tracks refitted by ITS should already be at the SPD vertex
1737 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1740 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1741 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1746 // Improve the reconstructed primary vertex position using the tracks
1748 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1749 if(fesd->GetPrimaryVertexSPD()) {
1750 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1751 if(vtitle.Contains("cosmics")) {
1752 runVertexFinderTracks=kFALSE;
1756 if (runVertexFinderTracks) {
1757 // TPC + ITS primary vertex
1758 ftVertexer->SetITSMode();
1759 ftVertexer->SetConstraintOff();
1760 // get cuts for vertexer from AliGRPRecoParam
1762 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1763 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1764 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1765 ftVertexer->SetCuts(cutsVertexer);
1766 delete [] cutsVertexer; cutsVertexer = NULL;
1767 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1768 ftVertexer->SetVtxStart(fDiamondProfile);
1770 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1772 if (pvtx->GetStatus()) {
1773 fesd->SetPrimaryVertexTracks(pvtx);
1774 for (Int_t i=0; i<ntracks; i++) {
1775 AliESDtrack *t = fesd->GetTrack(i);
1776 t->RelateToVertex(pvtx, kBz, kVeryBig);
1781 // TPC-only primary vertex
1782 ftVertexer->SetTPCMode();
1783 ftVertexer->SetConstraintOff();
1784 // get cuts for vertexer from AliGRPRecoParam
1786 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1787 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1788 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1789 ftVertexer->SetCuts(cutsVertexer);
1790 delete [] cutsVertexer; cutsVertexer = NULL;
1791 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1792 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1794 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1796 if (pvtx->GetStatus()) {
1797 fesd->SetPrimaryVertexTPC(pvtx);
1798 for (Int_t i=0; i<ntracks; i++) {
1799 AliESDtrack *t = fesd->GetTrack(i);
1800 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1806 delete[] selectedIdx;
1808 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1813 AliV0vertexer vtxer;
1814 vtxer.Tracks2V0vertices(fesd);
1816 if (fRunCascadeFinder) {
1818 AliCascadeVertexer cvtxer;
1819 cvtxer.V0sTracks2CascadeVertices(fesd);
1824 if (fCleanESD) CleanESD(fesd);
1827 fQAManager->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1828 fQAManager->RunOneEvent(fesd) ;
1831 AliQADataMaker *qadm = fQAManager->GetQADataMaker(AliQAv1::kGLOBAL);
1832 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1833 if (qadm && fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1834 qadm->Exec(AliQAv1::kESDS, fesd);
1837 if (fWriteESDfriend) {
1838 // fesdf->~AliESDfriend();
1839 // new (fesdf) AliESDfriend(); // Reset...
1840 fesd->GetESDfriend(fesdf);
1844 // Auto-save the ESD tree in case of prompt reco @P2
1845 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1846 ftree->AutoSave("SaveSelf");
1847 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1848 if (friendfile) friendfile->Save();
1855 if (fRunAliEVE) RunAliEVE();
1859 if (fWriteESDfriend) {
1860 fesdf->~AliESDfriend();
1861 new (fesdf) AliESDfriend(); // Reset...
1864 ProcInfo_t procInfo;
1865 gSystem->GetProcInfo(&procInfo);
1866 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1869 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1870 if (fReconstructor[iDet])
1871 fReconstructor[iDet]->SetRecoParam(NULL);
1874 if (fRunQA || fRunGlobalQA)
1875 fQAManager->Increment() ;
1880 //_____________________________________________________________________________
1881 void AliReconstruction::SlaveTerminate()
1883 // Finalize the run on the slave side
1884 // Called after the exit
1885 // from the event loop
1886 AliCodeTimerAuto("");
1888 if (fIsNewRunLoader) { // galice.root didn't exist
1889 fRunLoader->WriteHeader("OVERWRITE");
1890 fRunLoader->CdGAFile();
1891 fRunLoader->Write(0, TObject::kOverwrite);
1894 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1895 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1897 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1898 cdbMapCopy->SetOwner(1);
1899 cdbMapCopy->SetName("cdbMap");
1900 TIter iter(cdbMap->GetTable());
1903 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1904 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1905 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1906 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1909 TList *cdbListCopy = new TList();
1910 cdbListCopy->SetOwner(1);
1911 cdbListCopy->SetName("cdbList");
1913 TIter iter2(cdbList);
1916 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1917 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1920 ftree->GetUserInfo()->Add(cdbMapCopy);
1921 ftree->GetUserInfo()->Add(cdbListCopy);
1926 if (fWriteESDfriend)
1927 ftree->SetBranchStatus("ESDfriend*",0);
1928 // we want to have only one tree version number
1929 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1932 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1933 if (fRunPlaneEff && !FinishPlaneEff()) {
1934 AliWarning("Finish PlaneEff evaluation failed");
1937 // End of cycle for the in-loop
1939 fQAManager->EndOfCycle() ;
1942 AliQADataMaker *qadm = fQAManager->GetQADataMaker(AliQAv1::kGLOBAL);
1944 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS)))
1945 qadm->EndOfCycle(AliQAv1::kRECPOINTS);
1946 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1947 qadm->EndOfCycle(AliQAv1::kESDS);
1952 if (fRunQA || fRunGlobalQA) {
1954 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1955 TString qaOutputFile = outputFileName->GetTitle();
1956 qaOutputFile.ReplaceAll(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()),
1957 Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
1958 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
1959 qaProofFile->SetOutputFileName(qaOutputFile.Data());
1960 fOutput->Add(qaProofFile);
1961 MergeQA(qaProofFile->GetFileName());
1973 //_____________________________________________________________________________
1974 void AliReconstruction::Terminate()
1976 // Create tags for the events in the ESD tree (the ESD tree is always present)
1977 // In case of empty events the tags will contain dummy values
1978 AliCodeTimerAuto("");
1980 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1981 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
1983 // Cleanup of CDB manager: cache and active storages!
1984 AliCDBManager::Instance()->ClearCache();
1987 //_____________________________________________________________________________
1988 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1990 // run the local reconstruction
1992 static Int_t eventNr=0;
1993 AliCodeTimerAuto("")
1995 TString detStr = detectors;
1996 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1997 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1998 AliReconstructor* reconstructor = GetReconstructor(iDet);
1999 if (!reconstructor) continue;
2000 AliLoader* loader = fLoader[iDet];
2001 // Matthias April 2008: temporary fix to run HLT reconstruction
2002 // although the HLT loader is missing
2003 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2005 reconstructor->Reconstruct(fRawReader, NULL);
2008 reconstructor->Reconstruct(dummy, NULL);
2013 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2016 // conversion of digits
2017 if (fRawReader && reconstructor->HasDigitConversion()) {
2018 AliInfo(Form("converting raw data digits into root objects for %s",
2019 fgkDetectorName[iDet]));
2020 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2021 // fgkDetectorName[iDet]));
2022 loader->LoadDigits("update");
2023 loader->CleanDigits();
2024 loader->MakeDigitsContainer();
2025 TTree* digitsTree = loader->TreeD();
2026 reconstructor->ConvertDigits(fRawReader, digitsTree);
2027 loader->WriteDigits("OVERWRITE");
2028 loader->UnloadDigits();
2030 // local reconstruction
2031 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2032 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2033 loader->LoadRecPoints("update");
2034 loader->CleanRecPoints();
2035 loader->MakeRecPointsContainer();
2036 TTree* clustersTree = loader->TreeR();
2037 if (fRawReader && !reconstructor->HasDigitConversion()) {
2038 reconstructor->Reconstruct(fRawReader, clustersTree);
2040 loader->LoadDigits("read");
2041 TTree* digitsTree = loader->TreeD();
2043 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2044 if (fStopOnError) return kFALSE;
2046 reconstructor->Reconstruct(digitsTree, clustersTree);
2048 fQAManager->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2049 fQAManager->RunOneEventInOneDetector(iDet, digitsTree) ;
2052 loader->UnloadDigits();
2055 TString detQAStr(fQADetectors) ;
2057 fQAManager->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2058 fQAManager->RunOneEventInOneDetector(iDet, clustersTree) ;
2060 loader->WriteRecPoints("OVERWRITE");
2061 loader->UnloadRecPoints();
2062 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2064 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2065 AliError(Form("the following detectors were not found: %s",
2067 if (fStopOnError) return kFALSE;
2072 //_____________________________________________________________________________
2073 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2075 // run the SPD trackleting (for SPD efficiency purpouses)
2077 AliCodeTimerAuto("")
2079 Double_t vtxPos[3] = {0, 0, 0};
2080 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2082 TArrayF mcVertex(3);
2084 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2085 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2086 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2089 const AliESDVertex *vertex = esd->GetVertex();
2091 AliWarning("Vertex not found");
2094 vertex->GetXYZ(vtxPos);
2095 vertex->GetSigmaXYZ(vtxErr);
2096 if (fSPDTrackleter) {
2097 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2100 fLoader[0]->LoadRecPoints("read");
2101 TTree* tree = fLoader[0]->TreeR();
2103 AliError("Can't get the ITS cluster tree");
2106 fSPDTrackleter->LoadClusters(tree);
2107 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2109 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2110 AliError("AliITSTrackleterSPDEff Clusters2Tracks failed");
2111 // fLoader[0]->UnloadRecPoints();
2114 //fSPDTrackleter->UnloadRecPoints();
2116 AliWarning("SPDTrackleter not available");
2122 //_____________________________________________________________________________
2123 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2125 // run the barrel tracking
2127 AliCodeTimerAuto("")
2129 AliVertexer *vertexer = CreateVertexer();
2130 if (!vertexer) return kFALSE;
2132 AliInfo("running the ITS vertex finder");
2133 AliESDVertex* vertex = NULL;
2135 fLoader[0]->LoadRecPoints();
2136 TTree* cltree = fLoader[0]->TreeR();
2138 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2139 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2142 AliError("Can't get the ITS cluster tree");
2144 fLoader[0]->UnloadRecPoints();
2147 AliError("Can't get the ITS loader");
2150 AliWarning("Vertex not found");
2151 vertex = new AliESDVertex();
2152 vertex->SetName("default");
2155 vertex->SetName("reconstructed");
2160 vertex->GetXYZ(vtxPos);
2161 vertex->GetSigmaXYZ(vtxErr);
2163 esd->SetPrimaryVertexSPD(vertex);
2164 AliESDVertex *vpileup = NULL;
2165 Int_t novertices = 0;
2166 vpileup = vertexer->GetAllVertices(novertices);
2168 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2170 // if SPD multiplicity has been determined, it is stored in the ESD
2171 AliMultiplicity *mult = vertexer->GetMultiplicity();
2172 if(mult)esd->SetMultiplicity(mult);
2174 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2175 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2184 //_____________________________________________________________________________
2185 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2187 // run the HLT barrel tracking
2189 AliCodeTimerAuto("")
2192 AliError("Missing runLoader!");
2196 AliInfo("running HLT tracking");
2198 // Get a pointer to the HLT reconstructor
2199 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2200 if (!reconstructor) return kFALSE;
2203 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2204 TString detName = fgkDetectorName[iDet];
2205 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2206 reconstructor->SetOption(detName.Data());
2207 AliTracker *tracker = reconstructor->CreateTracker();
2209 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2210 if (fStopOnError) return kFALSE;
2214 Double_t vtxErr[3]={0.005,0.005,0.010};
2215 const AliESDVertex *vertex = esd->GetVertex();
2216 vertex->GetXYZ(vtxPos);
2217 tracker->SetVertex(vtxPos,vtxErr);
2219 fLoader[iDet]->LoadRecPoints("read");
2220 TTree* tree = fLoader[iDet]->TreeR();
2222 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2225 tracker->LoadClusters(tree);
2227 if (tracker->Clusters2Tracks(esd) != 0) {
2228 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2232 tracker->UnloadClusters();
2240 //_____________________________________________________________________________
2241 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2243 // run the muon spectrometer tracking
2245 AliCodeTimerAuto("")
2248 AliError("Missing runLoader!");
2251 Int_t iDet = 7; // for MUON
2253 AliInfo("is running...");
2255 // Get a pointer to the MUON reconstructor
2256 AliReconstructor *reconstructor = GetReconstructor(iDet);
2257 if (!reconstructor) return kFALSE;
2260 TString detName = fgkDetectorName[iDet];
2261 AliDebug(1, Form("%s tracking", detName.Data()));
2262 AliTracker *tracker = reconstructor->CreateTracker();
2264 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2269 fLoader[iDet]->LoadRecPoints("read");
2271 tracker->LoadClusters(fLoader[iDet]->TreeR());
2273 Int_t rv = tracker->Clusters2Tracks(esd);
2277 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2281 fLoader[iDet]->UnloadRecPoints();
2283 tracker->UnloadClusters();
2291 //_____________________________________________________________________________
2292 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2294 // run the barrel tracking
2295 static Int_t eventNr=0;
2296 AliCodeTimerAuto("")
2298 AliInfo("running tracking");
2300 //Fill the ESD with the T0 info (will be used by the TOF)
2301 if (fReconstructor[11] && fLoader[11]) {
2302 fLoader[11]->LoadRecPoints("READ");
2303 TTree *treeR = fLoader[11]->TreeR();
2305 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2309 // pass 1: TPC + ITS inwards
2310 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2311 if (!fTracker[iDet]) continue;
2312 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2315 fLoader[iDet]->LoadRecPoints("read");
2316 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2317 TTree* tree = fLoader[iDet]->TreeR();
2319 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2322 fTracker[iDet]->LoadClusters(tree);
2323 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2325 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2326 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2329 // preliminary PID in TPC needed by the ITS tracker
2331 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2332 AliESDpid::MakePID(esd);
2334 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2337 // pass 2: ALL backwards
2339 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2340 if (!fTracker[iDet]) continue;
2341 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2344 if (iDet > 1) { // all except ITS, TPC
2346 fLoader[iDet]->LoadRecPoints("read");
2347 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2348 tree = fLoader[iDet]->TreeR();
2350 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2353 fTracker[iDet]->LoadClusters(tree);
2354 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2358 if (iDet>1) // start filling residuals for the "outer" detectors
2359 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2361 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2362 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2367 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2368 fTracker[iDet]->UnloadClusters();
2369 fLoader[iDet]->UnloadRecPoints();
2371 // updated PID in TPC needed by the ITS tracker -MI
2373 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2374 AliESDpid::MakePID(esd);
2376 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2378 //stop filling residuals for the "outer" detectors
2379 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2381 // pass 3: TRD + TPC + ITS refit inwards
2383 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2384 if (!fTracker[iDet]) continue;
2385 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2388 if (iDet<2) // start filling residuals for TPC and ITS
2389 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2391 if (fTracker[iDet]->RefitInward(esd) != 0) {
2392 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2395 // run postprocessing
2396 if (fTracker[iDet]->PostProcess(esd) != 0) {
2397 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2400 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2403 // write space-points to the ESD in case alignment data output
2405 if (fWriteAlignmentData)
2406 WriteAlignmentData(esd);
2408 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2409 if (!fTracker[iDet]) continue;
2411 fTracker[iDet]->UnloadClusters();
2412 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2413 fLoader[iDet]->UnloadRecPoints();
2414 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2416 // stop filling residuals for TPC and ITS
2417 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2423 //_____________________________________________________________________________
2424 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2426 // Remove the data which are not needed for the physics analysis.
2429 Int_t nTracks=esd->GetNumberOfTracks();
2430 Int_t nV0s=esd->GetNumberOfV0s();
2432 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2434 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2435 Bool_t rc=esd->Clean(cleanPars);
2437 nTracks=esd->GetNumberOfTracks();
2438 nV0s=esd->GetNumberOfV0s();
2440 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2445 //_____________________________________________________________________________
2446 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2448 // fill the event summary data
2450 AliCodeTimerAuto("")
2451 static Int_t eventNr=0;
2452 TString detStr = detectors;
2454 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2455 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2456 AliReconstructor* reconstructor = GetReconstructor(iDet);
2457 if (!reconstructor) continue;
2458 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2459 TTree* clustersTree = NULL;
2460 if (fLoader[iDet]) {
2461 fLoader[iDet]->LoadRecPoints("read");
2462 clustersTree = fLoader[iDet]->TreeR();
2463 if (!clustersTree) {
2464 AliError(Form("Can't get the %s clusters tree",
2465 fgkDetectorName[iDet]));
2466 if (fStopOnError) return kFALSE;
2469 if (fRawReader && !reconstructor->HasDigitConversion()) {
2470 reconstructor->FillESD(fRawReader, clustersTree, esd);
2472 TTree* digitsTree = NULL;
2473 if (fLoader[iDet]) {
2474 fLoader[iDet]->LoadDigits("read");
2475 digitsTree = fLoader[iDet]->TreeD();
2477 AliError(Form("Can't get the %s digits tree",
2478 fgkDetectorName[iDet]));
2479 if (fStopOnError) return kFALSE;
2482 reconstructor->FillESD(digitsTree, clustersTree, esd);
2483 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2485 if (fLoader[iDet]) {
2486 fLoader[iDet]->UnloadRecPoints();
2490 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2491 AliError(Form("the following detectors were not found: %s",
2493 if (fStopOnError) return kFALSE;
2495 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2500 //_____________________________________________________________________________
2501 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2503 // Reads the trigger decision which is
2504 // stored in Trigger.root file and fills
2505 // the corresponding esd entries
2507 AliCodeTimerAuto("")
2509 AliInfo("Filling trigger information into the ESD");
2512 AliCTPRawStream input(fRawReader);
2513 if (!input.Next()) {
2514 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2517 if (esd->GetTriggerMask() != input.GetClassMask())
2518 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2519 input.GetClassMask(),esd->GetTriggerMask()));
2520 if (esd->GetOrbitNumber() != input.GetOrbitID())
2521 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2522 input.GetOrbitID(),esd->GetOrbitNumber()));
2523 if (esd->GetBunchCrossNumber() != input.GetBCID())
2524 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2525 input.GetBCID(),esd->GetBunchCrossNumber()));
2528 // Here one has to add the filling of trigger inputs and
2529 // interaction records
2539 //_____________________________________________________________________________
2540 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2543 // Filling information from RawReader Header
2546 if (!fRawReader) return kFALSE;
2548 AliInfo("Filling information from RawReader Header");
2550 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2551 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2552 esd->SetPeriodNumber(fRawReader->GetPeriod());
2554 esd->SetTimeStamp(fRawReader->GetTimestamp());
2555 esd->SetEventType(fRawReader->GetType());
2561 //_____________________________________________________________________________
2562 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2564 // check whether detName is contained in detectors
2565 // if yes, it is removed from detectors
2567 // check if all detectors are selected
2568 if ((detectors.CompareTo("ALL") == 0) ||
2569 detectors.BeginsWith("ALL ") ||
2570 detectors.EndsWith(" ALL") ||
2571 detectors.Contains(" ALL ")) {
2576 // search for the given detector
2577 Bool_t result = kFALSE;
2578 if ((detectors.CompareTo(detName) == 0) ||
2579 detectors.BeginsWith(detName+" ") ||
2580 detectors.EndsWith(" "+detName) ||
2581 detectors.Contains(" "+detName+" ")) {
2582 detectors.ReplaceAll(detName, "");
2586 // clean up the detectors string
2587 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2588 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2589 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2594 //_____________________________________________________________________________
2595 Bool_t AliReconstruction::InitRunLoader()
2597 // get or create the run loader
2599 if (gAlice) delete gAlice;
2602 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2603 // load all base libraries to get the loader classes
2604 TString libs = gSystem->GetLibraries();
2605 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2606 TString detName = fgkDetectorName[iDet];
2607 if (detName == "HLT") continue;
2608 if (libs.Contains("lib" + detName + "base.so")) continue;
2609 gSystem->Load("lib" + detName + "base.so");
2611 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2613 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2618 fRunLoader->CdGAFile();
2619 fRunLoader->LoadgAlice();
2621 //PH This is a temporary fix to give access to the kinematics
2622 //PH that is needed for the labels of ITS clusters
2623 fRunLoader->LoadHeader();
2624 fRunLoader->LoadKinematics();
2626 } else { // galice.root does not exist
2628 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2630 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2631 AliConfig::GetDefaultEventFolderName(),
2634 AliError(Form("could not create run loader in file %s",
2635 fGAliceFileName.Data()));
2639 fIsNewRunLoader = kTRUE;
2640 fRunLoader->MakeTree("E");
2642 if (fNumberOfEventsPerFile > 0)
2643 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2645 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2651 //_____________________________________________________________________________
2652 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2654 // get the reconstructor object and the loader for a detector
2656 if (fReconstructor[iDet]) {
2657 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2658 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2659 fReconstructor[iDet]->SetRecoParam(par);
2661 return fReconstructor[iDet];
2664 // load the reconstructor object
2665 TPluginManager* pluginManager = gROOT->GetPluginManager();
2666 TString detName = fgkDetectorName[iDet];
2667 TString recName = "Ali" + detName + "Reconstructor";
2669 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2671 AliReconstructor* reconstructor = NULL;
2672 // first check if a plugin is defined for the reconstructor
2673 TPluginHandler* pluginHandler =
2674 pluginManager->FindHandler("AliReconstructor", detName);
2675 // if not, add a plugin for it
2676 if (!pluginHandler) {
2677 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2678 TString libs = gSystem->GetLibraries();
2679 if (libs.Contains("lib" + detName + "base.so") ||
2680 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2681 pluginManager->AddHandler("AliReconstructor", detName,
2682 recName, detName + "rec", recName + "()");
2684 pluginManager->AddHandler("AliReconstructor", detName,
2685 recName, detName, recName + "()");
2687 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2689 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2690 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2692 if (reconstructor) {
2693 TObject* obj = fOptions.FindObject(detName.Data());
2694 if (obj) reconstructor->SetOption(obj->GetTitle());
2695 reconstructor->Init();
2696 fReconstructor[iDet] = reconstructor;
2699 // get or create the loader
2700 if (detName != "HLT") {
2701 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2702 if (!fLoader[iDet]) {
2703 AliConfig::Instance()
2704 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2706 // first check if a plugin is defined for the loader
2708 pluginManager->FindHandler("AliLoader", detName);
2709 // if not, add a plugin for it
2710 if (!pluginHandler) {
2711 TString loaderName = "Ali" + detName + "Loader";
2712 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2713 pluginManager->AddHandler("AliLoader", detName,
2714 loaderName, detName + "base",
2715 loaderName + "(const char*, TFolder*)");
2716 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2718 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2720 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2721 fRunLoader->GetEventFolder());
2723 if (!fLoader[iDet]) { // use default loader
2724 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2726 if (!fLoader[iDet]) {
2727 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2728 if (fStopOnError) return NULL;
2730 fRunLoader->AddLoader(fLoader[iDet]);
2731 fRunLoader->CdGAFile();
2732 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2733 fRunLoader->Write(0, TObject::kOverwrite);
2738 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2739 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2740 reconstructor->SetRecoParam(par);
2742 return reconstructor;
2745 //_____________________________________________________________________________
2746 AliVertexer* AliReconstruction::CreateVertexer()
2748 // create the vertexer
2749 // Please note that the caller is the owner of the
2752 AliVertexer* vertexer = NULL;
2753 AliReconstructor* itsReconstructor = GetReconstructor(0);
2754 if (itsReconstructor) {
2755 vertexer = itsReconstructor->CreateVertexer();
2758 AliWarning("couldn't create a vertexer for ITS");
2764 //_____________________________________________________________________________
2765 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2767 // create the trackers
2768 AliInfo("Creating trackers");
2770 TString detStr = detectors;
2771 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2772 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2773 AliReconstructor* reconstructor = GetReconstructor(iDet);
2774 if (!reconstructor) continue;
2775 TString detName = fgkDetectorName[iDet];
2776 if (detName == "HLT") {
2777 fRunHLTTracking = kTRUE;
2780 if (detName == "MUON") {
2781 fRunMuonTracking = kTRUE;
2786 fTracker[iDet] = reconstructor->CreateTracker();
2787 if (!fTracker[iDet] && (iDet < 7)) {
2788 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2789 if (fStopOnError) return kFALSE;
2791 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2797 //_____________________________________________________________________________
2798 void AliReconstruction::CleanUp()
2800 // delete trackers and the run loader and close and delete the file
2802 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2803 delete fReconstructor[iDet];
2804 fReconstructor[iDet] = NULL;
2805 fLoader[iDet] = NULL;
2806 delete fTracker[iDet];
2807 fTracker[iDet] = NULL;
2812 delete fSPDTrackleter;
2813 fSPDTrackleter = NULL;
2822 delete fParentRawReader;
2823 fParentRawReader=NULL;
2831 TGeoGlobalMagField::Instance()->SetField(NULL);
2834 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2836 // Write space-points which are then used in the alignment procedures
2837 // For the moment only ITS, TPC, TRD and TOF
2839 Int_t ntracks = esd->GetNumberOfTracks();
2840 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2842 AliESDtrack *track = esd->GetTrack(itrack);
2845 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
2846 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2847 nsp += track->GetNcls(iDet);
2849 if (iDet==0) { // ITS "extra" clusters
2850 track->GetClusters(iDet,idx);
2851 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2856 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2857 track->SetTrackPointArray(sp);
2859 for (Int_t iDet = 5; iDet >= 0; iDet--) {
2860 AliTracker *tracker = fTracker[iDet];
2861 if (!tracker) continue;
2862 Int_t nspdet = track->GetClusters(iDet,idx);
2864 if (iDet==0) // ITS "extra" clusters
2865 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2867 if (nspdet <= 0) continue;
2871 while (isp2 < nspdet) {
2872 Bool_t isvalid=kTRUE;
2874 Int_t index=idx[isp++];
2875 if (index < 0) continue;
2877 TString dets = fgkDetectorName[iDet];
2878 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2879 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2880 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2881 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2882 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2884 isvalid = tracker->GetTrackPoint(index,p);
2887 if (!isvalid) continue;
2888 if (iDet==0 && (isp-1)>=6) p.SetExtra();
2889 sp->AddPoint(isptrack,&p); isptrack++;
2896 //_____________________________________________________________________________
2897 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2899 // The method reads the raw-data error log
2900 // accumulated within the rawReader.
2901 // It extracts the raw-data errors related to
2902 // the current event and stores them into
2903 // a TClonesArray inside the esd object.
2905 if (!fRawReader) return;
2907 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2909 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2911 if (iEvent != log->GetEventNumber()) continue;
2913 esd->AddRawDataErrorLog(log);
2918 //_____________________________________________________________________________
2919 void AliReconstruction::CheckQA()
2921 // check the QA of SIM for this run and remove the detectors
2922 // with status Fatal
2924 // TString newRunLocalReconstruction ;
2925 // TString newRunTracking ;
2926 // TString newFillESD ;
2928 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
2929 // TString detName(AliQAv1::GetDetName(iDet)) ;
2930 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
2931 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
2932 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
2933 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
2935 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
2936 // fRunLocalReconstruction.Contains("ALL") ) {
2937 // newRunLocalReconstruction += detName ;
2938 // newRunLocalReconstruction += " " ;
2940 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
2941 // fRunTracking.Contains("ALL") ) {
2942 // newRunTracking += detName ;
2943 // newRunTracking += " " ;
2945 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
2946 // fFillESD.Contains("ALL") ) {
2947 // newFillESD += detName ;
2948 // newFillESD += " " ;
2952 // fRunLocalReconstruction = newRunLocalReconstruction ;
2953 // fRunTracking = newRunTracking ;
2954 // fFillESD = newFillESD ;
2957 //_____________________________________________________________________________
2958 Int_t AliReconstruction::GetDetIndex(const char* detector)
2960 // return the detector index corresponding to detector
2962 for (index = 0; index < kNDetectors ; index++) {
2963 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2968 //_____________________________________________________________________________
2969 Bool_t AliReconstruction::FinishPlaneEff() {
2971 // Here execute all the necessary operationis, at the end of the tracking phase,
2972 // in case that evaluation of PlaneEfficiencies was required for some detector.
2973 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2975 // This Preliminary version works only FOR ITS !!!!!
2976 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2979 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2982 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2983 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2984 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2985 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
2986 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2987 TString name=planeeff->GetName();
2989 TFile* pefile = TFile::Open(name, "RECREATE");
2990 ret=(Bool_t)planeeff->Write();
2992 if(planeeff->GetCreateHistos()) {
2993 TString hname=planeeff->GetName();
2994 hname+="Histo.root";
2995 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
2998 if(fSPDTrackleter) {
2999 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3000 TString name="AliITSPlaneEffSPDtracklet.root";
3001 TFile* pefile = TFile::Open(name, "RECREATE");
3002 ret=(Bool_t)planeeff->Write();
3004 AliESDEvent *dummy=NULL;
3005 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3010 //_____________________________________________________________________________
3011 Bool_t AliReconstruction::InitPlaneEff() {
3013 // Here execute all the necessary operations, before of the tracking phase,
3014 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3015 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3016 // which should be updated/recalculated.
3018 // This Preliminary version will work only FOR ITS !!!!!
3019 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3022 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3024 AliWarning(Form("Implementation of this method not yet completed !! Method return kTRUE"));
3026 fSPDTrackleter = NULL;
3027 AliReconstructor* itsReconstructor = GetReconstructor(0);
3028 if (itsReconstructor) {
3029 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3031 if (fSPDTrackleter) {
3032 AliInfo("Trackleter for SPD has been created");
3038 //_____________________________________________________________________________
3039 Bool_t AliReconstruction::InitAliEVE()
3041 // This method should be called only in case
3042 // AliReconstruction is run
3043 // within the alieve environment.
3044 // It will initialize AliEVE in a way
3045 // so that it can visualize event processed
3046 // by AliReconstruction.
3047 // The return flag shows whenever the
3048 // AliEVE initialization was successful or not.
3051 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3052 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3053 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3055 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3056 gROOT->ProcessLine("alieve_online_init()");
3061 //_____________________________________________________________________________
3062 void AliReconstruction::RunAliEVE()
3064 // Runs AliEVE visualisation of
3065 // the current event.
3066 // Should be executed only after
3067 // successful initialization of AliEVE.
3069 AliInfo("Running AliEVE...");
3070 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3074 //_____________________________________________________________________________
3075 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3077 // Allows to run QA for a selected set of detectors
3078 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3079 // all selected detectors run the same selected tasks
3081 if (!detAndAction.Contains(":")) {
3082 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3086 Int_t colon = detAndAction.Index(":") ;
3087 fQADetectors = detAndAction(0, colon) ;
3088 if (fQADetectors.Contains("ALL") )
3089 fQADetectors = fFillESD ;
3090 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3091 if (fQATasks.Contains("ALL") ) {
3092 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3094 fQATasks.ToUpper() ;
3096 if ( fQATasks.Contains("RAW") )
3097 tempo = Form("%d ", AliQAv1::kRAWS) ;
3098 if ( fQATasks.Contains("DIGIT") )
3099 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3100 if ( fQATasks.Contains("RECPOINT") )
3101 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3102 if ( fQATasks.Contains("ESD") )
3103 tempo += Form("%d ", AliQAv1::kESDS) ;
3105 if (fQATasks.IsNull()) {
3106 AliInfo("No QA requested\n") ;
3111 TString tempo(fQATasks) ;
3112 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3113 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3114 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3115 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3116 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3121 //_____________________________________________________________________________
3122 Bool_t AliReconstruction::InitRecoParams()
3124 // The method accesses OCDB and retrieves all
3125 // the available reco-param objects from there.
3127 Bool_t isOK = kTRUE;
3129 TString detStr = fLoadCDB;
3130 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3132 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3134 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3135 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3139 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3141 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3142 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3144 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3148 TObject *recoParamObj = entry->GetObject();
3149 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3150 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3151 // Registering them in AliRecoParam
3152 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3154 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3155 // The detector has only onse set of reco parameters
3156 // Registering it in AliRecoParam
3157 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3158 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3159 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3162 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3166 AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3170 if (AliDebugLevel() > 0) fRecoParam.Print();
3175 //_____________________________________________________________________________
3176 Bool_t AliReconstruction::GetEventInfo()
3178 // Fill the event info object
3180 AliCodeTimerAuto("")
3182 AliCentralTrigger *aCTP = NULL;
3184 fEventInfo.SetEventType(fRawReader->GetType());
3186 ULong64_t mask = fRawReader->GetClassMask();
3187 fEventInfo.SetTriggerMask(mask);
3188 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3189 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3191 aCTP = new AliCentralTrigger();
3192 TString configstr("");
3193 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3194 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3198 aCTP->SetClassMask(mask);
3199 aCTP->SetClusterMask(clmask);
3202 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3204 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3205 aCTP = fRunLoader->GetTrigger();
3206 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3207 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3210 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3215 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3217 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3218 if (fRawReader) delete aCTP;
3222 UChar_t clustmask = 0;
3224 ULong64_t trmask = fEventInfo.GetTriggerMask();
3225 const TObjArray& classesArray = config->GetClasses();
3226 Int_t nclasses = classesArray.GetEntriesFast();
3227 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3228 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3230 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3231 fesd->SetTriggerClass(trclass->GetName(),trindex);
3232 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3233 if (trmask & (1 << trindex)) {
3235 trclasses += trclass->GetName();
3237 clustmask |= trclass->GetCluster()->GetClusterMask();
3241 fEventInfo.SetTriggerClasses(trclasses);
3243 // Set the information in ESD
3244 fesd->SetTriggerMask(trmask);
3245 fesd->SetTriggerCluster(clustmask);
3247 if (!aCTP->CheckTriggeredDetectors()) {
3248 if (fRawReader) delete aCTP;
3252 if (fRawReader) delete aCTP;
3254 // We have to fill also the HLT decision here!!
3260 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3262 // Match the detector list found in the rec.C or the default 'ALL'
3263 // to the list found in the GRP (stored there by the shuttle PP which
3264 // gets the information from ECS)
3265 static TString resultList;
3266 TString detList = detectorList;
3270 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3271 if ((detectorMask >> iDet) & 0x1) {
3272 TString det = AliDAQ::OfflineModuleName(iDet);
3273 if ((detList.CompareTo("ALL") == 0) ||
3274 ((detList.BeginsWith("ALL ") ||
3275 detList.EndsWith(" ALL") ||
3276 detList.Contains(" ALL ")) &&
3277 !(detList.BeginsWith("-"+det+" ") ||
3278 detList.EndsWith(" -"+det) ||
3279 detList.Contains(" -"+det+" "))) ||
3280 (detList.CompareTo(det) == 0) ||
3281 detList.BeginsWith(det+" ") ||
3282 detList.EndsWith(" "+det) ||
3283 detList.Contains( " "+det+" " )) {
3284 if (!resultList.EndsWith(det + " ")) {
3293 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3294 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3295 if ((detList.CompareTo("ALL") == 0) ||
3296 ((detList.BeginsWith("ALL ") ||
3297 detList.EndsWith(" ALL") ||
3298 detList.Contains(" ALL ")) &&
3299 !(detList.BeginsWith("-"+hltDet+" ") ||
3300 detList.EndsWith(" -"+hltDet) ||
3301 detList.Contains(" -"+hltDet+" "))) ||
3302 (detList.CompareTo(hltDet) == 0) ||
3303 detList.BeginsWith(hltDet+" ") ||
3304 detList.EndsWith(" "+hltDet) ||
3305 detList.Contains( " "+hltDet+" " )) {
3306 resultList += hltDet;
3310 return resultList.Data();
3314 //______________________________________________________________________________
3315 void AliReconstruction::Abort(const char *method, EAbort what)
3317 // Abort processing. If what = kAbortProcess, the Process() loop will be
3318 // aborted. If what = kAbortFile, the current file in a chain will be
3319 // aborted and the processing will continue with the next file, if there
3320 // is no next file then Process() will be aborted. Abort() can also be
3321 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3322 // the SlaveTerminate() and Terminate() are always called. The abort flag
3323 // can be checked in these methods using GetAbort().
3325 // The method is overwritten in AliReconstruction for better handling of
3326 // reco specific errors
3328 if (!fStopOnError) return;
3332 TString whyMess = method;
3333 whyMess += " failed! Aborting...";
3335 AliError(whyMess.Data());
3338 TString mess = "Abort";
3339 if (fAbort == kAbortProcess)
3340 mess = "AbortProcess";
3341 else if (fAbort == kAbortFile)
3344 Info(mess, whyMess.Data());