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),
258 fSameQACycle(kFALSE),
259 fInitQACalled(kFALSE),
260 fWriteQAExpertData(kTRUE),
261 fRunPlaneEff(kFALSE),
270 fIsNewRunLoader(kFALSE),
274 // create reconstruction object with default parameters
277 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
278 fReconstructor[iDet] = NULL;
279 fLoader[iDet] = NULL;
280 fTracker[iDet] = NULL;
282 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
283 fQACycles[iDet] = 999999 ;
284 fQAWriteExpert[iDet] = kFALSE ;
290 //_____________________________________________________________________________
291 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
293 fUniformField(rec.fUniformField),
294 fRunVertexFinder(rec.fRunVertexFinder),
295 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
296 fRunHLTTracking(rec.fRunHLTTracking),
297 fRunMuonTracking(rec.fRunMuonTracking),
298 fRunV0Finder(rec.fRunV0Finder),
299 fRunCascadeFinder(rec.fRunCascadeFinder),
300 fStopOnError(rec.fStopOnError),
301 fWriteAlignmentData(rec.fWriteAlignmentData),
302 fWriteESDfriend(rec.fWriteESDfriend),
303 fFillTriggerESD(rec.fFillTriggerESD),
305 fCleanESD(rec.fCleanESD),
306 fV0DCAmax(rec.fV0DCAmax),
307 fV0CsPmin(rec.fV0CsPmin),
311 fRunLocalReconstruction(rec.fRunLocalReconstruction),
312 fRunTracking(rec.fRunTracking),
313 fFillESD(rec.fFillESD),
314 fLoadCDB(rec.fLoadCDB),
315 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
316 fGAliceFileName(rec.fGAliceFileName),
317 fRawInput(rec.fRawInput),
318 fESDOutput(rec.fESDOutput),
319 fEquipIdMap(rec.fEquipIdMap),
320 fFirstEvent(rec.fFirstEvent),
321 fLastEvent(rec.fLastEvent),
322 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
324 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
325 fLoadAlignData(rec.fLoadAlignData),
326 fUseHLTData(rec.fUseHLTData),
332 fParentRawReader(NULL),
334 fRecoParam(rec.fRecoParam),
336 fSPDTrackleter(NULL),
338 fDiamondProfileSPD(rec.fDiamondProfileSPD),
339 fDiamondProfile(rec.fDiamondProfile),
340 fDiamondProfileTPC(rec.fDiamondProfileTPC),
344 fAlignObjArray(rec.fAlignObjArray),
345 fCDBUri(rec.fCDBUri),
346 fQARefUri(rec.fQARefUri),
348 fInitCDBCalled(rec.fInitCDBCalled),
349 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
350 fQADetectors(rec.fQADetectors),
351 fQATasks(rec.fQATasks),
353 fRunGlobalQA(rec.fRunGlobalQA),
354 fSameQACycle(rec.fSameQACycle),
355 fInitQACalled(rec.fInitQACalled),
356 fWriteQAExpertData(rec.fWriteQAExpertData),
357 fRunPlaneEff(rec.fRunPlaneEff),
366 fIsNewRunLoader(rec.fIsNewRunLoader),
372 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
373 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
375 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
376 fReconstructor[iDet] = NULL;
377 fLoader[iDet] = NULL;
378 fTracker[iDet] = NULL;
381 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
382 fQACycles[iDet] = rec.fQACycles[iDet];
383 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
386 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
387 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
391 //_____________________________________________________________________________
392 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
394 // assignment operator
395 // Used in PROOF mode
396 // Be very careful while modifing it!
397 // Simple rules to follow:
398 // for persistent data members - use their assignment operators
399 // for non-persistent ones - do nothing or take the default values from constructor
400 // TSelector members should not be touched
401 if(&rec == this) return *this;
403 fUniformField = rec.fUniformField;
404 fRunVertexFinder = rec.fRunVertexFinder;
405 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
406 fRunHLTTracking = rec.fRunHLTTracking;
407 fRunMuonTracking = rec.fRunMuonTracking;
408 fRunV0Finder = rec.fRunV0Finder;
409 fRunCascadeFinder = rec.fRunCascadeFinder;
410 fStopOnError = rec.fStopOnError;
411 fWriteAlignmentData = rec.fWriteAlignmentData;
412 fWriteESDfriend = rec.fWriteESDfriend;
413 fFillTriggerESD = rec.fFillTriggerESD;
415 fCleanESD = rec.fCleanESD;
416 fV0DCAmax = rec.fV0DCAmax;
417 fV0CsPmin = rec.fV0CsPmin;
421 fRunLocalReconstruction = rec.fRunLocalReconstruction;
422 fRunTracking = rec.fRunTracking;
423 fFillESD = rec.fFillESD;
424 fLoadCDB = rec.fLoadCDB;
425 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
426 fGAliceFileName = rec.fGAliceFileName;
427 fRawInput = rec.fRawInput;
428 fESDOutput = rec.fESDOutput;
429 fEquipIdMap = rec.fEquipIdMap;
430 fFirstEvent = rec.fFirstEvent;
431 fLastEvent = rec.fLastEvent;
432 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
434 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
435 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
438 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
439 fLoadAlignData = rec.fLoadAlignData;
440 fUseHLTData = rec.fUseHLTData;
442 delete fRunInfo; fRunInfo = NULL;
443 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
445 fEventInfo = rec.fEventInfo;
449 fParentRawReader = NULL;
451 fRecoParam = rec.fRecoParam;
453 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
454 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
455 delete fLoader[iDet]; fLoader[iDet] = NULL;
456 delete fTracker[iDet]; fTracker[iDet] = NULL;
459 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
460 fQACycles[iDet] = rec.fQACycles[iDet];
461 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
464 delete fSPDTrackleter; fSPDTrackleter = NULL;
466 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
467 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
468 delete fDiamondProfile; fDiamondProfile = NULL;
469 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
470 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
471 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
473 delete fGRPData; fGRPData = NULL;
474 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
475 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
477 delete fAlignObjArray; fAlignObjArray = NULL;
480 fQARefUri = rec.fQARefUri;
481 fSpecCDBUri.Delete();
482 fInitCDBCalled = rec.fInitCDBCalled;
483 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
484 fQADetectors = rec.fQADetectors;
485 fQATasks = rec.fQATasks;
487 fRunGlobalQA = rec.fRunGlobalQA;
488 fSameQACycle = rec.fSameQACycle;
489 fInitQACalled = rec.fInitQACalled;
490 fWriteQAExpertData = rec.fWriteQAExpertData;
491 fRunPlaneEff = rec.fRunPlaneEff;
500 fIsNewRunLoader = rec.fIsNewRunLoader;
507 //_____________________________________________________________________________
508 AliReconstruction::~AliReconstruction()
515 if (fAlignObjArray) {
516 fAlignObjArray->Delete();
517 delete fAlignObjArray;
519 fSpecCDBUri.Delete();
521 AliCodeTimer::Instance()->Print();
524 //_____________________________________________________________________________
525 void AliReconstruction::InitQA()
527 //Initialize the QA and start of cycle
528 AliCodeTimerAuto("");
530 if (fInitQACalled) return;
531 fInitQACalled = kTRUE;
533 AliQAManager * qam = AliQAManager::QAManager("rec") ;
534 if (fWriteQAExpertData)
535 qam->SetWriteExpert() ;
537 if (qam->IsDefaultStorageSet()) {
538 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
539 AliWarning("Default QA reference storage has been already set !");
540 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fQARefUri.Data()));
541 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
542 fQARefUri = qam->GetDefaultStorage()->GetURI();
544 if (fQARefUri.Length() > 0) {
545 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
546 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
547 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
549 fQARefUri="local://$ALICE_ROOT/QAref";
550 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
551 AliWarning("Default QA refeference storage not yet set !!!!");
552 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
553 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
556 qam->SetDefaultStorage(fQARefUri);
560 qam->SetActiveDetectors(fQADetectors) ;
561 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) {
562 qam->SetCycleLength(AliQAv1::DETECTORINDEX_t(det), fQACycles[det]) ;
563 qam->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
565 if (!fRawReader && !fInput && fQATasks.Contains(AliQAv1::kRAWS))
566 fQATasks.ReplaceAll(Form("%d",AliQAv1::kRAWS), "") ;
567 qam->SetTasks(fQATasks) ;
568 qam->InitQADataMaker(AliCDBManager::Instance()->GetRun()) ;
571 Bool_t sameCycle = kFALSE ;
572 AliQADataMaker *qadm = qam->GetQADataMaker(AliQAv1::kGLOBAL);
573 AliInfo(Form("Initializing the global QA data maker"));
574 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS))) {
575 qadm->StartOfCycle(AliQAv1::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
576 TObjArray **arr=qadm->Init(AliQAv1::kRECPOINTS);
577 AliTracker::SetResidualsArray(arr);
580 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS))) {
581 qadm->StartOfCycle(AliQAv1::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
582 qadm->Init(AliQAv1::kESDS);
585 AliSysInfo::AddStamp("InitQA") ;
588 //_____________________________________________________________________________
589 void AliReconstruction::MergeQA(const char *fileName)
591 //Initialize the QA and start of cycle
592 AliCodeTimerAuto("") ;
593 AliQAManager::QAManager()->Merge(AliCDBManager::Instance()->GetRun(),fileName) ;
594 AliSysInfo::AddStamp("MergeQA") ;
597 //_____________________________________________________________________________
598 void AliReconstruction::InitCDB()
600 // activate a default CDB storage
601 // First check if we have any CDB storage set, because it is used
602 // to retrieve the calibration and alignment constants
603 AliCodeTimerAuto("");
605 if (fInitCDBCalled) return;
606 fInitCDBCalled = kTRUE;
608 AliCDBManager* man = AliCDBManager::Instance();
609 if (man->IsDefaultStorageSet())
611 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
612 AliWarning("Default CDB storage has been already set !");
613 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
614 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
615 fCDBUri = man->GetDefaultStorage()->GetURI();
618 if (fCDBUri.Length() > 0)
620 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
621 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
622 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
624 fCDBUri="local://$ALICE_ROOT/OCDB";
625 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
626 AliWarning("Default CDB storage not yet set !!!!");
627 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
628 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
631 man->SetDefaultStorage(fCDBUri);
634 // Now activate the detector specific CDB storage locations
635 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
636 TObject* obj = fSpecCDBUri[i];
638 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
639 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
640 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
641 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
643 AliSysInfo::AddStamp("InitCDB");
646 //_____________________________________________________________________________
647 void AliReconstruction::SetDefaultStorage(const char* uri) {
648 // Store the desired default CDB storage location
649 // Activate it later within the Run() method
655 //_____________________________________________________________________________
656 void AliReconstruction::SetQARefDefaultStorage(const char* uri) {
657 // Store the desired default CDB storage location
658 // Activate it later within the Run() method
661 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
664 //_____________________________________________________________________________
665 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
666 // Store a detector-specific CDB storage location
667 // Activate it later within the Run() method
669 AliCDBPath aPath(calibType);
670 if(!aPath.IsValid()){
671 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
672 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
673 if(!strcmp(calibType, fgkDetectorName[iDet])) {
674 aPath.SetPath(Form("%s/*", calibType));
675 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
679 if(!aPath.IsValid()){
680 AliError(Form("Not a valid path or detector: %s", calibType));
685 // // check that calibType refers to a "valid" detector name
686 // Bool_t isDetector = kFALSE;
687 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
688 // TString detName = fgkDetectorName[iDet];
689 // if(aPath.GetLevel0() == detName) {
690 // isDetector = kTRUE;
696 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
700 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
701 if (obj) fSpecCDBUri.Remove(obj);
702 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
706 //_____________________________________________________________________________
707 Bool_t AliReconstruction::SetRunNumberFromData()
709 // The method is called in Run() in order
710 // to set a correct run number.
711 // In case of raw data reconstruction the
712 // run number is taken from the raw data header
714 if (fSetRunNumberFromDataCalled) return kTRUE;
715 fSetRunNumberFromDataCalled = kTRUE;
717 AliCDBManager* man = AliCDBManager::Instance();
720 if(fRawReader->NextEvent()) {
721 if(man->GetRun() > 0) {
722 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
724 man->SetRun(fRawReader->GetRunNumber());
725 fRawReader->RewindEvents();
728 if(man->GetRun() > 0) {
729 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
732 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
738 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
740 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
745 // read run number from gAlice
746 if(rl->GetHeader()) {
747 man->SetRun(rl->GetHeader()->GetRun());
752 AliError("Neither run-loader header nor RawReader objects are found !");
764 //_____________________________________________________________________________
765 void AliReconstruction::SetCDBLock() {
766 // Set CDB lock: from now on it is forbidden to reset the run number
767 // or the default storage or to activate any further storage!
769 AliCDBManager::Instance()->SetLock(1);
772 //_____________________________________________________________________________
773 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
775 // Read the alignment objects from CDB.
776 // Each detector is supposed to have the
777 // alignment objects in DET/Align/Data CDB path.
778 // All the detector objects are then collected,
779 // sorted by geometry level (starting from ALIC) and
780 // then applied to the TGeo geometry.
781 // Finally an overlaps check is performed.
783 // Load alignment data from CDB and fill fAlignObjArray
784 if(fLoadAlignFromCDB){
786 TString detStr = detectors;
787 TString loadAlObjsListOfDets = "";
789 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
790 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
791 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
793 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
795 loadAlObjsListOfDets += fgkDetectorName[iDet];
796 loadAlObjsListOfDets += " ";
798 } // end loop over detectors
800 if(AliGeomManager::GetNalignable("GRP") != 0)
801 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
802 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
803 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
805 // Check if the array with alignment objects was
806 // provided by the user. If yes, apply the objects
807 // to the present TGeo geometry
808 if (fAlignObjArray) {
809 if (gGeoManager && gGeoManager->IsClosed()) {
810 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
811 AliError("The misalignment of one or more volumes failed!"
812 "Compare the list of simulated detectors and the list of detector alignment data!");
817 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
823 if (fAlignObjArray) {
824 fAlignObjArray->Delete();
825 delete fAlignObjArray; fAlignObjArray=NULL;
831 //_____________________________________________________________________________
832 void AliReconstruction::SetGAliceFile(const char* fileName)
834 // set the name of the galice file
836 fGAliceFileName = fileName;
839 //_____________________________________________________________________________
840 void AliReconstruction::SetInput(const char* input)
842 // In case the input string starts with 'mem://', we run in an online mode
843 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
844 // file is assumed. One can give as an input:
845 // mem://: - events taken from DAQ monitoring libs online
847 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
848 if (input) fRawInput = input;
851 //_____________________________________________________________________________
852 void AliReconstruction::SetOutput(const char* output)
854 // Set the output ESD filename
855 // 'output' is a normalt ROOT url
856 // The method is used in case of raw-data reco with PROOF
857 if (output) fESDOutput.SetUrl(output);
860 //_____________________________________________________________________________
861 void AliReconstruction::SetOption(const char* detector, const char* option)
863 // set options for the reconstruction of a detector
865 TObject* obj = fOptions.FindObject(detector);
866 if (obj) fOptions.Remove(obj);
867 fOptions.Add(new TNamed(detector, option));
870 //_____________________________________________________________________________
871 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
873 // Set custom reconstruction parameters for a given detector
874 // Single set of parameters for all the events
876 // First check if the reco-params are global
877 if(!strcmp(detector, "GRP")) {
879 fRecoParam.AddDetRecoParam(kNDetectors,par);
883 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
884 if(!strcmp(detector, fgkDetectorName[iDet])) {
886 fRecoParam.AddDetRecoParam(iDet,par);
893 //_____________________________________________________________________________
894 Bool_t AliReconstruction::SetFieldMap(Float_t l3Cur, Float_t diCur, Float_t l3Pol,
895 Float_t diPol, Float_t beamenergy,
896 const Char_t *beamtype, const Char_t *path)
898 //------------------------------------------------
899 // The magnetic field map, defined externally...
900 // L3 current 30000 A -> 0.5 T
901 // L3 current 12000 A -> 0.2 T
902 // dipole current 6000 A
903 // The polarities must be the same
904 //------------------------------------------------
905 const Float_t l3NominalCurrent1=30000.; // (A)
906 const Float_t l3NominalCurrent2=12000.; // (A)
907 const Float_t diNominalCurrent =6000. ; // (A)
909 const Float_t tolerance=0.03; // relative current tolerance
910 const Float_t zero=77.; // "zero" current (A)
912 TString s=(l3Pol < 0) ? "L3: -" : "L3: +";
914 AliMagF::BMap_t map = AliMagF::k5kG;
918 l3Cur = TMath::Abs(l3Cur);
919 if (TMath::Abs(l3Cur-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
920 fcL3 = l3Cur/l3NominalCurrent1;
923 } else if (TMath::Abs(l3Cur-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
924 fcL3 = l3Cur/l3NominalCurrent2;
927 } else if (l3Cur <= zero) {
929 map = AliMagF::k5kGUniform;
931 fUniformField=kTRUE; // track with the uniform (zero) B field
933 AliError(Form("Wrong L3 current (%f A)!",l3Cur));
937 diCur = TMath::Abs(diCur);
938 if (TMath::Abs(diCur-diNominalCurrent)/diNominalCurrent < tolerance) {
939 // 3% current tolerance...
940 fcDip = diCur/diNominalCurrent;
942 } else if (diCur <= zero) { // some small current..
946 AliError(Form("Wrong dipole current (%f A)!",diCur));
950 if (l3Pol!=diPol && (map==AliMagF::k5kG || map==AliMagF::k2kG) && fcDip!=0) {
951 AliError("L3 and Dipole polarities must be the same");
955 if (l3Pol<0) fcL3 = -fcL3;
956 if (diPol<0) fcDip = -fcDip;
958 AliMagF::BeamType_t btype = AliMagF::kNoBeamField;
959 TString btypestr = beamtype;
961 TPRegexp protonBeam("(proton|p)\\s*-?\\s*\\1");
962 TPRegexp ionBeam("(lead|pb|ion|a)\\s*-?\\s*\\1");
963 if (btypestr.Contains(ionBeam)) btype = AliMagF::kBeamTypeAA;
964 else if (btypestr.Contains(protonBeam)) btype = AliMagF::kBeamTypepp;
966 AliInfo(Form("Cannot determine the beam type from %s, assume no LHC magnet field",beamtype));
969 AliMagF* fld = new AliMagF("MagneticFieldMap", s.Data(), 2, fcL3, fcDip, 10., map, path,
971 TGeoGlobalMagField::Instance()->SetField( fld );
972 TGeoGlobalMagField::Instance()->Lock();
978 Bool_t AliReconstruction::InitGRP() {
979 //------------------------------------
980 // Initialization of the GRP entry
981 //------------------------------------
982 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
986 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
989 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
991 fGRPData = new AliGRPObject();
992 fGRPData->ReadValuesFromMap(m);
996 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
997 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
1001 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
1005 AliError("No GRP entry found in OCDB!");
1009 TString lhcState = fGRPData->GetLHCState();
1010 if (lhcState==AliGRPObject::GetInvalidString()) {
1011 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
1012 lhcState = "UNKNOWN";
1015 TString beamType = fGRPData->GetBeamType();
1016 if (beamType==AliGRPObject::GetInvalidString()) {
1017 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
1018 beamType = "UNKNOWN";
1021 Float_t beamEnergy = fGRPData->GetBeamEnergy();
1022 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
1023 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
1026 // energy is provided in MeV*120
1027 beamEnergy /= 120E3;
1029 TString runType = fGRPData->GetRunType();
1030 if (runType==AliGRPObject::GetInvalidString()) {
1031 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
1032 runType = "UNKNOWN";
1035 Int_t activeDetectors = fGRPData->GetDetectorMask();
1036 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
1037 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
1038 activeDetectors = 1074790399;
1041 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1045 // Process the list of active detectors
1046 if (activeDetectors) {
1047 UInt_t detMask = activeDetectors;
1048 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1049 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1050 fFillESD = MatchDetectorList(fFillESD,detMask);
1051 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1052 fLoadCDB.Form("%s %s %s %s",
1053 fRunLocalReconstruction.Data(),
1054 fRunTracking.Data(),
1056 fQADetectors.Data());
1057 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1058 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1)) {
1059 // switch off the vertexer
1060 AliInfo("SPD is not in the list of active detectors. Vertexer switched off.");
1061 fRunVertexFinder = kFALSE;
1063 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1064 // switch off the reading of CTP raw-data payload
1065 if (fFillTriggerESD) {
1066 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1067 fFillTriggerESD = kFALSE;
1072 AliInfo("===================================================================================");
1073 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1074 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1075 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1076 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1077 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1078 AliInfo("===================================================================================");
1080 //*** Dealing with the magnetic field map
1081 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {AliInfo("Running with the externally locked B field !");}
1083 // Construct the field map out of the information retrieved from GRP.
1086 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1087 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1088 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1092 Char_t l3Polarity = fGRPData->GetL3Polarity();
1093 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1094 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1099 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1100 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1101 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1105 Char_t diPolarity = fGRPData->GetDipolePolarity();
1106 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1107 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1112 TObjString *l3Current=
1113 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
1115 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1118 TObjString *l3Polarity=
1119 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
1121 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1126 TObjString *diCurrent=
1127 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
1129 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1132 TObjString *diPolarity=
1133 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
1135 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1141 if ( !SetFieldMap(l3Current, diCurrent, l3Polarity ? -1:1, diPolarity ? -1:1) )
1142 AliFatal("Failed to creat a B field map ! Exiting...");
1143 AliInfo("Running with the B field constructed out of GRP !");
1145 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1149 //*** Get the diamond profiles from OCDB
1150 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1152 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1154 AliError("No SPD diamond profile found in OCDB!");
1157 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1159 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1161 AliError("No diamond profile found in OCDB!");
1164 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1166 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1168 AliError("No TPC diamond profile found in OCDB!");
1174 //_____________________________________________________________________________
1175 Bool_t AliReconstruction::LoadCDB()
1177 AliCodeTimerAuto("");
1179 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1181 TString detStr = fLoadCDB;
1182 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1183 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1184 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1189 //_____________________________________________________________________________
1190 Bool_t AliReconstruction::Run(const char* input)
1193 AliCodeTimerAuto("");
1196 if (GetAbort() != TSelector::kContinue) return kFALSE;
1198 TChain *chain = NULL;
1199 if (fRawReader && (chain = fRawReader->GetChain())) {
1202 gProof->AddInput(this);
1203 if (!fESDOutput.IsValid()) {
1204 fESDOutput.SetProtocol("root",kTRUE);
1205 fESDOutput.SetHost(gSystem->HostName());
1206 fESDOutput.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1208 AliInfo(Form("Output file with ESDs is %s",fESDOutput.GetUrl()));
1209 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",fESDOutput.GetUrl()));
1210 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1212 chain->Process("AliReconstruction");
1215 chain->Process(this);
1220 if (GetAbort() != TSelector::kContinue) return kFALSE;
1222 if (GetAbort() != TSelector::kContinue) return kFALSE;
1223 //******* The loop over events
1224 AliInfo("Starting looping over events");
1226 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1227 (fRawReader && fRawReader->NextEvent())) {
1228 if (!ProcessEvent(iEvent)) {
1229 Abort("ProcessEvent",TSelector::kAbortFile);
1235 if (GetAbort() != TSelector::kContinue) return kFALSE;
1237 if (GetAbort() != TSelector::kContinue) return kFALSE;
1243 //_____________________________________________________________________________
1244 void AliReconstruction::InitRawReader(const char* input)
1246 AliCodeTimerAuto("");
1248 // Init raw-reader and
1249 // set the input in case of raw data
1250 if (input) fRawInput = input;
1251 fRawReader = AliRawReader::Create(fRawInput.Data());
1253 AliInfo("Reconstruction will run over digits");
1255 if (!fEquipIdMap.IsNull() && fRawReader)
1256 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1258 if (!fUseHLTData.IsNull()) {
1259 // create the RawReaderHLT which performs redirection of HLT input data for
1260 // the specified detectors
1261 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1263 fParentRawReader=fRawReader;
1264 fRawReader=pRawReader;
1266 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1269 AliSysInfo::AddStamp("CreateRawReader");
1272 //_____________________________________________________________________________
1273 void AliReconstruction::InitRun(const char* input)
1275 // Initialization of raw-reader,
1276 // run number, CDB etc.
1277 AliCodeTimerAuto("");
1278 AliSysInfo::AddStamp("Start");
1280 // Initialize raw-reader if any
1281 InitRawReader(input);
1283 // Initialize the CDB storage
1286 // Set run number in CDBManager (if it is not already set by the user)
1287 if (!SetRunNumberFromData()) {
1288 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1292 // Set CDB lock: from now on it is forbidden to reset the run number
1293 // or the default storage or to activate any further storage!
1298 //_____________________________________________________________________________
1299 void AliReconstruction::Begin(TTree *)
1301 // Initialize AlReconstruction before
1302 // going into the event loop
1303 // Should follow the TSelector convention
1304 // i.e. initialize only the object on the client side
1305 AliCodeTimerAuto("");
1307 AliReconstruction *reco = NULL;
1309 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1312 AliSysInfo::AddStamp("ReadInputInBegin");
1315 // Import ideal TGeo geometry and apply misalignment
1317 TString geom(gSystem->DirName(fGAliceFileName));
1318 geom += "/geometry.root";
1319 AliGeomManager::LoadGeometry(geom.Data());
1321 Abort("LoadGeometry", TSelector::kAbortProcess);
1324 AliSysInfo::AddStamp("LoadGeom");
1325 TString detsToCheck=fRunLocalReconstruction;
1326 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1327 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1330 AliSysInfo::AddStamp("CheckGeom");
1333 if (!MisalignGeometry(fLoadAlignData)) {
1334 Abort("MisalignGeometry", TSelector::kAbortProcess);
1337 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1338 AliSysInfo::AddStamp("MisalignGeom");
1341 Abort("InitGRP", TSelector::kAbortProcess);
1344 AliSysInfo::AddStamp("InitGRP");
1347 Abort("LoadCDB", TSelector::kAbortProcess);
1350 AliSysInfo::AddStamp("LoadCDB");
1352 // Read the reconstruction parameters from OCDB
1353 if (!InitRecoParams()) {
1354 AliWarning("Not all detectors have correct RecoParam objects initialized");
1356 AliSysInfo::AddStamp("InitRecoParams");
1358 if (fInput && gProof) {
1359 if (reco) *reco = *this;
1361 gProof->AddInputData(gGeoManager,kTRUE);
1363 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1364 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1365 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1366 magFieldMap->SetName("MagneticFieldMap");
1367 gProof->AddInputData(magFieldMap,kTRUE);
1372 //_____________________________________________________________________________
1373 void AliReconstruction::SlaveBegin(TTree*)
1375 // Initialization related to run-loader,
1376 // vertexer, trackers, recontructors
1377 // In proof mode it is executed on the slave
1378 AliCodeTimerAuto("");
1380 TProofOutputFile *outProofFile = NULL;
1382 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1385 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1387 AliGeomManager::SetGeometry(tgeo);
1389 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1390 Int_t runNumber = -1;
1391 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1392 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1393 man->SetCacheFlag(kTRUE);
1394 man->SetLock(kTRUE);
1398 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1399 TGeoGlobalMagField::Instance()->SetField(map);
1401 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1402 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1403 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1404 fOutput->Add(outProofFile);
1406 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1409 // get the run loader
1410 if (!InitRunLoader()) {
1411 Abort("InitRunLoader", TSelector::kAbortProcess);
1414 AliSysInfo::AddStamp("LoadLoader");
1416 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1419 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1420 Abort("CreateTrackers", TSelector::kAbortProcess);
1423 AliSysInfo::AddStamp("CreateTrackers");
1425 // create the ESD output file and tree
1426 if (!outProofFile) {
1427 ffile = TFile::Open("AliESDs.root", "RECREATE");
1428 ffile->SetCompressionLevel(2);
1429 if (!ffile->IsOpen()) {
1430 Abort("OpenESDFile", TSelector::kAbortProcess);
1435 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1436 Abort(Form("Problems opening output PROOF file: %s/%s",
1437 outProofFile->GetDir(), outProofFile->GetFileName()),
1438 TSelector::kAbortProcess);
1443 ftree = new TTree("esdTree", "Tree with ESD objects");
1444 fesd = new AliESDEvent();
1445 fesd->CreateStdContent();
1447 fesd->WriteToTree(ftree);
1448 if (fWriteESDfriend) {
1450 // Since we add the branch manually we must
1451 // book and add it after WriteToTree
1452 // otherwise it is created twice,
1453 // once via writetotree and once here.
1454 // The case for AliESDfriend is now
1455 // caught also in AlIESDEvent::WriteToTree but
1456 // be careful when changing the name (AliESDfriend is not
1457 // a TNamed so we had to hardwire it)
1458 fesdf = new AliESDfriend();
1459 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1460 br->SetFile("AliESDfriends.root");
1461 fesd->AddObject(fesdf);
1463 ftree->GetUserInfo()->Add(fesd);
1465 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1466 fhltesd = new AliESDEvent();
1467 fhltesd->CreateStdContent();
1469 // read the ESD template from CDB
1470 // HLT is allowed to put non-std content to its ESD, the non-std
1471 // objects need to be created before invocation of WriteToTree in
1472 // order to create all branches. Initialization is done from an
1473 // ESD layout template in CDB
1474 AliCDBManager* man = AliCDBManager::Instance();
1475 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1476 AliCDBEntry* hltESDConfig=NULL;
1477 if (man->GetId(hltESDConfigPath)!=NULL &&
1478 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1479 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1481 // init all internal variables from the list of objects
1482 pESDLayout->GetStdContent();
1484 // copy content and create non-std objects
1485 *fhltesd=*pESDLayout;
1488 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1489 hltESDConfigPath.GetPath().Data()));
1493 fhltesd->WriteToTree(fhlttree);
1494 fhlttree->GetUserInfo()->Add(fhltesd);
1496 ProcInfo_t procInfo;
1497 gSystem->GetProcInfo(&procInfo);
1498 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1501 //Initialize the QA and start of cycle
1502 if (fRunQA || fRunGlobalQA)
1505 //Initialize the Plane Efficiency framework
1506 if (fRunPlaneEff && !InitPlaneEff()) {
1507 Abort("InitPlaneEff", TSelector::kAbortProcess);
1511 if (strcmp(gProgName,"alieve") == 0)
1512 fRunAliEVE = InitAliEVE();
1517 //_____________________________________________________________________________
1518 Bool_t AliReconstruction::Process(Long64_t entry)
1520 // run the reconstruction over a single entry
1521 // from the chain with raw data
1522 AliCodeTimerAuto("");
1524 TTree *currTree = fChain->GetTree();
1525 AliRawVEvent *event = NULL;
1526 currTree->SetBranchAddress("rawevent",&event);
1527 currTree->GetEntry(entry);
1528 fRawReader = new AliRawReaderRoot(event);
1529 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1537 //_____________________________________________________________________________
1538 void AliReconstruction::Init(TTree *tree)
1541 AliError("The input tree is not found!");
1547 //_____________________________________________________________________________
1548 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1550 // run the reconstruction over a single event
1551 // The event loop is steered in Run method
1553 AliCodeTimerAuto("");
1555 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1556 fRunLoader->SetEventNumber(iEvent);
1557 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1559 fRunLoader->TreeE()->Fill();
1560 if (fRawReader && fRawReader->UseAutoSaveESD())
1561 fRunLoader->TreeE()->AutoSave("SaveSelf");
1564 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1568 AliInfo(Form("processing event %d", iEvent));
1570 fRunLoader->GetEvent(iEvent);
1572 // Fill Event-info object
1574 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo);
1575 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1577 // Set the reco-params
1579 TString detStr = fLoadCDB;
1580 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1581 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1582 AliReconstructor *reconstructor = GetReconstructor(iDet);
1583 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1584 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1585 reconstructor->SetRecoParam(par);
1587 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1595 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1596 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1598 // local single event reconstruction
1599 if (!fRunLocalReconstruction.IsNull()) {
1600 TString detectors=fRunLocalReconstruction;
1601 // run HLT event reconstruction first
1602 // ;-( IsSelected changes the string
1603 if (IsSelected("HLT", detectors) &&
1604 !RunLocalEventReconstruction("HLT")) {
1605 if (fStopOnError) {CleanUp(); return kFALSE;}
1607 detectors=fRunLocalReconstruction;
1608 detectors.ReplaceAll("HLT", "");
1609 if (!RunLocalEventReconstruction(detectors)) {
1610 if (fStopOnError) {CleanUp(); return kFALSE;}
1614 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1615 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1616 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1617 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1619 // Set magnetic field from the tracker
1620 fesd->SetMagneticField(AliTracker::GetBz());
1621 fhltesd->SetMagneticField(AliTracker::GetBz());
1623 // Set most probable pt, for B=0 tracking
1624 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1625 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1626 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1628 // Fill raw-data error log into the ESD
1629 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1632 if (fRunVertexFinder) {
1633 if (!RunVertexFinder(fesd)) {
1634 if (fStopOnError) {CleanUp(); return kFALSE;}
1638 // For Plane Efficiency: run the SPD trackleter
1639 if (fRunPlaneEff && fSPDTrackleter) {
1640 if (!RunSPDTrackleting(fesd)) {
1641 if (fStopOnError) {CleanUp(); return kFALSE;}
1646 if (!fRunTracking.IsNull()) {
1647 if (fRunMuonTracking) {
1648 if (!RunMuonTracking(fesd)) {
1649 if (fStopOnError) {CleanUp(); return kFALSE;}
1655 if (!fRunTracking.IsNull()) {
1656 if (!RunTracking(fesd)) {
1657 if (fStopOnError) {CleanUp(); return kFALSE;}
1662 if (!fFillESD.IsNull()) {
1663 TString detectors=fFillESD;
1664 // run HLT first and on hltesd
1665 // ;-( IsSelected changes the string
1666 if (IsSelected("HLT", detectors) &&
1667 !FillESD(fhltesd, "HLT")) {
1668 if (fStopOnError) {CleanUp(); return kFALSE;}
1671 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1672 if (detectors.Contains("ALL")) {
1674 for (Int_t idet=0; idet<kNDetectors; ++idet){
1675 detectors += fgkDetectorName[idet];
1679 detectors.ReplaceAll("HLT", "");
1680 if (!FillESD(fesd, detectors)) {
1681 if (fStopOnError) {CleanUp(); return kFALSE;}
1685 // fill Event header information from the RawEventHeader
1686 if (fRawReader){FillRawEventHeaderESD(fesd);}
1689 AliESDpid::MakePID(fesd);
1691 if (fFillTriggerESD) {
1692 if (!FillTriggerESD(fesd)) {
1693 if (fStopOnError) {CleanUp(); return kFALSE;}
1700 // Propagate track to the beam pipe (if not already done by ITS)
1702 const Int_t ntracks = fesd->GetNumberOfTracks();
1703 const Double_t kBz = fesd->GetMagneticField();
1704 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1707 UShort_t *selectedIdx=new UShort_t[ntracks];
1709 for (Int_t itrack=0; itrack<ntracks; itrack++){
1710 const Double_t kMaxStep = 1; //max step over the material
1713 AliESDtrack *track = fesd->GetTrack(itrack);
1714 if (!track) continue;
1716 AliExternalTrackParam *tpcTrack =
1717 (AliExternalTrackParam *)track->GetTPCInnerParam();
1721 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1724 Int_t n=trkArray.GetEntriesFast();
1725 selectedIdx[n]=track->GetID();
1726 trkArray.AddLast(tpcTrack);
1729 //Tracks refitted by ITS should already be at the SPD vertex
1730 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1733 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1734 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1739 // Improve the reconstructed primary vertex position using the tracks
1741 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1742 if(fesd->GetPrimaryVertexSPD()) {
1743 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1744 if(vtitle.Contains("cosmics")) {
1745 runVertexFinderTracks=kFALSE;
1749 if (runVertexFinderTracks) {
1750 // TPC + ITS primary vertex
1751 ftVertexer->SetITSMode();
1752 ftVertexer->SetConstraintOff();
1753 // get cuts for vertexer from AliGRPRecoParam
1755 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1756 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1757 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1758 ftVertexer->SetCuts(cutsVertexer);
1759 delete [] cutsVertexer; cutsVertexer = NULL;
1760 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1761 ftVertexer->SetVtxStart(fDiamondProfile);
1763 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1765 if (pvtx->GetStatus()) {
1766 fesd->SetPrimaryVertexTracks(pvtx);
1767 for (Int_t i=0; i<ntracks; i++) {
1768 AliESDtrack *t = fesd->GetTrack(i);
1769 t->RelateToVertex(pvtx, kBz, kVeryBig);
1774 // TPC-only primary vertex
1775 ftVertexer->SetTPCMode();
1776 ftVertexer->SetConstraintOff();
1777 // get cuts for vertexer from AliGRPRecoParam
1779 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1780 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1781 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1782 ftVertexer->SetCuts(cutsVertexer);
1783 delete [] cutsVertexer; cutsVertexer = NULL;
1784 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1785 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1787 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1789 if (pvtx->GetStatus()) {
1790 fesd->SetPrimaryVertexTPC(pvtx);
1791 for (Int_t i=0; i<ntracks; i++) {
1792 AliESDtrack *t = fesd->GetTrack(i);
1793 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1799 delete[] selectedIdx;
1801 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1806 AliV0vertexer vtxer;
1807 vtxer.Tracks2V0vertices(fesd);
1809 if (fRunCascadeFinder) {
1811 AliCascadeVertexer cvtxer;
1812 cvtxer.V0sTracks2CascadeVertices(fesd);
1817 if (fCleanESD) CleanESD(fesd);
1820 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1821 AliQAManager::QAManager()->RunOneEvent(fesd) ;
1824 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1825 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1826 if (qadm && fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1827 qadm->Exec(AliQAv1::kESDS, fesd);
1830 if (fWriteESDfriend) {
1831 // fesdf->~AliESDfriend();
1832 // new (fesdf) AliESDfriend(); // Reset...
1833 fesd->GetESDfriend(fesdf);
1837 // Auto-save the ESD tree in case of prompt reco @P2
1838 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1839 ftree->AutoSave("SaveSelf");
1840 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1841 if (friendfile) friendfile->Save();
1848 if (fRunAliEVE) RunAliEVE();
1852 if (fWriteESDfriend) {
1853 fesdf->~AliESDfriend();
1854 new (fesdf) AliESDfriend(); // Reset...
1857 ProcInfo_t procInfo;
1858 gSystem->GetProcInfo(&procInfo);
1859 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1862 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1863 if (fReconstructor[iDet])
1864 fReconstructor[iDet]->SetRecoParam(NULL);
1867 if (fRunQA || fRunGlobalQA)
1868 AliQAManager::QAManager()->Increment() ;
1873 //_____________________________________________________________________________
1874 void AliReconstruction::SlaveTerminate()
1876 // Finalize the run on the slave side
1877 // Called after the exit
1878 // from the event loop
1879 AliCodeTimerAuto("");
1881 if (fIsNewRunLoader) { // galice.root didn't exist
1882 fRunLoader->WriteHeader("OVERWRITE");
1883 fRunLoader->CdGAFile();
1884 fRunLoader->Write(0, TObject::kOverwrite);
1887 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1888 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1890 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1891 cdbMapCopy->SetOwner(1);
1892 cdbMapCopy->SetName("cdbMap");
1893 TIter iter(cdbMap->GetTable());
1896 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1897 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1898 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1899 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1902 TList *cdbListCopy = new TList();
1903 cdbListCopy->SetOwner(1);
1904 cdbListCopy->SetName("cdbList");
1906 TIter iter2(cdbList);
1909 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1910 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1913 ftree->GetUserInfo()->Add(cdbMapCopy);
1914 ftree->GetUserInfo()->Add(cdbListCopy);
1919 if (fWriteESDfriend)
1920 ftree->SetBranchStatus("ESDfriend*",0);
1921 // we want to have only one tree version number
1922 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1925 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1926 if (fRunPlaneEff && !FinishPlaneEff()) {
1927 AliWarning("Finish PlaneEff evaluation failed");
1930 // End of cycle for the in-loop
1932 AliQAManager::QAManager()->EndOfCycle() ;
1935 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1937 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS)))
1938 qadm->EndOfCycle(AliQAv1::kRECPOINTS);
1939 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1940 qadm->EndOfCycle(AliQAv1::kESDS);
1945 if (fRunQA || fRunGlobalQA) {
1947 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1948 TString qaOutputFile = outputFileName->GetTitle();
1949 qaOutputFile.ReplaceAll(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()),
1950 Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
1951 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",AliQAv1::GetQADataFileName()));
1952 qaProofFile->SetOutputFileName(qaOutputFile.Data());
1953 fOutput->Add(qaProofFile);
1954 MergeQA(qaProofFile->GetFileName());
1966 //_____________________________________________________________________________
1967 void AliReconstruction::Terminate()
1969 // Create tags for the events in the ESD tree (the ESD tree is always present)
1970 // In case of empty events the tags will contain dummy values
1971 AliCodeTimerAuto("");
1973 // Do not call the ESD tag creator in case of PROOF-based reconstruction
1975 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1976 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
1979 // Cleanup of CDB manager: cache and active storages!
1980 AliCDBManager::Instance()->ClearCache();
1983 //_____________________________________________________________________________
1984 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1986 // run the local reconstruction
1988 static Int_t eventNr=0;
1989 AliCodeTimerAuto("")
1991 TString detStr = detectors;
1992 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1993 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1994 AliReconstructor* reconstructor = GetReconstructor(iDet);
1995 if (!reconstructor) continue;
1996 AliLoader* loader = fLoader[iDet];
1997 // Matthias April 2008: temporary fix to run HLT reconstruction
1998 // although the HLT loader is missing
1999 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2001 reconstructor->Reconstruct(fRawReader, NULL);
2004 reconstructor->Reconstruct(dummy, NULL);
2009 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2012 // conversion of digits
2013 if (fRawReader && reconstructor->HasDigitConversion()) {
2014 AliInfo(Form("converting raw data digits into root objects for %s",
2015 fgkDetectorName[iDet]));
2016 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2017 // fgkDetectorName[iDet]));
2018 loader->LoadDigits("update");
2019 loader->CleanDigits();
2020 loader->MakeDigitsContainer();
2021 TTree* digitsTree = loader->TreeD();
2022 reconstructor->ConvertDigits(fRawReader, digitsTree);
2023 loader->WriteDigits("OVERWRITE");
2024 loader->UnloadDigits();
2026 // local reconstruction
2027 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2028 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2029 loader->LoadRecPoints("update");
2030 loader->CleanRecPoints();
2031 loader->MakeRecPointsContainer();
2032 TTree* clustersTree = loader->TreeR();
2033 if (fRawReader && !reconstructor->HasDigitConversion()) {
2034 reconstructor->Reconstruct(fRawReader, clustersTree);
2036 loader->LoadDigits("read");
2037 TTree* digitsTree = loader->TreeD();
2039 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2040 if (fStopOnError) return kFALSE;
2042 reconstructor->Reconstruct(digitsTree, clustersTree);
2044 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2045 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2048 loader->UnloadDigits();
2051 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2052 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2054 loader->WriteRecPoints("OVERWRITE");
2055 loader->UnloadRecPoints();
2056 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2058 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2059 AliError(Form("the following detectors were not found: %s",
2061 if (fStopOnError) return kFALSE;
2066 //_____________________________________________________________________________
2067 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2069 // run the SPD trackleting (for SPD efficiency purpouses)
2071 AliCodeTimerAuto("")
2073 Double_t vtxPos[3] = {0, 0, 0};
2074 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2076 TArrayF mcVertex(3);
2078 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2079 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2080 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2083 const AliESDVertex *vertex = esd->GetVertex();
2085 AliWarning("Vertex not found");
2088 vertex->GetXYZ(vtxPos);
2089 vertex->GetSigmaXYZ(vtxErr);
2090 if (fSPDTrackleter) {
2091 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2094 fLoader[0]->LoadRecPoints("read");
2095 TTree* tree = fLoader[0]->TreeR();
2097 AliError("Can't get the ITS cluster tree");
2100 fSPDTrackleter->LoadClusters(tree);
2101 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2103 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2104 AliError("AliITSTrackleterSPDEff Clusters2Tracks failed");
2105 // fLoader[0]->UnloadRecPoints();
2108 //fSPDTrackleter->UnloadRecPoints();
2110 AliWarning("SPDTrackleter not available");
2116 //_____________________________________________________________________________
2117 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2119 // run the barrel tracking
2121 AliCodeTimerAuto("")
2123 AliVertexer *vertexer = CreateVertexer();
2124 if (!vertexer) return kFALSE;
2126 AliInfo("running the ITS vertex finder");
2127 AliESDVertex* vertex = NULL;
2129 fLoader[0]->LoadRecPoints();
2130 TTree* cltree = fLoader[0]->TreeR();
2132 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2133 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2136 AliError("Can't get the ITS cluster tree");
2138 fLoader[0]->UnloadRecPoints();
2141 AliError("Can't get the ITS loader");
2144 AliWarning("Vertex not found");
2145 vertex = new AliESDVertex();
2146 vertex->SetName("default");
2149 vertex->SetName("reconstructed");
2154 vertex->GetXYZ(vtxPos);
2155 vertex->GetSigmaXYZ(vtxErr);
2157 esd->SetPrimaryVertexSPD(vertex);
2158 AliESDVertex *vpileup = NULL;
2159 Int_t novertices = 0;
2160 vpileup = vertexer->GetAllVertices(novertices);
2162 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2164 // if SPD multiplicity has been determined, it is stored in the ESD
2165 AliMultiplicity *mult = vertexer->GetMultiplicity();
2166 if(mult)esd->SetMultiplicity(mult);
2168 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2169 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2178 //_____________________________________________________________________________
2179 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2181 // run the HLT barrel tracking
2183 AliCodeTimerAuto("")
2186 AliError("Missing runLoader!");
2190 AliInfo("running HLT tracking");
2192 // Get a pointer to the HLT reconstructor
2193 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2194 if (!reconstructor) return kFALSE;
2197 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2198 TString detName = fgkDetectorName[iDet];
2199 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2200 reconstructor->SetOption(detName.Data());
2201 AliTracker *tracker = reconstructor->CreateTracker();
2203 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2204 if (fStopOnError) return kFALSE;
2208 Double_t vtxErr[3]={0.005,0.005,0.010};
2209 const AliESDVertex *vertex = esd->GetVertex();
2210 vertex->GetXYZ(vtxPos);
2211 tracker->SetVertex(vtxPos,vtxErr);
2213 fLoader[iDet]->LoadRecPoints("read");
2214 TTree* tree = fLoader[iDet]->TreeR();
2216 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2219 tracker->LoadClusters(tree);
2221 if (tracker->Clusters2Tracks(esd) != 0) {
2222 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2226 tracker->UnloadClusters();
2234 //_____________________________________________________________________________
2235 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2237 // run the muon spectrometer tracking
2239 AliCodeTimerAuto("")
2242 AliError("Missing runLoader!");
2245 Int_t iDet = 7; // for MUON
2247 AliInfo("is running...");
2249 // Get a pointer to the MUON reconstructor
2250 AliReconstructor *reconstructor = GetReconstructor(iDet);
2251 if (!reconstructor) return kFALSE;
2254 TString detName = fgkDetectorName[iDet];
2255 AliDebug(1, Form("%s tracking", detName.Data()));
2256 AliTracker *tracker = reconstructor->CreateTracker();
2258 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2263 fLoader[iDet]->LoadRecPoints("read");
2265 tracker->LoadClusters(fLoader[iDet]->TreeR());
2267 Int_t rv = tracker->Clusters2Tracks(esd);
2271 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2275 fLoader[iDet]->UnloadRecPoints();
2277 tracker->UnloadClusters();
2285 //_____________________________________________________________________________
2286 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2288 // run the barrel tracking
2289 static Int_t eventNr=0;
2290 AliCodeTimerAuto("")
2292 AliInfo("running tracking");
2295 //Fill the ESD with the T0 info (will be used by the TOF)
2296 if (fReconstructor[11] && fLoader[11]) {
2297 fLoader[11]->LoadRecPoints("READ");
2298 TTree *treeR = fLoader[11]->TreeR();
2300 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2304 // pass 1: TPC + ITS inwards
2305 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2306 if (!fTracker[iDet]) continue;
2307 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2310 fLoader[iDet]->LoadRecPoints("read");
2311 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2312 TTree* tree = fLoader[iDet]->TreeR();
2314 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2317 fTracker[iDet]->LoadClusters(tree);
2318 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2320 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2321 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2324 // preliminary PID in TPC needed by the ITS tracker
2326 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2327 AliESDpid::MakePID(esd);
2329 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2332 // pass 2: ALL backwards
2334 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2335 if (!fTracker[iDet]) continue;
2336 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2339 if (iDet > 1) { // all except ITS, TPC
2341 fLoader[iDet]->LoadRecPoints("read");
2342 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2343 tree = fLoader[iDet]->TreeR();
2345 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2348 fTracker[iDet]->LoadClusters(tree);
2349 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2353 if (iDet>1) // start filling residuals for the "outer" detectors
2355 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2356 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2357 if ( ! arr[fRecoParam.GetEventSpecie()]->At(0) ) {
2358 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2359 qadm->InitRecPointsForTracker() ;
2362 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2363 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2368 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2369 fTracker[iDet]->UnloadClusters();
2370 fLoader[iDet]->UnloadRecPoints();
2372 // updated PID in TPC needed by the ITS tracker -MI
2374 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2375 AliESDpid::MakePID(esd);
2377 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2379 //stop filling residuals for the "outer" detectors
2380 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2382 // pass 3: TRD + TPC + ITS refit inwards
2384 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2385 if (!fTracker[iDet]) continue;
2386 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2389 if (iDet<2) // start filling residuals for TPC and ITS
2391 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2392 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2393 if ( ! arr[fRecoParam.GetEventSpecie()]->At(0) ) {
2394 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2395 qadm->InitRecPointsForTracker() ;
2399 if (fTracker[iDet]->RefitInward(esd) != 0) {
2400 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2403 // run postprocessing
2404 if (fTracker[iDet]->PostProcess(esd) != 0) {
2405 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2408 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2411 // write space-points to the ESD in case alignment data output
2413 if (fWriteAlignmentData)
2414 WriteAlignmentData(esd);
2416 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2417 if (!fTracker[iDet]) continue;
2419 fTracker[iDet]->UnloadClusters();
2420 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2421 fLoader[iDet]->UnloadRecPoints();
2422 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2424 // stop filling residuals for TPC and ITS
2425 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2431 //_____________________________________________________________________________
2432 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2434 // Remove the data which are not needed for the physics analysis.
2437 Int_t nTracks=esd->GetNumberOfTracks();
2438 Int_t nV0s=esd->GetNumberOfV0s();
2440 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2442 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2443 Bool_t rc=esd->Clean(cleanPars);
2445 nTracks=esd->GetNumberOfTracks();
2446 nV0s=esd->GetNumberOfV0s();
2448 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2453 //_____________________________________________________________________________
2454 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2456 // fill the event summary data
2458 AliCodeTimerAuto("")
2459 static Int_t eventNr=0;
2460 TString detStr = detectors;
2462 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2463 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2464 AliReconstructor* reconstructor = GetReconstructor(iDet);
2465 if (!reconstructor) continue;
2466 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2467 TTree* clustersTree = NULL;
2468 if (fLoader[iDet]) {
2469 fLoader[iDet]->LoadRecPoints("read");
2470 clustersTree = fLoader[iDet]->TreeR();
2471 if (!clustersTree) {
2472 AliError(Form("Can't get the %s clusters tree",
2473 fgkDetectorName[iDet]));
2474 if (fStopOnError) return kFALSE;
2477 if (fRawReader && !reconstructor->HasDigitConversion()) {
2478 reconstructor->FillESD(fRawReader, clustersTree, esd);
2480 TTree* digitsTree = NULL;
2481 if (fLoader[iDet]) {
2482 fLoader[iDet]->LoadDigits("read");
2483 digitsTree = fLoader[iDet]->TreeD();
2485 AliError(Form("Can't get the %s digits tree",
2486 fgkDetectorName[iDet]));
2487 if (fStopOnError) return kFALSE;
2490 reconstructor->FillESD(digitsTree, clustersTree, esd);
2491 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2493 if (fLoader[iDet]) {
2494 fLoader[iDet]->UnloadRecPoints();
2498 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2499 AliError(Form("the following detectors were not found: %s",
2501 if (fStopOnError) return kFALSE;
2503 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2508 //_____________________________________________________________________________
2509 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2511 // Reads the trigger decision which is
2512 // stored in Trigger.root file and fills
2513 // the corresponding esd entries
2515 AliCodeTimerAuto("")
2517 AliInfo("Filling trigger information into the ESD");
2520 AliCTPRawStream input(fRawReader);
2521 if (!input.Next()) {
2522 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2525 if (esd->GetTriggerMask() != input.GetClassMask())
2526 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2527 input.GetClassMask(),esd->GetTriggerMask()));
2528 if (esd->GetOrbitNumber() != input.GetOrbitID())
2529 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2530 input.GetOrbitID(),esd->GetOrbitNumber()));
2531 if (esd->GetBunchCrossNumber() != input.GetBCID())
2532 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2533 input.GetBCID(),esd->GetBunchCrossNumber()));
2536 // Here one has to add the filling of trigger inputs and
2537 // interaction records
2547 //_____________________________________________________________________________
2548 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2551 // Filling information from RawReader Header
2554 if (!fRawReader) return kFALSE;
2556 AliInfo("Filling information from RawReader Header");
2558 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2559 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2560 esd->SetPeriodNumber(fRawReader->GetPeriod());
2562 esd->SetTimeStamp(fRawReader->GetTimestamp());
2563 esd->SetEventType(fRawReader->GetType());
2569 //_____________________________________________________________________________
2570 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2572 // check whether detName is contained in detectors
2573 // if yes, it is removed from detectors
2575 // check if all detectors are selected
2576 if ((detectors.CompareTo("ALL") == 0) ||
2577 detectors.BeginsWith("ALL ") ||
2578 detectors.EndsWith(" ALL") ||
2579 detectors.Contains(" ALL ")) {
2584 // search for the given detector
2585 Bool_t result = kFALSE;
2586 if ((detectors.CompareTo(detName) == 0) ||
2587 detectors.BeginsWith(detName+" ") ||
2588 detectors.EndsWith(" "+detName) ||
2589 detectors.Contains(" "+detName+" ")) {
2590 detectors.ReplaceAll(detName, "");
2594 // clean up the detectors string
2595 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2596 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2597 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2602 //_____________________________________________________________________________
2603 Bool_t AliReconstruction::InitRunLoader()
2605 // get or create the run loader
2607 if (gAlice) delete gAlice;
2610 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2611 // load all base libraries to get the loader classes
2612 TString libs = gSystem->GetLibraries();
2613 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2614 TString detName = fgkDetectorName[iDet];
2615 if (detName == "HLT") continue;
2616 if (libs.Contains("lib" + detName + "base.so")) continue;
2617 gSystem->Load("lib" + detName + "base.so");
2619 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2621 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2626 fRunLoader->CdGAFile();
2627 fRunLoader->LoadgAlice();
2629 //PH This is a temporary fix to give access to the kinematics
2630 //PH that is needed for the labels of ITS clusters
2631 fRunLoader->LoadHeader();
2632 fRunLoader->LoadKinematics();
2634 } else { // galice.root does not exist
2636 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2638 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2639 AliConfig::GetDefaultEventFolderName(),
2642 AliError(Form("could not create run loader in file %s",
2643 fGAliceFileName.Data()));
2647 fIsNewRunLoader = kTRUE;
2648 fRunLoader->MakeTree("E");
2650 if (fNumberOfEventsPerFile > 0)
2651 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2653 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2659 //_____________________________________________________________________________
2660 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2662 // get the reconstructor object and the loader for a detector
2664 if (fReconstructor[iDet]) {
2665 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2666 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2667 fReconstructor[iDet]->SetRecoParam(par);
2669 return fReconstructor[iDet];
2672 // load the reconstructor object
2673 TPluginManager* pluginManager = gROOT->GetPluginManager();
2674 TString detName = fgkDetectorName[iDet];
2675 TString recName = "Ali" + detName + "Reconstructor";
2677 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2679 AliReconstructor* reconstructor = NULL;
2680 // first check if a plugin is defined for the reconstructor
2681 TPluginHandler* pluginHandler =
2682 pluginManager->FindHandler("AliReconstructor", detName);
2683 // if not, add a plugin for it
2684 if (!pluginHandler) {
2685 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2686 TString libs = gSystem->GetLibraries();
2687 if (libs.Contains("lib" + detName + "base.so") ||
2688 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2689 pluginManager->AddHandler("AliReconstructor", detName,
2690 recName, detName + "rec", recName + "()");
2692 pluginManager->AddHandler("AliReconstructor", detName,
2693 recName, detName, recName + "()");
2695 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2697 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2698 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2700 if (reconstructor) {
2701 TObject* obj = fOptions.FindObject(detName.Data());
2702 if (obj) reconstructor->SetOption(obj->GetTitle());
2703 reconstructor->Init();
2704 fReconstructor[iDet] = reconstructor;
2707 // get or create the loader
2708 if (detName != "HLT") {
2709 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2710 if (!fLoader[iDet]) {
2711 AliConfig::Instance()
2712 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2714 // first check if a plugin is defined for the loader
2716 pluginManager->FindHandler("AliLoader", detName);
2717 // if not, add a plugin for it
2718 if (!pluginHandler) {
2719 TString loaderName = "Ali" + detName + "Loader";
2720 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2721 pluginManager->AddHandler("AliLoader", detName,
2722 loaderName, detName + "base",
2723 loaderName + "(const char*, TFolder*)");
2724 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2726 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2728 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2729 fRunLoader->GetEventFolder());
2731 if (!fLoader[iDet]) { // use default loader
2732 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2734 if (!fLoader[iDet]) {
2735 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2736 if (fStopOnError) return NULL;
2738 fRunLoader->AddLoader(fLoader[iDet]);
2739 fRunLoader->CdGAFile();
2740 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2741 fRunLoader->Write(0, TObject::kOverwrite);
2746 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2747 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2748 reconstructor->SetRecoParam(par);
2750 return reconstructor;
2753 //_____________________________________________________________________________
2754 AliVertexer* AliReconstruction::CreateVertexer()
2756 // create the vertexer
2757 // Please note that the caller is the owner of the
2760 AliVertexer* vertexer = NULL;
2761 AliReconstructor* itsReconstructor = GetReconstructor(0);
2762 if (itsReconstructor && (fRunLocalReconstruction.Contains("ITS"))) {
2763 vertexer = itsReconstructor->CreateVertexer();
2766 AliWarning("couldn't create a vertexer for ITS");
2772 //_____________________________________________________________________________
2773 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2775 // create the trackers
2776 AliInfo("Creating trackers");
2778 TString detStr = detectors;
2779 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2780 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2781 AliReconstructor* reconstructor = GetReconstructor(iDet);
2782 if (!reconstructor) continue;
2783 TString detName = fgkDetectorName[iDet];
2784 if (detName == "HLT") {
2785 fRunHLTTracking = kTRUE;
2788 if (detName == "MUON") {
2789 fRunMuonTracking = kTRUE;
2794 fTracker[iDet] = reconstructor->CreateTracker();
2795 if (!fTracker[iDet] && (iDet < 7)) {
2796 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2797 if (fStopOnError) return kFALSE;
2799 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2805 //_____________________________________________________________________________
2806 void AliReconstruction::CleanUp()
2808 // delete trackers and the run loader and close and delete the file
2810 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2811 delete fReconstructor[iDet];
2812 fReconstructor[iDet] = NULL;
2813 fLoader[iDet] = NULL;
2814 delete fTracker[iDet];
2815 fTracker[iDet] = NULL;
2820 delete fSPDTrackleter;
2821 fSPDTrackleter = NULL;
2830 delete fParentRawReader;
2831 fParentRawReader=NULL;
2839 AliQAManager::Destroy() ;
2841 TGeoGlobalMagField::Instance()->SetField(NULL);
2844 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2846 // Write space-points which are then used in the alignment procedures
2847 // For the moment only ITS, TPC, TRD and TOF
2849 Int_t ntracks = esd->GetNumberOfTracks();
2850 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2852 AliESDtrack *track = esd->GetTrack(itrack);
2855 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
2856 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2857 nsp += track->GetNcls(iDet);
2859 if (iDet==0) { // ITS "extra" clusters
2860 track->GetClusters(iDet,idx);
2861 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2866 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2867 track->SetTrackPointArray(sp);
2869 for (Int_t iDet = 5; iDet >= 0; iDet--) {
2870 AliTracker *tracker = fTracker[iDet];
2871 if (!tracker) continue;
2872 Int_t nspdet = track->GetClusters(iDet,idx);
2874 if (iDet==0) // ITS "extra" clusters
2875 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2877 if (nspdet <= 0) continue;
2881 while (isp2 < nspdet) {
2882 Bool_t isvalid=kTRUE;
2884 Int_t index=idx[isp++];
2885 if (index < 0) continue;
2887 TString dets = fgkDetectorName[iDet];
2888 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2889 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2890 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2891 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2892 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2894 isvalid = tracker->GetTrackPoint(index,p);
2897 if (!isvalid) continue;
2898 if (iDet==0 && (isp-1)>=6) p.SetExtra();
2899 sp->AddPoint(isptrack,&p); isptrack++;
2906 //_____________________________________________________________________________
2907 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2909 // The method reads the raw-data error log
2910 // accumulated within the rawReader.
2911 // It extracts the raw-data errors related to
2912 // the current event and stores them into
2913 // a TClonesArray inside the esd object.
2915 if (!fRawReader) return;
2917 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2919 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2921 if (iEvent != log->GetEventNumber()) continue;
2923 esd->AddRawDataErrorLog(log);
2928 //_____________________________________________________________________________
2929 void AliReconstruction::CheckQA()
2931 // check the QA of SIM for this run and remove the detectors
2932 // with status Fatal
2934 // TString newRunLocalReconstruction ;
2935 // TString newRunTracking ;
2936 // TString newFillESD ;
2938 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
2939 // TString detName(AliQAv1::GetDetName(iDet)) ;
2940 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
2941 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
2942 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
2943 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
2945 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
2946 // fRunLocalReconstruction.Contains("ALL") ) {
2947 // newRunLocalReconstruction += detName ;
2948 // newRunLocalReconstruction += " " ;
2950 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
2951 // fRunTracking.Contains("ALL") ) {
2952 // newRunTracking += detName ;
2953 // newRunTracking += " " ;
2955 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
2956 // fFillESD.Contains("ALL") ) {
2957 // newFillESD += detName ;
2958 // newFillESD += " " ;
2962 // fRunLocalReconstruction = newRunLocalReconstruction ;
2963 // fRunTracking = newRunTracking ;
2964 // fFillESD = newFillESD ;
2967 //_____________________________________________________________________________
2968 Int_t AliReconstruction::GetDetIndex(const char* detector)
2970 // return the detector index corresponding to detector
2972 for (index = 0; index < kNDetectors ; index++) {
2973 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2978 //_____________________________________________________________________________
2979 Bool_t AliReconstruction::FinishPlaneEff() {
2981 // Here execute all the necessary operationis, at the end of the tracking phase,
2982 // in case that evaluation of PlaneEfficiencies was required for some detector.
2983 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2985 // This Preliminary version works only FOR ITS !!!!!
2986 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2989 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2992 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2993 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2994 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2995 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
2996 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2997 TString name=planeeff->GetName();
2999 TFile* pefile = TFile::Open(name, "RECREATE");
3000 ret=(Bool_t)planeeff->Write();
3002 if(planeeff->GetCreateHistos()) {
3003 TString hname=planeeff->GetName();
3004 hname+="Histo.root";
3005 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3008 if(fSPDTrackleter) {
3009 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3010 TString name="AliITSPlaneEffSPDtracklet.root";
3011 TFile* pefile = TFile::Open(name, "RECREATE");
3012 ret=(Bool_t)planeeff->Write();
3014 AliESDEvent *dummy=NULL;
3015 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3020 //_____________________________________________________________________________
3021 Bool_t AliReconstruction::InitPlaneEff() {
3023 // Here execute all the necessary operations, before of the tracking phase,
3024 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3025 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3026 // which should be updated/recalculated.
3028 // This Preliminary version will work only FOR ITS !!!!!
3029 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3032 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3034 AliWarning(Form("Implementation of this method not yet completed !! Method return kTRUE"));
3036 fSPDTrackleter = NULL;
3037 AliReconstructor* itsReconstructor = GetReconstructor(0);
3038 if (itsReconstructor) {
3039 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3041 if (fSPDTrackleter) {
3042 AliInfo("Trackleter for SPD has been created");
3048 //_____________________________________________________________________________
3049 Bool_t AliReconstruction::InitAliEVE()
3051 // This method should be called only in case
3052 // AliReconstruction is run
3053 // within the alieve environment.
3054 // It will initialize AliEVE in a way
3055 // so that it can visualize event processed
3056 // by AliReconstruction.
3057 // The return flag shows whenever the
3058 // AliEVE initialization was successful or not.
3061 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3062 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3063 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3065 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3066 gROOT->ProcessLine("alieve_online_init()");
3071 //_____________________________________________________________________________
3072 void AliReconstruction::RunAliEVE()
3074 // Runs AliEVE visualisation of
3075 // the current event.
3076 // Should be executed only after
3077 // successful initialization of AliEVE.
3079 AliInfo("Running AliEVE...");
3080 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3084 //_____________________________________________________________________________
3085 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3087 // Allows to run QA for a selected set of detectors
3088 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3089 // all selected detectors run the same selected tasks
3091 if (!detAndAction.Contains(":")) {
3092 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3096 Int_t colon = detAndAction.Index(":") ;
3097 fQADetectors = detAndAction(0, colon) ;
3098 if (fQADetectors.Contains("ALL") )
3099 fQADetectors = fFillESD ;
3100 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3101 if (fQATasks.Contains("ALL") ) {
3102 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3104 fQATasks.ToUpper() ;
3106 if ( fQATasks.Contains("RAW") )
3107 tempo = Form("%d ", AliQAv1::kRAWS) ;
3108 if ( fQATasks.Contains("DIGIT") )
3109 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3110 if ( fQATasks.Contains("RECPOINT") )
3111 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3112 if ( fQATasks.Contains("ESD") )
3113 tempo += Form("%d ", AliQAv1::kESDS) ;
3115 if (fQATasks.IsNull()) {
3116 AliInfo("No QA requested\n") ;
3121 TString tempo(fQATasks) ;
3122 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3123 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3124 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3125 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3126 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3131 //_____________________________________________________________________________
3132 Bool_t AliReconstruction::InitRecoParams()
3134 // The method accesses OCDB and retrieves all
3135 // the available reco-param objects from there.
3137 Bool_t isOK = kTRUE;
3139 TString detStr = fLoadCDB;
3140 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3142 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3144 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3145 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3149 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3151 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3152 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3154 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3158 TObject *recoParamObj = entry->GetObject();
3159 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3160 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3161 // Registering them in AliRecoParam
3162 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3164 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3165 // The detector has only onse set of reco parameters
3166 // Registering it in AliRecoParam
3167 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3168 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3169 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3172 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3176 AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3180 if (AliDebugLevel() > 0) fRecoParam.Print();
3185 //_____________________________________________________________________________
3186 Bool_t AliReconstruction::GetEventInfo()
3188 // Fill the event info object
3190 AliCodeTimerAuto("")
3192 AliCentralTrigger *aCTP = NULL;
3194 fEventInfo.SetEventType(fRawReader->GetType());
3196 ULong64_t mask = fRawReader->GetClassMask();
3197 fEventInfo.SetTriggerMask(mask);
3198 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3199 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3201 aCTP = new AliCentralTrigger();
3202 TString configstr("");
3203 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3204 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3208 aCTP->SetClassMask(mask);
3209 aCTP->SetClusterMask(clmask);
3212 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3214 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3215 aCTP = fRunLoader->GetTrigger();
3216 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3217 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3220 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3225 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3227 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3228 if (fRawReader) delete aCTP;
3232 UChar_t clustmask = 0;
3234 ULong64_t trmask = fEventInfo.GetTriggerMask();
3235 const TObjArray& classesArray = config->GetClasses();
3236 Int_t nclasses = classesArray.GetEntriesFast();
3237 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3238 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3240 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3241 fesd->SetTriggerClass(trclass->GetName(),trindex);
3242 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3243 if (trmask & (1 << trindex)) {
3245 trclasses += trclass->GetName();
3247 clustmask |= trclass->GetCluster()->GetClusterMask();
3251 fEventInfo.SetTriggerClasses(trclasses);
3253 // Set the information in ESD
3254 fesd->SetTriggerMask(trmask);
3255 fesd->SetTriggerCluster(clustmask);
3257 if (!aCTP->CheckTriggeredDetectors()) {
3258 if (fRawReader) delete aCTP;
3262 if (fRawReader) delete aCTP;
3264 // We have to fill also the HLT decision here!!
3270 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3272 // Match the detector list found in the rec.C or the default 'ALL'
3273 // to the list found in the GRP (stored there by the shuttle PP which
3274 // gets the information from ECS)
3275 static TString resultList;
3276 TString detList = detectorList;
3280 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3281 if ((detectorMask >> iDet) & 0x1) {
3282 TString det = AliDAQ::OfflineModuleName(iDet);
3283 if ((detList.CompareTo("ALL") == 0) ||
3284 ((detList.BeginsWith("ALL ") ||
3285 detList.EndsWith(" ALL") ||
3286 detList.Contains(" ALL ")) &&
3287 !(detList.BeginsWith("-"+det+" ") ||
3288 detList.EndsWith(" -"+det) ||
3289 detList.Contains(" -"+det+" "))) ||
3290 (detList.CompareTo(det) == 0) ||
3291 detList.BeginsWith(det+" ") ||
3292 detList.EndsWith(" "+det) ||
3293 detList.Contains( " "+det+" " )) {
3294 if (!resultList.EndsWith(det + " ")) {
3303 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3304 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3305 if ((detList.CompareTo("ALL") == 0) ||
3306 ((detList.BeginsWith("ALL ") ||
3307 detList.EndsWith(" ALL") ||
3308 detList.Contains(" ALL ")) &&
3309 !(detList.BeginsWith("-"+hltDet+" ") ||
3310 detList.EndsWith(" -"+hltDet) ||
3311 detList.Contains(" -"+hltDet+" "))) ||
3312 (detList.CompareTo(hltDet) == 0) ||
3313 detList.BeginsWith(hltDet+" ") ||
3314 detList.EndsWith(" "+hltDet) ||
3315 detList.Contains( " "+hltDet+" " )) {
3316 resultList += hltDet;
3320 return resultList.Data();
3324 //______________________________________________________________________________
3325 void AliReconstruction::Abort(const char *method, EAbort what)
3327 // Abort processing. If what = kAbortProcess, the Process() loop will be
3328 // aborted. If what = kAbortFile, the current file in a chain will be
3329 // aborted and the processing will continue with the next file, if there
3330 // is no next file then Process() will be aborted. Abort() can also be
3331 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3332 // the SlaveTerminate() and Terminate() are always called. The abort flag
3333 // can be checked in these methods using GetAbort().
3335 // The method is overwritten in AliReconstruction for better handling of
3336 // reco specific errors
3338 if (!fStopOnError) return;
3342 TString whyMess = method;
3343 whyMess += " failed! Aborting...";
3345 AliError(whyMess.Data());
3348 TString mess = "Abort";
3349 if (fAbort == kAbortProcess)
3350 mess = "AbortProcess";
3351 else if (fAbort == kAbortFile)
3354 Info(mess, whyMess.Data());