1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for running the reconstruction //
22 // Clusters and tracks are created for all detectors and all events by //
25 // AliReconstruction rec; //
28 // The Run method returns kTRUE in case of successful execution. //
30 // If the input to the reconstruction are not simulated digits but raw data, //
31 // this can be specified by an argument of the Run method or by the method //
33 // rec.SetInput("..."); //
35 // The input formats and the corresponding argument are: //
36 // - DDL raw data files: directory name, ends with "/" //
37 // - raw data root file: root file name, extension ".root" //
38 // - raw data DATE file: DATE file name, any other non-empty string //
39 // - MC root files : empty string, default //
41 // By default all events are reconstructed. The reconstruction can be //
42 // limited to a range of events by giving the index of the first and the //
43 // last event as an argument to the Run method or by calling //
45 // rec.SetEventRange(..., ...); //
47 // The index -1 (default) can be used for the last event to indicate no //
48 // upper limit of the event range. //
50 // In case of raw-data reconstruction the user can modify the default //
51 // number of events per digits/clusters/tracks file. In case the option //
52 // is not used the number is set 1. In case the user provides 0, than //
53 // the number of events is equal to the number of events inside the //
54 // raw-data file (i.e. one digits/clusters/tracks file): //
56 // rec.SetNumberOfEventsPerFile(...); //
59 // The name of the galice file can be changed from the default //
60 // "galice.root" by passing it as argument to the AliReconstruction //
61 // constructor or by //
63 // rec.SetGAliceFile("..."); //
65 // The local reconstruction can be switched on or off for individual //
68 // rec.SetRunLocalReconstruction("..."); //
70 // The argument is a (case sensitive) string with the names of the //
71 // detectors separated by a space. The special string "ALL" selects all //
72 // available detectors. This is the default. //
74 // The reconstruction of the primary vertex position can be switched off by //
76 // rec.SetRunVertexFinder(kFALSE); //
78 // The tracking and the creation of ESD tracks can be switched on for //
79 // selected detectors by //
81 // rec.SetRunTracking("..."); //
83 // Uniform/nonuniform field tracking switches (default: uniform field) //
85 // rec.SetUniformFieldTracking(); ( rec.SetUniformFieldTracking(kFALSE); ) //
87 // The filling of additional ESD information can be steered by //
89 // rec.SetFillESD("..."); //
91 // Again, for both methods the string specifies the list of detectors. //
92 // The default is "ALL". //
94 // The call of the shortcut method //
96 // rec.SetRunReconstruction("..."); //
98 // is equivalent to calling SetRunLocalReconstruction, SetRunTracking and //
99 // SetFillESD with the same detector selecting string as argument. //
101 // The reconstruction requires digits or raw data as input. For the creation //
102 // of digits and raw data have a look at the class AliSimulation. //
104 // The input data of a detector can be replaced by the corresponding HLT //
105 // data by calling (usual detector string) //
106 // SetUseHLTData("..."); //
109 ///////////////////////////////////////////////////////////////////////////////
116 #include <TGeoGlobalMagField.h>
117 #include <TGeoManager.h>
119 #include <TLorentzVector.h>
121 #include <TObjArray.h>
122 #include <TPRegexp.h>
123 #include <TParameter.h>
124 #include <TPluginManager.h>
126 #include <TProofOutputFile.h>
129 #include <THashTable.h>
131 #include <TMessage.h>
133 #include "AliAlignObj.h"
134 #include "AliCDBEntry.h"
135 #include "AliCDBManager.h"
136 #include "AliCDBStorage.h"
137 #include "AliCTPRawStream.h"
138 #include "AliCascadeVertexer.h"
139 #include "AliCentralTrigger.h"
140 #include "AliCodeTimer.h"
142 #include "AliDetectorRecoParam.h"
143 #include "AliESDCaloCells.h"
144 #include "AliESDCaloCluster.h"
145 #include "AliESDEvent.h"
146 #include "AliESDMuonTrack.h"
147 #include "AliESDPmdTrack.h"
148 #include "AliESDTagCreator.h"
149 #include "AliESDVertex.h"
150 #include "AliESDcascade.h"
151 #include "AliESDfriend.h"
152 #include "AliESDkink.h"
153 #include "AliESDpid.h"
154 #include "AliESDtrack.h"
155 #include "AliESDtrack.h"
156 #include "AliEventInfo.h"
157 #include "AliGRPObject.h"
158 #include "AliGRPRecoParam.h"
159 #include "AliGenEventHeader.h"
160 #include "AliGeomManager.h"
161 #include "AliGlobalQADataMaker.h"
162 #include "AliHeader.h"
165 #include "AliMultiplicity.h"
167 #include "AliPlaneEff.h"
169 #include "AliQADataMakerRec.h"
170 #include "AliQAManager.h"
171 #include "AliRawVEvent.h"
172 #include "AliRawEventHeaderBase.h"
173 #include "AliRawHLTManager.h"
174 #include "AliRawReaderDate.h"
175 #include "AliRawReaderFile.h"
176 #include "AliRawReaderRoot.h"
177 #include "AliReconstruction.h"
178 #include "AliReconstructor.h"
180 #include "AliRunInfo.h"
181 #include "AliRunLoader.h"
182 #include "AliSysInfo.h" // memory snapshots
183 #include "AliTrackPointArray.h"
184 #include "AliTracker.h"
185 #include "AliTriggerClass.h"
186 #include "AliTriggerCluster.h"
187 #include "AliTriggerIR.h"
188 #include "AliTriggerConfiguration.h"
189 #include "AliV0vertexer.h"
190 #include "AliVertexer.h"
191 #include "AliVertexerTracks.h"
193 ClassImp(AliReconstruction)
195 //_____________________________________________________________________________
196 const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
198 //_____________________________________________________________________________
199 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
201 fRunVertexFinder(kTRUE),
202 fRunVertexFinderTracks(kTRUE),
203 fRunHLTTracking(kFALSE),
204 fRunMuonTracking(kFALSE),
206 fRunCascadeFinder(kTRUE),
207 fStopOnError(kFALSE),
208 fWriteAlignmentData(kFALSE),
209 fWriteESDfriend(kFALSE),
210 fFillTriggerESD(kTRUE),
218 fRunLocalReconstruction("ALL"),
222 fUseTrackingErrorsForAlignment(""),
223 fGAliceFileName(gAliceFilename),
226 fProofOutputFileName(""),
227 fProofOutputLocation(""),
228 fProofOutputDataset(kFALSE),
229 fProofOutputArchive(""),
233 fNumberOfEventsPerFile((UInt_t)-1),
235 fLoadAlignFromCDB(kTRUE),
236 fLoadAlignData("ALL"),
244 fParentRawReader(NULL),
248 fSPDTrackleter(NULL),
250 fDiamondProfileSPD(NULL),
251 fDiamondProfile(NULL),
252 fDiamondProfileTPC(NULL),
253 fListOfCosmicTriggers(NULL),
257 fAlignObjArray(NULL),
261 fInitCDBCalled(kFALSE),
262 fSetRunNumberFromDataCalled(kFALSE),
267 fSameQACycle(kFALSE),
268 fInitQACalled(kFALSE),
269 fWriteQAExpertData(kTRUE),
270 fRunPlaneEff(kFALSE),
279 fIsNewRunLoader(kFALSE),
283 // create reconstruction object with default parameters
286 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
287 fReconstructor[iDet] = NULL;
288 fLoader[iDet] = NULL;
289 fTracker[iDet] = NULL;
291 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
292 fQACycles[iDet] = 999999 ;
293 fQAWriteExpert[iDet] = kFALSE ;
299 //_____________________________________________________________________________
300 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
302 fRunVertexFinder(rec.fRunVertexFinder),
303 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
304 fRunHLTTracking(rec.fRunHLTTracking),
305 fRunMuonTracking(rec.fRunMuonTracking),
306 fRunV0Finder(rec.fRunV0Finder),
307 fRunCascadeFinder(rec.fRunCascadeFinder),
308 fStopOnError(rec.fStopOnError),
309 fWriteAlignmentData(rec.fWriteAlignmentData),
310 fWriteESDfriend(rec.fWriteESDfriend),
311 fFillTriggerESD(rec.fFillTriggerESD),
313 fCleanESD(rec.fCleanESD),
314 fV0DCAmax(rec.fV0DCAmax),
315 fV0CsPmin(rec.fV0CsPmin),
319 fRunLocalReconstruction(rec.fRunLocalReconstruction),
320 fRunTracking(rec.fRunTracking),
321 fFillESD(rec.fFillESD),
322 fLoadCDB(rec.fLoadCDB),
323 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
324 fGAliceFileName(rec.fGAliceFileName),
325 fRawInput(rec.fRawInput),
326 fESDOutput(rec.fESDOutput),
327 fProofOutputFileName(rec.fProofOutputFileName),
328 fProofOutputLocation(rec.fProofOutputLocation),
329 fProofOutputDataset(rec.fProofOutputDataset),
330 fProofOutputArchive(rec.fProofOutputArchive),
331 fEquipIdMap(rec.fEquipIdMap),
332 fFirstEvent(rec.fFirstEvent),
333 fLastEvent(rec.fLastEvent),
334 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
336 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
337 fLoadAlignData(rec.fLoadAlignData),
338 fUseHLTData(rec.fUseHLTData),
345 fParentRawReader(NULL),
347 fRecoParam(rec.fRecoParam),
349 fSPDTrackleter(NULL),
351 fDiamondProfileSPD(rec.fDiamondProfileSPD),
352 fDiamondProfile(rec.fDiamondProfile),
353 fDiamondProfileTPC(rec.fDiamondProfileTPC),
354 fListOfCosmicTriggers(NULL),
358 fAlignObjArray(rec.fAlignObjArray),
359 fCDBUri(rec.fCDBUri),
360 fQARefUri(rec.fQARefUri),
362 fInitCDBCalled(rec.fInitCDBCalled),
363 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
364 fQADetectors(rec.fQADetectors),
365 fQATasks(rec.fQATasks),
367 fRunGlobalQA(rec.fRunGlobalQA),
368 fSameQACycle(rec.fSameQACycle),
369 fInitQACalled(rec.fInitQACalled),
370 fWriteQAExpertData(rec.fWriteQAExpertData),
371 fRunPlaneEff(rec.fRunPlaneEff),
380 fIsNewRunLoader(rec.fIsNewRunLoader),
386 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
387 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
389 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
390 fReconstructor[iDet] = NULL;
391 fLoader[iDet] = NULL;
392 fTracker[iDet] = NULL;
395 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
396 fQACycles[iDet] = rec.fQACycles[iDet];
397 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
400 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
401 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
405 //_____________________________________________________________________________
406 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
408 // assignment operator
409 // Used in PROOF mode
410 // Be very careful while modifing it!
411 // Simple rules to follow:
412 // for persistent data members - use their assignment operators
413 // for non-persistent ones - do nothing or take the default values from constructor
414 // TSelector members should not be touched
415 if(&rec == this) return *this;
417 fRunVertexFinder = rec.fRunVertexFinder;
418 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
419 fRunHLTTracking = rec.fRunHLTTracking;
420 fRunMuonTracking = rec.fRunMuonTracking;
421 fRunV0Finder = rec.fRunV0Finder;
422 fRunCascadeFinder = rec.fRunCascadeFinder;
423 fStopOnError = rec.fStopOnError;
424 fWriteAlignmentData = rec.fWriteAlignmentData;
425 fWriteESDfriend = rec.fWriteESDfriend;
426 fFillTriggerESD = rec.fFillTriggerESD;
428 fCleanESD = rec.fCleanESD;
429 fV0DCAmax = rec.fV0DCAmax;
430 fV0CsPmin = rec.fV0CsPmin;
434 fRunLocalReconstruction = rec.fRunLocalReconstruction;
435 fRunTracking = rec.fRunTracking;
436 fFillESD = rec.fFillESD;
437 fLoadCDB = rec.fLoadCDB;
438 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
439 fGAliceFileName = rec.fGAliceFileName;
440 fRawInput = rec.fRawInput;
441 fESDOutput = rec.fESDOutput;
442 fProofOutputFileName = rec.fProofOutputFileName;
443 fProofOutputLocation = rec.fProofOutputLocation;
444 fProofOutputDataset = rec.fProofOutputDataset;
445 fProofOutputArchive = rec.fProofOutputArchive;
446 fEquipIdMap = rec.fEquipIdMap;
447 fFirstEvent = rec.fFirstEvent;
448 fLastEvent = rec.fLastEvent;
449 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
451 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
452 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
455 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
456 fLoadAlignData = rec.fLoadAlignData;
457 fUseHLTData = rec.fUseHLTData;
459 delete fRunInfo; fRunInfo = NULL;
460 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
462 fEventInfo = rec.fEventInfo;
464 delete fRunScalers; fRunScalers = NULL;
465 if (rec.fRunScalers) fRunScalers = new AliTriggerRunScalers(*rec.fRunScalers);
470 fParentRawReader = NULL;
472 fRecoParam = rec.fRecoParam;
474 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
475 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
476 delete fLoader[iDet]; fLoader[iDet] = NULL;
477 delete fTracker[iDet]; fTracker[iDet] = NULL;
480 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
481 fQACycles[iDet] = rec.fQACycles[iDet];
482 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
485 delete fSPDTrackleter; fSPDTrackleter = NULL;
487 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
488 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
489 delete fDiamondProfile; fDiamondProfile = NULL;
490 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
491 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
492 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
494 delete fListOfCosmicTriggers; fListOfCosmicTriggers = NULL;
495 if (rec.fListOfCosmicTriggers) fListOfCosmicTriggers = (THashTable*)((rec.fListOfCosmicTriggers)->Clone());
497 delete fGRPData; fGRPData = NULL;
498 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
499 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
501 delete fAlignObjArray; fAlignObjArray = NULL;
504 fQARefUri = rec.fQARefUri;
505 fSpecCDBUri.Delete();
506 fInitCDBCalled = rec.fInitCDBCalled;
507 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
508 fQADetectors = rec.fQADetectors;
509 fQATasks = rec.fQATasks;
511 fRunGlobalQA = rec.fRunGlobalQA;
512 fSameQACycle = rec.fSameQACycle;
513 fInitQACalled = rec.fInitQACalled;
514 fWriteQAExpertData = rec.fWriteQAExpertData;
515 fRunPlaneEff = rec.fRunPlaneEff;
524 fIsNewRunLoader = rec.fIsNewRunLoader;
531 //_____________________________________________________________________________
532 AliReconstruction::~AliReconstruction()
537 if (fListOfCosmicTriggers) {
538 fListOfCosmicTriggers->Delete();
539 delete fListOfCosmicTriggers;
544 if (fAlignObjArray) {
545 fAlignObjArray->Delete();
546 delete fAlignObjArray;
548 fSpecCDBUri.Delete();
550 AliCodeTimer::Instance()->Print();
553 //_____________________________________________________________________________
554 void AliReconstruction::InitQA()
556 //Initialize the QA and start of cycle
557 AliCodeTimerAuto("");
559 if (fInitQACalled) return;
560 fInitQACalled = kTRUE;
562 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kRECMODE) ;
563 if (fWriteQAExpertData)
564 qam->SetWriteExpert() ;
566 if (qam->IsDefaultStorageSet()) {
567 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
568 AliWarning("Default QA reference storage has been already set !");
569 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fQARefUri.Data()));
570 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
571 fQARefUri = qam->GetDefaultStorage()->GetURI();
573 if (fQARefUri.Length() > 0) {
574 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
575 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
576 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
578 fQARefUri="local://$ALICE_ROOT/QAref";
579 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
580 AliWarning("Default QA refeference storage not yet set !!!!");
581 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
582 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
585 qam->SetDefaultStorage(fQARefUri);
589 qam->SetActiveDetectors(fQADetectors) ;
590 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) {
591 qam->SetCycleLength(AliQAv1::DETECTORINDEX_t(det), fQACycles[det]) ;
592 qam->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
594 if (!fRawReader && !fInput && fQATasks.Contains(AliQAv1::kRAWS))
595 fQATasks.ReplaceAll(Form("%d",AliQAv1::kRAWS), "") ;
596 qam->SetTasks(fQATasks) ;
597 qam->InitQADataMaker(AliCDBManager::Instance()->GetRun()) ;
600 Bool_t sameCycle = kFALSE ;
601 AliQADataMaker *qadm = qam->GetQADataMaker(AliQAv1::kGLOBAL);
602 AliInfo(Form("Initializing the global QA data maker"));
603 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS))) {
604 qadm->StartOfCycle(AliQAv1::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
605 TObjArray **arr=qadm->Init(AliQAv1::kRECPOINTS);
606 AliTracker::SetResidualsArray(arr);
609 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS))) {
610 qadm->StartOfCycle(AliQAv1::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
611 qadm->Init(AliQAv1::kESDS);
614 AliSysInfo::AddStamp("InitQA") ;
617 //_____________________________________________________________________________
618 void AliReconstruction::MergeQA(const char *fileName)
620 //Initialize the QA and start of cycle
621 AliCodeTimerAuto("") ;
622 AliQAManager::QAManager()->Merge(AliCDBManager::Instance()->GetRun(),fileName) ;
623 AliSysInfo::AddStamp("MergeQA") ;
626 //_____________________________________________________________________________
627 void AliReconstruction::InitCDB()
629 // activate a default CDB storage
630 // First check if we have any CDB storage set, because it is used
631 // to retrieve the calibration and alignment constants
632 AliCodeTimerAuto("");
634 if (fInitCDBCalled) return;
635 fInitCDBCalled = kTRUE;
637 AliCDBManager* man = AliCDBManager::Instance();
638 if (man->IsDefaultStorageSet())
640 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
641 AliWarning("Default CDB storage has been already set !");
642 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
643 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
644 fCDBUri = man->GetDefaultStorage()->GetURI();
647 if (fCDBUri.Length() > 0)
649 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
650 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
651 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
652 man->SetDefaultStorage(fCDBUri);
654 else if (!man->GetRaw()){
655 fCDBUri="local://$ALICE_ROOT/OCDB";
656 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
657 AliWarning("Default CDB storage not yet set !!!!");
658 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
659 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
660 man->SetDefaultStorage(fCDBUri);
663 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
664 AliWarning("Default storage will be set after setting the Run Number!!!");
665 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
669 // Now activate the detector specific CDB storage locations
670 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
671 TObject* obj = fSpecCDBUri[i];
673 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
674 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
675 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
676 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
678 AliSysInfo::AddStamp("InitCDB");
681 //_____________________________________________________________________________
682 void AliReconstruction::SetDefaultStorage(const char* uri) {
683 // Store the desired default CDB storage location
684 // Activate it later within the Run() method
690 //_____________________________________________________________________________
691 void AliReconstruction::SetQARefDefaultStorage(const char* uri) {
692 // Store the desired default CDB storage location
693 // Activate it later within the Run() method
696 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
699 //_____________________________________________________________________________
700 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
701 // Store a detector-specific CDB storage location
702 // Activate it later within the Run() method
704 AliCDBPath aPath(calibType);
705 if(!aPath.IsValid()){
706 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
707 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
708 if(!strcmp(calibType, fgkDetectorName[iDet])) {
709 aPath.SetPath(Form("%s/*", calibType));
710 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
714 if(!aPath.IsValid()){
715 AliError(Form("Not a valid path or detector: %s", calibType));
720 // // check that calibType refers to a "valid" detector name
721 // Bool_t isDetector = kFALSE;
722 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
723 // TString detName = fgkDetectorName[iDet];
724 // if(aPath.GetLevel0() == detName) {
725 // isDetector = kTRUE;
731 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
735 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
736 if (obj) fSpecCDBUri.Remove(obj);
737 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
741 //_____________________________________________________________________________
742 Bool_t AliReconstruction::SetRunNumberFromData()
744 // The method is called in Run() in order
745 // to set a correct run number.
746 // In case of raw data reconstruction the
747 // run number is taken from the raw data header
749 if (fSetRunNumberFromDataCalled) return kTRUE;
750 fSetRunNumberFromDataCalled = kTRUE;
752 AliCDBManager* man = AliCDBManager::Instance();
755 if(fRawReader->NextEvent()) {
756 if(man->GetRun() > 0) {
757 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
759 man->SetRun(fRawReader->GetRunNumber());
760 fRawReader->RewindEvents();
763 if(man->GetRun() > 0) {
764 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
767 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
773 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
775 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
780 // read run number from gAlice
781 if(rl->GetHeader()) {
782 man->SetRun(rl->GetHeader()->GetRun());
787 AliError("Neither run-loader header nor RawReader objects are found !");
799 //_____________________________________________________________________________
800 void AliReconstruction::SetCDBLock() {
801 // Set CDB lock: from now on it is forbidden to reset the run number
802 // or the default storage or to activate any further storage!
804 AliCDBManager::Instance()->SetLock(1);
807 //_____________________________________________________________________________
808 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
810 // Read the alignment objects from CDB.
811 // Each detector is supposed to have the
812 // alignment objects in DET/Align/Data CDB path.
813 // All the detector objects are then collected,
814 // sorted by geometry level (starting from ALIC) and
815 // then applied to the TGeo geometry.
816 // Finally an overlaps check is performed.
818 // Load alignment data from CDB and fill fAlignObjArray
819 if(fLoadAlignFromCDB){
821 TString detStr = detectors;
822 TString loadAlObjsListOfDets = "";
824 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
825 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
826 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
828 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
830 loadAlObjsListOfDets += fgkDetectorName[iDet];
831 loadAlObjsListOfDets += " ";
833 } // end loop over detectors
835 if(AliGeomManager::GetNalignable("GRP") != 0)
836 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
837 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
838 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
840 // Check if the array with alignment objects was
841 // provided by the user. If yes, apply the objects
842 // to the present TGeo geometry
843 if (fAlignObjArray) {
844 if (gGeoManager && gGeoManager->IsClosed()) {
845 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
846 AliError("The misalignment of one or more volumes failed!"
847 "Compare the list of simulated detectors and the list of detector alignment data!");
852 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
858 if (fAlignObjArray) {
859 fAlignObjArray->Delete();
860 delete fAlignObjArray; fAlignObjArray=NULL;
866 //_____________________________________________________________________________
867 void AliReconstruction::SetGAliceFile(const char* fileName)
869 // set the name of the galice file
871 fGAliceFileName = fileName;
874 //_____________________________________________________________________________
875 void AliReconstruction::SetInput(const char* input)
877 // In case the input string starts with 'mem://', we run in an online mode
878 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
879 // file is assumed. One can give as an input:
880 // mem://: - events taken from DAQ monitoring libs online
882 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
883 if (input) fRawInput = input;
886 //_____________________________________________________________________________
887 void AliReconstruction::SetOutput(const char* output)
889 // Set the output ESD filename
890 // 'output' is a normalt ROOT url
891 // The method is used in case of raw-data reco with PROOF
892 if (output) fESDOutput = output;
895 //_____________________________________________________________________________
896 void AliReconstruction::SetOption(const char* detector, const char* option)
898 // set options for the reconstruction of a detector
900 TObject* obj = fOptions.FindObject(detector);
901 if (obj) fOptions.Remove(obj);
902 fOptions.Add(new TNamed(detector, option));
905 //_____________________________________________________________________________
906 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
908 // Set custom reconstruction parameters for a given detector
909 // Single set of parameters for all the events
911 // First check if the reco-params are global
912 if(!strcmp(detector, "GRP")) {
914 fRecoParam.AddDetRecoParam(kNDetectors,par);
918 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
919 if(!strcmp(detector, fgkDetectorName[iDet])) {
921 fRecoParam.AddDetRecoParam(iDet,par);
928 //_____________________________________________________________________________
929 Bool_t AliReconstruction::InitGRP() {
930 //------------------------------------
931 // Initialization of the GRP entry
932 //------------------------------------
933 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
937 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
940 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
942 fGRPData = new AliGRPObject();
943 fGRPData->ReadValuesFromMap(m);
947 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
948 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
952 // FIX ME: The unloading of GRP entry is temporarily disabled
953 // because ZDC and VZERO are using it in order to initialize
954 // their reconstructor objects. In the future one has to think
955 // of propagating AliRunInfo to the reconstructors.
956 // AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
960 AliError("No GRP entry found in OCDB!");
964 TString lhcState = fGRPData->GetLHCState();
965 if (lhcState==AliGRPObject::GetInvalidString()) {
966 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
967 lhcState = "UNKNOWN";
970 TString beamType = fGRPData->GetBeamType();
971 if (beamType==AliGRPObject::GetInvalidString()) {
972 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
973 beamType = "UNKNOWN";
976 Float_t beamEnergy = fGRPData->GetBeamEnergy();
977 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
978 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
981 // LHC: "multiply by 120 to get the energy in MeV"
984 TString runType = fGRPData->GetRunType();
985 if (runType==AliGRPObject::GetInvalidString()) {
986 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
990 Int_t activeDetectors = fGRPData->GetDetectorMask();
991 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
992 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
993 activeDetectors = 1074790399;
996 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1000 // Process the list of active detectors
1001 if (activeDetectors) {
1002 UInt_t detMask = activeDetectors;
1003 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1004 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1005 fFillESD = MatchDetectorList(fFillESD,detMask);
1006 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1007 fLoadCDB.Form("%s %s %s %s",
1008 fRunLocalReconstruction.Data(),
1009 fRunTracking.Data(),
1011 fQADetectors.Data());
1012 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1013 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1)) {
1014 // switch off the vertexer
1015 AliInfo("SPD is not in the list of active detectors. Vertexer switched off.");
1016 fRunVertexFinder = kFALSE;
1018 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1019 // switch off the reading of CTP raw-data payload
1020 if (fFillTriggerESD) {
1021 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1022 fFillTriggerESD = kFALSE;
1027 AliInfo("===================================================================================");
1028 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1029 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1030 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1031 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1032 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1033 AliInfo("===================================================================================");
1035 //*** Dealing with the magnetic field map
1036 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {
1037 if (TGeoGlobalMagField::Instance()->GetField()->TestBit(AliMagF::kOverrideGRP)) {
1038 AliInfo("ExpertMode!!! GRP information will be ignored !");
1039 AliInfo("ExpertMode!!! Running with the externally locked B field !");
1042 AliInfo("Destroying existing B field instance!");
1043 delete TGeoGlobalMagField::Instance();
1046 if ( !TGeoGlobalMagField::Instance()->IsLocked() ) {
1047 // Construct the field map out of the information retrieved from GRP.
1050 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1051 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1052 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1056 Char_t l3Polarity = fGRPData->GetL3Polarity();
1057 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1058 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1063 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1064 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1065 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1069 Char_t diPolarity = fGRPData->GetDipolePolarity();
1070 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1071 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1075 // read special bits for the polarity convention and map type
1076 Int_t polConvention = fGRPData->IsPolarityConventionLHC() ? AliMagF::kConvLHC : AliMagF::kConvDCS2008;
1077 Bool_t uniformB = fGRPData->IsUniformBMap();
1080 AliMagF* fld = AliMagF::CreateFieldMap(TMath::Abs(l3Current) * (l3Polarity ? -1:1),
1081 TMath::Abs(diCurrent) * (diPolarity ? -1:1),
1082 polConvention,uniformB,beamEnergy, beamType.Data());
1084 TGeoGlobalMagField::Instance()->SetField( fld );
1085 TGeoGlobalMagField::Instance()->Lock();
1086 AliInfo("Running with the B field constructed out of GRP !");
1088 else AliFatal("Failed to create a B field map !");
1090 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1093 //*** Get the diamond profiles from OCDB
1094 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1096 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1098 AliError("No SPD diamond profile found in OCDB!");
1101 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1103 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1105 AliError("No diamond profile found in OCDB!");
1108 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1110 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1112 AliError("No TPC diamond profile found in OCDB!");
1115 entry = AliCDBManager::Instance()->Get("GRP/Calib/CosmicTriggers");
1117 fListOfCosmicTriggers = dynamic_cast<THashTable*>(entry->GetObject());
1119 AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
1122 if (!fListOfCosmicTriggers) {
1123 AliWarning("Can not get list of cosmic triggers from OCDB! Cosmic event specie will be effectively disabled!");
1129 //_____________________________________________________________________________
1130 Bool_t AliReconstruction::LoadCDB()
1132 AliCodeTimerAuto("");
1134 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1136 TString detStr = fLoadCDB;
1137 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1138 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1139 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1143 //_____________________________________________________________________________
1144 Bool_t AliReconstruction::LoadTriggerScalersCDB()
1146 AliCodeTimerAuto("");
1148 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
1152 AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
1153 fRunScalers = dynamic_cast<AliTriggerRunScalers*> (entry->GetObject());
1155 if (fRunScalers->CorrectScalersOverflow() == 0) AliInfo("32bit Trigger counters corrected for overflow");
1160 //_____________________________________________________________________________
1161 Bool_t AliReconstruction::Run(const char* input)
1164 AliCodeTimerAuto("");
1167 if (GetAbort() != TSelector::kContinue) return kFALSE;
1169 TChain *chain = NULL;
1170 if (fRawReader && (chain = fRawReader->GetChain())) {
1171 Long64_t nEntries = (fLastEvent < 0) ? (TChain::kBigNumber) : (fLastEvent - fFirstEvent + 1);
1176 gProof->Exec("TGrid::Connect(\"alien://\")",kTRUE);
1178 TMessage::EnableSchemaEvolutionForAll(kTRUE);
1179 gProof->Exec("TMessage::EnableSchemaEvolutionForAll(kTRUE)",kTRUE);
1181 gProof->AddInput(this);
1183 if (!ParseOutput()) return kFALSE;
1185 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1187 chain->Process("AliReconstruction","",nEntries,fFirstEvent);
1190 chain->Process(this,"",nEntries,fFirstEvent);
1195 if (GetAbort() != TSelector::kContinue) return kFALSE;
1197 if (GetAbort() != TSelector::kContinue) return kFALSE;
1198 //******* The loop over events
1199 AliInfo("Starting looping over events");
1201 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1202 (fRawReader && fRawReader->NextEvent())) {
1203 if (!ProcessEvent(iEvent)) {
1204 Abort("ProcessEvent",TSelector::kAbortFile);
1210 if (GetAbort() != TSelector::kContinue) return kFALSE;
1212 if (GetAbort() != TSelector::kContinue) return kFALSE;
1218 //_____________________________________________________________________________
1219 void AliReconstruction::InitRawReader(const char* input)
1221 AliCodeTimerAuto("");
1223 // Init raw-reader and
1224 // set the input in case of raw data
1225 if (input) fRawInput = input;
1226 fRawReader = AliRawReader::Create(fRawInput.Data());
1228 AliInfo("Reconstruction will run over digits");
1230 if (!fEquipIdMap.IsNull() && fRawReader)
1231 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1233 if (!fUseHLTData.IsNull()) {
1234 // create the RawReaderHLT which performs redirection of HLT input data for
1235 // the specified detectors
1236 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1238 fParentRawReader=fRawReader;
1239 fRawReader=pRawReader;
1241 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1244 AliSysInfo::AddStamp("CreateRawReader");
1247 //_____________________________________________________________________________
1248 void AliReconstruction::InitRun(const char* input)
1250 // Initialization of raw-reader,
1251 // run number, CDB etc.
1252 AliCodeTimerAuto("");
1253 AliSysInfo::AddStamp("Start");
1255 // Initialize raw-reader if any
1256 InitRawReader(input);
1258 // Initialize the CDB storage
1261 // Set run number in CDBManager (if it is not already set by the user)
1262 if (!SetRunNumberFromData()) {
1263 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1267 // Set CDB lock: from now on it is forbidden to reset the run number
1268 // or the default storage or to activate any further storage!
1273 //_____________________________________________________________________________
1274 void AliReconstruction::Begin(TTree *)
1276 // Initialize AlReconstruction before
1277 // going into the event loop
1278 // Should follow the TSelector convention
1279 // i.e. initialize only the object on the client side
1280 AliCodeTimerAuto("");
1282 AliReconstruction *reco = NULL;
1284 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1287 AliSysInfo::AddStamp("ReadInputInBegin");
1290 // Import ideal TGeo geometry and apply misalignment
1292 TString geom(gSystem->DirName(fGAliceFileName));
1293 geom += "/geometry.root";
1294 AliGeomManager::LoadGeometry(geom.Data());
1296 Abort("LoadGeometry", TSelector::kAbortProcess);
1299 AliSysInfo::AddStamp("LoadGeom");
1300 TString detsToCheck=fRunLocalReconstruction;
1301 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1302 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1305 AliSysInfo::AddStamp("CheckGeom");
1308 if (!MisalignGeometry(fLoadAlignData)) {
1309 Abort("MisalignGeometry", TSelector::kAbortProcess);
1312 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1313 AliSysInfo::AddStamp("MisalignGeom");
1316 Abort("InitGRP", TSelector::kAbortProcess);
1319 AliSysInfo::AddStamp("InitGRP");
1322 Abort("LoadCDB", TSelector::kAbortProcess);
1325 AliSysInfo::AddStamp("LoadCDB");
1327 if (!LoadTriggerScalersCDB()) {
1328 Abort("LoadTriggerScalersCDB", TSelector::kAbortProcess);
1331 AliSysInfo::AddStamp("LoadTriggerScalersCDB");
1334 // Read the reconstruction parameters from OCDB
1335 if (!InitRecoParams()) {
1336 AliWarning("Not all detectors have correct RecoParam objects initialized");
1338 AliSysInfo::AddStamp("InitRecoParams");
1340 if (fInput && gProof) {
1341 if (reco) *reco = *this;
1343 gGeoManager->SetName("Geometry");
1344 gProof->AddInputData(gGeoManager,kTRUE);
1346 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1347 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1348 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1349 magFieldMap->SetName("MagneticFieldMap");
1350 gProof->AddInputData(magFieldMap,kTRUE);
1355 //_____________________________________________________________________________
1356 void AliReconstruction::SlaveBegin(TTree*)
1358 // Initialization related to run-loader,
1359 // vertexer, trackers, recontructors
1360 // In proof mode it is executed on the slave
1361 AliCodeTimerAuto("");
1363 TProofOutputFile *outProofFile = NULL;
1365 if (AliDebugLevel() > 0) fInput->Print();
1366 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1369 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1371 AliGeomManager::SetGeometry(tgeo);
1373 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1374 Int_t runNumber = -1;
1375 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1376 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1377 man->SetCacheFlag(kTRUE);
1378 man->SetLock(kTRUE);
1382 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1383 TGeoGlobalMagField::Instance()->SetField(map);
1385 if (TNamed *outputFileName = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE"))
1386 fProofOutputFileName = outputFileName->GetTitle();
1387 if (TNamed *outputLocation = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_LOCATION"))
1388 fProofOutputLocation = outputLocation->GetTitle();
1389 if (fInput->FindObject("PROOF_OUTPUTFILE_DATASET"))
1390 fProofOutputDataset = kTRUE;
1391 if (TNamed *archiveList = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_ARCHIVE"))
1392 fProofOutputArchive = archiveList->GetTitle();
1393 if (!fProofOutputFileName.IsNull() &&
1394 !fProofOutputLocation.IsNull() &&
1395 fProofOutputArchive.IsNull()) {
1396 if (!fProofOutputDataset) {
1397 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"M");
1398 outProofFile->SetOutputFileName(Form("%s%s",fProofOutputLocation.Data(),fProofOutputFileName.Data()));
1401 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"DROV",fProofOutputLocation.Data());
1403 if (AliDebugLevel() > 0) outProofFile->Dump();
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 AliInfo(Form("Opening output PROOF file: %s/%s",
1436 outProofFile->GetDir(), outProofFile->GetFileName()));
1437 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1438 Abort(Form("Problems opening output PROOF file: %s/%s",
1439 outProofFile->GetDir(), outProofFile->GetFileName()),
1440 TSelector::kAbortProcess);
1445 ftree = new TTree("esdTree", "Tree with ESD objects");
1446 fesd = new AliESDEvent();
1447 fesd->CreateStdContent();
1449 fesd->WriteToTree(ftree);
1450 if (fWriteESDfriend) {
1452 // Since we add the branch manually we must
1453 // book and add it after WriteToTree
1454 // otherwise it is created twice,
1455 // once via writetotree and once here.
1456 // The case for AliESDfriend is now
1457 // caught also in AlIESDEvent::WriteToTree but
1458 // be careful when changing the name (AliESDfriend is not
1459 // a TNamed so we had to hardwire it)
1460 fesdf = new AliESDfriend();
1461 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1462 br->SetFile("AliESDfriends.root");
1463 fesd->AddObject(fesdf);
1465 ftree->GetUserInfo()->Add(fesd);
1467 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1468 fhltesd = new AliESDEvent();
1469 fhltesd->CreateStdContent();
1471 // read the ESD template from CDB
1472 // HLT is allowed to put non-std content to its ESD, the non-std
1473 // objects need to be created before invocation of WriteToTree in
1474 // order to create all branches. Initialization is done from an
1475 // ESD layout template in CDB
1476 AliCDBManager* man = AliCDBManager::Instance();
1477 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1478 AliCDBEntry* hltESDConfig=NULL;
1479 if (man->GetId(hltESDConfigPath)!=NULL &&
1480 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1481 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1483 // init all internal variables from the list of objects
1484 pESDLayout->GetStdContent();
1486 // copy content and create non-std objects
1487 *fhltesd=*pESDLayout;
1490 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1491 hltESDConfigPath.GetPath().Data()));
1495 fhltesd->WriteToTree(fhlttree);
1496 fhlttree->GetUserInfo()->Add(fhltesd);
1498 ProcInfo_t procInfo;
1499 gSystem->GetProcInfo(&procInfo);
1500 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1503 //Initialize the QA and start of cycle
1504 if (fRunQA || fRunGlobalQA)
1507 //Initialize the Plane Efficiency framework
1508 if (fRunPlaneEff && !InitPlaneEff()) {
1509 Abort("InitPlaneEff", TSelector::kAbortProcess);
1513 if (strcmp(gProgName,"alieve") == 0)
1514 fRunAliEVE = InitAliEVE();
1519 //_____________________________________________________________________________
1520 Bool_t AliReconstruction::Process(Long64_t entry)
1522 // run the reconstruction over a single entry
1523 // from the chain with raw data
1524 AliCodeTimerAuto("");
1526 TTree *currTree = fChain->GetTree();
1527 AliRawVEvent *event = NULL;
1528 currTree->SetBranchAddress("rawevent",&event);
1529 currTree->GetEntry(entry);
1530 fRawReader = new AliRawReaderRoot(event);
1531 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1539 //_____________________________________________________________________________
1540 void AliReconstruction::Init(TTree *tree)
1543 AliError("The input tree is not found!");
1549 //_____________________________________________________________________________
1550 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1552 // run the reconstruction over a single event
1553 // The event loop is steered in Run method
1555 AliCodeTimerAuto("");
1557 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1558 fRunLoader->SetEventNumber(iEvent);
1559 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1561 fRunLoader->TreeE()->Fill();
1562 if (fRawReader && fRawReader->UseAutoSaveESD())
1563 fRunLoader->TreeE()->AutoSave("SaveSelf");
1566 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1570 AliInfo(Form("processing event %d", iEvent));
1572 fRunLoader->GetEvent(iEvent);
1574 // Fill Event-info object
1576 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo,fListOfCosmicTriggers);
1577 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1579 // Set the reco-params
1581 TString detStr = fLoadCDB;
1582 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1583 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1584 AliReconstructor *reconstructor = GetReconstructor(iDet);
1585 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1586 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1587 reconstructor->SetRecoParam(par);
1588 reconstructor->SetEventInfo(&fEventInfo);
1590 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1591 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(par->GetEventSpecie())) ;
1599 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1600 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1602 // local single event reconstruction
1603 if (!fRunLocalReconstruction.IsNull()) {
1604 TString detectors=fRunLocalReconstruction;
1605 // run HLT event reconstruction first
1606 // ;-( IsSelected changes the string
1607 if (IsSelected("HLT", detectors) &&
1608 !RunLocalEventReconstruction("HLT")) {
1609 if (fStopOnError) {CleanUp(); return kFALSE;}
1611 detectors=fRunLocalReconstruction;
1612 detectors.ReplaceAll("HLT", "");
1613 if (!RunLocalEventReconstruction(detectors)) {
1614 if (fStopOnError) {CleanUp(); return kFALSE;}
1618 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1619 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1620 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1621 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1623 // Set magnetic field from the tracker
1624 fesd->SetMagneticField(AliTracker::GetBz());
1625 fhltesd->SetMagneticField(AliTracker::GetBz());
1627 AliMagF* fld = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1628 if (fld) { // set info needed for field initialization
1629 fesd->SetCurrentL3(fld->GetCurrentSol());
1630 fesd->SetCurrentDip(fld->GetCurrentDip());
1631 fesd->SetBeamEnergy(fld->GetBeamEnergy());
1632 fesd->SetBeamType(fld->GetBeamTypeText());
1633 fesd->SetUniformBMap(fld->IsUniform());
1634 fesd->SetBInfoStored();
1636 fhltesd->SetCurrentL3(fld->GetCurrentSol());
1637 fhltesd->SetCurrentDip(fld->GetCurrentDip());
1638 fhltesd->SetBeamEnergy(fld->GetBeamEnergy());
1639 fhltesd->SetBeamType(fld->GetBeamTypeText());
1640 fhltesd->SetUniformBMap(fld->IsUniform());
1641 fhltesd->SetBInfoStored();
1644 // Set most probable pt, for B=0 tracking
1645 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1646 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1647 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1649 // Fill raw-data error log into the ESD
1650 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1653 if (fRunVertexFinder) {
1654 if (!RunVertexFinder(fesd)) {
1655 if (fStopOnError) {CleanUp(); return kFALSE;}
1659 // For Plane Efficiency: run the SPD trackleter
1660 if (fRunPlaneEff && fSPDTrackleter) {
1661 if (!RunSPDTrackleting(fesd)) {
1662 if (fStopOnError) {CleanUp(); return kFALSE;}
1667 if (!fRunTracking.IsNull()) {
1668 if (fRunMuonTracking) {
1669 if (!RunMuonTracking(fesd)) {
1670 if (fStopOnError) {CleanUp(); return kFALSE;}
1676 if (!fRunTracking.IsNull()) {
1677 if (!RunTracking(fesd)) {
1678 if (fStopOnError) {CleanUp(); return kFALSE;}
1683 if (!fFillESD.IsNull()) {
1684 TString detectors=fFillESD;
1685 // run HLT first and on hltesd
1686 // ;-( IsSelected changes the string
1687 if (IsSelected("HLT", detectors) &&
1688 !FillESD(fhltesd, "HLT")) {
1689 if (fStopOnError) {CleanUp(); return kFALSE;}
1692 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1693 if (detectors.Contains("ALL")) {
1695 for (Int_t idet=0; idet<kNDetectors; ++idet){
1696 detectors += fgkDetectorName[idet];
1700 detectors.ReplaceAll("HLT", "");
1701 if (!FillESD(fesd, detectors)) {
1702 if (fStopOnError) {CleanUp(); return kFALSE;}
1706 // fill Event header information from the RawEventHeader
1707 if (fRawReader){FillRawEventHeaderESD(fesd);}
1708 if (fRawReader){FillRawEventHeaderESD(fhltesd);}
1711 AliESDpid::MakePID(fesd);
1713 if (fFillTriggerESD) {
1714 if (!FillTriggerESD(fesd)) {
1715 if (fStopOnError) {CleanUp(); return kFALSE;}
1718 // Always fill scalers
1719 if (!FillTriggerScalers(fesd)) {
1720 if (fStopOnError) {CleanUp(); return kFALSE;}
1727 // Propagate track to the beam pipe (if not already done by ITS)
1729 const Int_t ntracks = fesd->GetNumberOfTracks();
1730 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1733 UShort_t *selectedIdx=new UShort_t[ntracks];
1735 for (Int_t itrack=0; itrack<ntracks; itrack++){
1736 const Double_t kMaxStep = 1; //max step over the material
1739 AliESDtrack *track = fesd->GetTrack(itrack);
1740 if (!track) continue;
1742 AliExternalTrackParam *tpcTrack =
1743 (AliExternalTrackParam *)track->GetTPCInnerParam();
1747 PropagateTrackToBxByBz(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1750 Int_t n=trkArray.GetEntriesFast();
1751 selectedIdx[n]=track->GetID();
1752 trkArray.AddLast(tpcTrack);
1755 //Tracks refitted by ITS should already be at the SPD vertex
1756 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1759 PropagateTrackToBxByBz(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1760 Double_t x[3]; track->GetXYZ(x);
1761 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1762 track->RelateToVertexBxByBz(fesd->GetPrimaryVertexSPD(), b, kVeryBig);
1767 // Improve the reconstructed primary vertex position using the tracks
1769 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1770 if(fesd->GetPrimaryVertexSPD()) {
1771 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1772 if(vtitle.Contains("cosmics")) {
1773 runVertexFinderTracks=kFALSE;
1777 if (runVertexFinderTracks) {
1778 // TPC + ITS primary vertex
1779 ftVertexer->SetITSMode();
1780 ftVertexer->SetConstraintOff();
1781 // get cuts for vertexer from AliGRPRecoParam
1783 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1784 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1785 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1786 ftVertexer->SetCuts(cutsVertexer);
1787 delete [] cutsVertexer; cutsVertexer = NULL;
1788 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1789 ftVertexer->SetVtxStart(fDiamondProfile);
1791 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1793 if (pvtx->GetStatus()) {
1794 fesd->SetPrimaryVertexTracks(pvtx);
1795 for (Int_t i=0; i<ntracks; i++) {
1796 AliESDtrack *t = fesd->GetTrack(i);
1797 Double_t x[3]; t->GetXYZ(x);
1798 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1799 t->RelateToVertexBxByBz(pvtx, b, kVeryBig);
1804 // TPC-only primary vertex
1805 ftVertexer->SetTPCMode();
1806 ftVertexer->SetConstraintOff();
1807 // get cuts for vertexer from AliGRPRecoParam
1809 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1810 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1811 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1812 ftVertexer->SetCuts(cutsVertexer);
1813 delete [] cutsVertexer; cutsVertexer = NULL;
1814 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1815 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1817 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1819 if (pvtx->GetStatus()) {
1820 fesd->SetPrimaryVertexTPC(pvtx);
1821 for (Int_t i=0; i<ntracks; i++) {
1822 AliESDtrack *t = fesd->GetTrack(i);
1823 Double_t x[3]; t->GetXYZ(x);
1824 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1825 t->RelateToVertexTPCBxByBz(pvtx, b, kVeryBig);
1831 delete[] selectedIdx;
1833 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1838 AliV0vertexer vtxer;
1839 vtxer.Tracks2V0vertices(fesd);
1841 if (fRunCascadeFinder) {
1843 AliCascadeVertexer cvtxer;
1844 cvtxer.V0sTracks2CascadeVertices(fesd);
1849 if (fCleanESD) CleanESD(fesd);
1852 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1853 AliQAManager::QAManager()->RunOneEvent(fesd) ;
1856 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1857 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1858 if (qadm && fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1859 qadm->Exec(AliQAv1::kESDS, fesd);
1862 if (fWriteESDfriend) {
1863 // fesdf->~AliESDfriend();
1864 // new (fesdf) AliESDfriend(); // Reset...
1865 fesd->GetESDfriend(fesdf);
1869 // Auto-save the ESD tree in case of prompt reco @P2
1870 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1871 ftree->AutoSave("SaveSelf");
1872 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1873 if (friendfile) friendfile->Save();
1880 if (fRunAliEVE) RunAliEVE();
1884 if (fWriteESDfriend) {
1885 fesdf->~AliESDfriend();
1886 new (fesdf) AliESDfriend(); // Reset...
1889 ProcInfo_t procInfo;
1890 gSystem->GetProcInfo(&procInfo);
1891 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1894 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1895 if (fReconstructor[iDet]) {
1896 fReconstructor[iDet]->SetRecoParam(NULL);
1897 fReconstructor[iDet]->SetEventInfo(NULL);
1899 if (fTracker[iDet]) fTracker[iDet]->SetEventInfo(NULL);
1902 if (fRunQA || fRunGlobalQA)
1903 AliQAManager::QAManager()->Increment() ;
1908 //_____________________________________________________________________________
1909 void AliReconstruction::SlaveTerminate()
1911 // Finalize the run on the slave side
1912 // Called after the exit
1913 // from the event loop
1914 AliCodeTimerAuto("");
1916 if (fIsNewRunLoader) { // galice.root didn't exist
1917 fRunLoader->WriteHeader("OVERWRITE");
1918 fRunLoader->CdGAFile();
1919 fRunLoader->Write(0, TObject::kOverwrite);
1922 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1923 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1925 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1926 cdbMapCopy->SetOwner(1);
1927 cdbMapCopy->SetName("cdbMap");
1928 TIter iter(cdbMap->GetTable());
1931 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1932 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1933 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1934 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1937 TList *cdbListCopy = new TList();
1938 cdbListCopy->SetOwner(1);
1939 cdbListCopy->SetName("cdbList");
1941 TIter iter2(cdbList);
1944 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1945 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1948 ftree->GetUserInfo()->Add(cdbMapCopy);
1949 ftree->GetUserInfo()->Add(cdbListCopy);
1954 if (fWriteESDfriend)
1955 ftree->SetBranchStatus("ESDfriend*",0);
1956 // we want to have only one tree version number
1957 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1958 fhlttree->Write(fhlttree->GetName(),TObject::kOverwrite);
1960 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1961 if (fRunPlaneEff && !FinishPlaneEff()) {
1962 AliWarning("Finish PlaneEff evaluation failed");
1965 // End of cycle for the in-loop
1967 AliQAManager::QAManager()->EndOfCycle() ;
1970 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1972 if (fQATasks.Contains(Form("%d", AliQAv1::kRECPOINTS)))
1973 qadm->EndOfCycle(AliQAv1::kRECPOINTS);
1974 if (fQATasks.Contains(Form("%d", AliQAv1::kESDS)))
1975 qadm->EndOfCycle(AliQAv1::kESDS);
1980 if (fRunQA || fRunGlobalQA) {
1982 !fProofOutputLocation.IsNull() &&
1983 fProofOutputArchive.IsNull() &&
1984 !fProofOutputDataset) {
1985 TString qaOutputFile(Form("%sMerged.%s.Data.root",
1986 fProofOutputLocation.Data(),
1987 AliQAv1::GetQADataFileName()));
1988 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",
1989 AliQAv1::GetQADataFileName()));
1990 qaProofFile->SetOutputFileName(qaOutputFile.Data());
1991 if (AliDebugLevel() > 0) qaProofFile->Dump();
1992 fOutput->Add(qaProofFile);
1993 MergeQA(qaProofFile->GetFileName());
2004 if (!fProofOutputFileName.IsNull() &&
2005 !fProofOutputLocation.IsNull() &&
2006 fProofOutputDataset &&
2007 !fProofOutputArchive.IsNull()) {
2008 TProofOutputFile *zipProofFile = new TProofOutputFile(fProofOutputFileName.Data(),
2010 fProofOutputLocation.Data());
2011 if (AliDebugLevel() > 0) zipProofFile->Dump();
2012 fOutput->Add(zipProofFile);
2013 TString fileList(fProofOutputArchive.Data());
2014 fileList.ReplaceAll(","," ");
2015 AliInfo(Form("Executing: zip -n root %s %s",zipProofFile->GetFileName(),fileList.Data()));
2016 gSystem->Exec(Form("zip -n root %s %s",zipProofFile->GetFileName(),fileList.Data()));
2021 //_____________________________________________________________________________
2022 void AliReconstruction::Terminate()
2024 // Create tags for the events in the ESD tree (the ESD tree is always present)
2025 // In case of empty events the tags will contain dummy values
2026 AliCodeTimerAuto("");
2028 // Do not call the ESD tag creator in case of PROOF-based reconstruction
2030 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
2031 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
2034 // Cleanup of CDB manager: cache and active storages!
2035 AliCDBManager::Instance()->ClearCache();
2038 //_____________________________________________________________________________
2039 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
2041 // run the local reconstruction
2043 static Int_t eventNr=0;
2044 AliCodeTimerAuto("")
2046 TString detStr = detectors;
2047 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2048 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2049 AliReconstructor* reconstructor = GetReconstructor(iDet);
2050 if (!reconstructor) continue;
2051 AliLoader* loader = fLoader[iDet];
2052 // Matthias April 2008: temporary fix to run HLT reconstruction
2053 // although the HLT loader is missing
2054 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2056 reconstructor->Reconstruct(fRawReader, NULL);
2059 reconstructor->Reconstruct(dummy, NULL);
2064 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2067 // conversion of digits
2068 if (fRawReader && reconstructor->HasDigitConversion()) {
2069 AliInfo(Form("converting raw data digits into root objects for %s",
2070 fgkDetectorName[iDet]));
2071 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2072 // fgkDetectorName[iDet]));
2073 loader->LoadDigits("update");
2074 loader->CleanDigits();
2075 loader->MakeDigitsContainer();
2076 TTree* digitsTree = loader->TreeD();
2077 reconstructor->ConvertDigits(fRawReader, digitsTree);
2078 loader->WriteDigits("OVERWRITE");
2079 loader->UnloadDigits();
2081 // local reconstruction
2082 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2083 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2084 loader->LoadRecPoints("update");
2085 loader->CleanRecPoints();
2086 loader->MakeRecPointsContainer();
2087 TTree* clustersTree = loader->TreeR();
2088 if (fRawReader && !reconstructor->HasDigitConversion()) {
2089 reconstructor->Reconstruct(fRawReader, clustersTree);
2091 loader->LoadDigits("read");
2092 TTree* digitsTree = loader->TreeD();
2094 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2095 if (fStopOnError) return kFALSE;
2097 reconstructor->Reconstruct(digitsTree, clustersTree);
2099 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2100 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2103 loader->UnloadDigits();
2106 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2107 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2109 loader->WriteRecPoints("OVERWRITE");
2110 loader->UnloadRecPoints();
2111 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2113 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2114 AliError(Form("the following detectors were not found: %s",
2116 if (fStopOnError) return kFALSE;
2121 //_____________________________________________________________________________
2122 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2124 // run the SPD trackleting (for SPD efficiency purpouses)
2126 AliCodeTimerAuto("")
2128 Double_t vtxPos[3] = {0, 0, 0};
2129 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2131 TArrayF mcVertex(3);
2133 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2134 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2135 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2138 const AliESDVertex *vertex = esd->GetVertex();
2140 AliWarning("Vertex not found");
2143 vertex->GetXYZ(vtxPos);
2144 vertex->GetSigmaXYZ(vtxErr);
2145 if (fSPDTrackleter) {
2146 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2149 fLoader[0]->LoadRecPoints("read");
2150 TTree* tree = fLoader[0]->TreeR();
2152 AliError("Can't get the ITS cluster tree");
2155 fSPDTrackleter->LoadClusters(tree);
2156 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2158 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2159 AliError("AliITSTrackleterSPDEff Clusters2Tracks failed");
2160 // fLoader[0]->UnloadRecPoints();
2163 //fSPDTrackleter->UnloadRecPoints();
2165 AliWarning("SPDTrackleter not available");
2171 //_____________________________________________________________________________
2172 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2174 // run the barrel tracking
2176 AliCodeTimerAuto("")
2178 AliVertexer *vertexer = CreateVertexer();
2179 if (!vertexer) return kFALSE;
2181 AliInfo("running the ITS vertex finder");
2182 AliESDVertex* vertex = NULL;
2184 fLoader[0]->LoadRecPoints();
2185 TTree* cltree = fLoader[0]->TreeR();
2187 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2188 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2191 AliError("Can't get the ITS cluster tree");
2193 fLoader[0]->UnloadRecPoints();
2196 AliError("Can't get the ITS loader");
2199 AliWarning("Vertex not found");
2200 vertex = new AliESDVertex();
2201 vertex->SetName("default");
2204 vertex->SetName("reconstructed");
2209 vertex->GetXYZ(vtxPos);
2210 vertex->GetSigmaXYZ(vtxErr);
2212 esd->SetPrimaryVertexSPD(vertex);
2213 AliESDVertex *vpileup = NULL;
2214 Int_t novertices = 0;
2215 vpileup = vertexer->GetAllVertices(novertices);
2217 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2219 // if SPD multiplicity has been determined, it is stored in the ESD
2220 AliMultiplicity *mult = vertexer->GetMultiplicity();
2221 if(mult)esd->SetMultiplicity(mult);
2223 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2224 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2233 //_____________________________________________________________________________
2234 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2236 // run the HLT barrel tracking
2238 AliCodeTimerAuto("")
2241 AliError("Missing runLoader!");
2245 AliInfo("running HLT tracking");
2247 // Get a pointer to the HLT reconstructor
2248 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2249 if (!reconstructor) return kFALSE;
2252 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2253 TString detName = fgkDetectorName[iDet];
2254 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2255 reconstructor->SetOption(detName.Data());
2256 AliTracker *tracker = reconstructor->CreateTracker();
2258 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2259 if (fStopOnError) return kFALSE;
2263 Double_t vtxErr[3]={0.005,0.005,0.010};
2264 const AliESDVertex *vertex = esd->GetVertex();
2265 vertex->GetXYZ(vtxPos);
2266 tracker->SetVertex(vtxPos,vtxErr);
2268 fLoader[iDet]->LoadRecPoints("read");
2269 TTree* tree = fLoader[iDet]->TreeR();
2271 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2274 tracker->LoadClusters(tree);
2276 if (tracker->Clusters2Tracks(esd) != 0) {
2277 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2281 tracker->UnloadClusters();
2289 //_____________________________________________________________________________
2290 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2292 // run the muon spectrometer tracking
2294 AliCodeTimerAuto("")
2297 AliError("Missing runLoader!");
2300 Int_t iDet = 7; // for MUON
2302 AliInfo("is running...");
2304 // Get a pointer to the MUON reconstructor
2305 AliReconstructor *reconstructor = GetReconstructor(iDet);
2306 if (!reconstructor) return kFALSE;
2309 TString detName = fgkDetectorName[iDet];
2310 AliDebug(1, Form("%s tracking", detName.Data()));
2311 AliTracker *tracker = reconstructor->CreateTracker();
2313 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2318 fLoader[iDet]->LoadRecPoints("read");
2320 tracker->LoadClusters(fLoader[iDet]->TreeR());
2322 Int_t rv = tracker->Clusters2Tracks(esd);
2326 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2330 fLoader[iDet]->UnloadRecPoints();
2332 tracker->UnloadClusters();
2340 //_____________________________________________________________________________
2341 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2343 // run the barrel tracking
2344 static Int_t eventNr=0;
2345 AliCodeTimerAuto("")
2347 AliInfo("running tracking");
2349 // Set the event info which is used
2350 // by the trackers in order to obtain
2351 // information about read-out detectors,
2353 AliDebug(1, "Setting event info");
2354 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2355 if (!fTracker[iDet]) continue;
2356 fTracker[iDet]->SetEventInfo(&fEventInfo);
2359 //Fill the ESD with the T0 info (will be used by the TOF)
2360 if (fReconstructor[11] && fLoader[11]) {
2361 fLoader[11]->LoadRecPoints("READ");
2362 TTree *treeR = fLoader[11]->TreeR();
2364 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2368 // pass 1: TPC + ITS inwards
2369 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2370 if (!fTracker[iDet]) continue;
2371 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2374 fLoader[iDet]->LoadRecPoints("read");
2375 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2376 TTree* tree = fLoader[iDet]->TreeR();
2378 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2381 fTracker[iDet]->LoadClusters(tree);
2382 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2384 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2385 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2388 // preliminary PID in TPC needed by the ITS tracker
2390 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2391 AliESDpid::MakePID(esd);
2393 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2396 // pass 2: ALL backwards
2398 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2399 if (!fTracker[iDet]) continue;
2400 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2403 if (iDet > 1) { // all except ITS, TPC
2405 fLoader[iDet]->LoadRecPoints("read");
2406 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2407 tree = fLoader[iDet]->TreeR();
2409 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2412 fTracker[iDet]->LoadClusters(tree);
2413 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2417 if (iDet>1) // start filling residuals for the "outer" detectors
2419 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2420 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2422 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2423 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2424 if ( elem && (! elem->At(0)) ) {
2425 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2426 if (qadm) qadm->InitRecPointsForTracker() ;
2430 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2431 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2436 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2437 fTracker[iDet]->UnloadClusters();
2438 fLoader[iDet]->UnloadRecPoints();
2440 // updated PID in TPC needed by the ITS tracker -MI
2442 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2443 AliESDpid::MakePID(esd);
2445 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2447 //stop filling residuals for the "outer" detectors
2448 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2450 // pass 3: TRD + TPC + ITS refit inwards
2452 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2453 if (!fTracker[iDet]) continue;
2454 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2457 if (iDet<2) // start filling residuals for TPC and ITS
2459 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2460 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2462 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2463 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2464 if ( elem && (! elem->At(0)) ) {
2465 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2466 if (qadm) qadm->InitRecPointsForTracker() ;
2471 if (fTracker[iDet]->RefitInward(esd) != 0) {
2472 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2475 // run postprocessing
2476 if (fTracker[iDet]->PostProcess(esd) != 0) {
2477 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2480 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2483 // write space-points to the ESD in case alignment data output
2485 if (fWriteAlignmentData)
2486 WriteAlignmentData(esd);
2488 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2489 if (!fTracker[iDet]) continue;
2491 fTracker[iDet]->UnloadClusters();
2492 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2493 fLoader[iDet]->UnloadRecPoints();
2494 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2496 // stop filling residuals for TPC and ITS
2497 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2503 //_____________________________________________________________________________
2504 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2506 // Remove the data which are not needed for the physics analysis.
2509 Int_t nTracks=esd->GetNumberOfTracks();
2510 Int_t nV0s=esd->GetNumberOfV0s();
2512 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2514 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2515 Bool_t rc=esd->Clean(cleanPars);
2517 nTracks=esd->GetNumberOfTracks();
2518 nV0s=esd->GetNumberOfV0s();
2520 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2525 //_____________________________________________________________________________
2526 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2528 // fill the event summary data
2530 AliCodeTimerAuto("")
2531 static Int_t eventNr=0;
2532 TString detStr = detectors;
2534 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2535 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2536 AliReconstructor* reconstructor = GetReconstructor(iDet);
2537 if (!reconstructor) continue;
2538 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2539 TTree* clustersTree = NULL;
2540 if (fLoader[iDet]) {
2541 fLoader[iDet]->LoadRecPoints("read");
2542 clustersTree = fLoader[iDet]->TreeR();
2543 if (!clustersTree) {
2544 AliError(Form("Can't get the %s clusters tree",
2545 fgkDetectorName[iDet]));
2546 if (fStopOnError) return kFALSE;
2549 if (fRawReader && !reconstructor->HasDigitConversion()) {
2550 reconstructor->FillESD(fRawReader, clustersTree, esd);
2552 TTree* digitsTree = NULL;
2553 if (fLoader[iDet]) {
2554 fLoader[iDet]->LoadDigits("read");
2555 digitsTree = fLoader[iDet]->TreeD();
2557 AliError(Form("Can't get the %s digits tree",
2558 fgkDetectorName[iDet]));
2559 if (fStopOnError) return kFALSE;
2562 reconstructor->FillESD(digitsTree, clustersTree, esd);
2563 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2565 if (fLoader[iDet]) {
2566 fLoader[iDet]->UnloadRecPoints();
2570 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2571 AliError(Form("the following detectors were not found: %s",
2573 if (fStopOnError) return kFALSE;
2575 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2580 //_____________________________________________________________________________
2581 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2583 // Reads the trigger decision which is
2584 // stored in Trigger.root file and fills
2585 // the corresponding esd entries
2587 AliCodeTimerAuto("")
2589 AliInfo("Filling trigger information into the ESD");
2592 AliCTPRawStream input(fRawReader);
2593 if (!input.Next()) {
2594 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2597 if (esd->GetTriggerMask() != input.GetClassMask())
2598 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2599 input.GetClassMask(),esd->GetTriggerMask()));
2600 if (esd->GetOrbitNumber() != input.GetOrbitID())
2601 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2602 input.GetOrbitID(),esd->GetOrbitNumber()));
2603 if (esd->GetBunchCrossNumber() != input.GetBCID())
2604 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2605 input.GetBCID(),esd->GetBunchCrossNumber()));
2606 AliESDHeader* esdheader = esd->GetHeader();
2607 esdheader->SetL0TriggerInputs(input.GetL0Inputs());
2608 esdheader->SetL1TriggerInputs(input.GetL1Inputs());
2609 esdheader->SetL2TriggerInputs(input.GetL2Inputs());
2611 UInt_t orbit=input.GetOrbitID();
2612 for(Int_t i=0 ; i<input.GetNIRs() ; i++ )
2613 if(TMath::Abs(Int_t(orbit-(input.GetIR(i))->GetOrbit()))<=1){
2614 esdheader->AddTriggerIR(input.GetIR(i));
2620 //_____________________________________________________________________________
2621 Bool_t AliReconstruction::FillTriggerScalers(AliESDEvent*& esd)
2624 //fRunScalers->Print();
2625 if(fRunScalers && fRunScalers->CheckRunScalers()){
2626 AliTimeStamp* timestamp = new AliTimeStamp(esd->GetOrbitNumber(), esd->GetPeriodNumber(), esd->GetBunchCrossNumber());
2627 //AliTimeStamp* timestamp = new AliTimeStamp(10308000, 0, (ULong64_t)486238);
2628 AliESDHeader* esdheader = fesd->GetHeader();
2629 for(Int_t i=0;i<50;i++){
2630 if((1<<i) & esd->GetTriggerMask()){
2631 AliTriggerScalersESD* scalesd = fRunScalers->GetScalersForEventClass( timestamp, i+1);
2632 if(scalesd)esdheader->SetTriggerScalersRecord(scalesd);
2638 //_____________________________________________________________________________
2639 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2642 // Filling information from RawReader Header
2645 if (!fRawReader) return kFALSE;
2647 AliInfo("Filling information from RawReader Header");
2649 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2650 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2651 esd->SetPeriodNumber(fRawReader->GetPeriod());
2653 esd->SetTimeStamp(fRawReader->GetTimestamp());
2654 esd->SetEventType(fRawReader->GetType());
2660 //_____________________________________________________________________________
2661 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2663 // check whether detName is contained in detectors
2664 // if yes, it is removed from detectors
2666 // check if all detectors are selected
2667 if ((detectors.CompareTo("ALL") == 0) ||
2668 detectors.BeginsWith("ALL ") ||
2669 detectors.EndsWith(" ALL") ||
2670 detectors.Contains(" ALL ")) {
2675 // search for the given detector
2676 Bool_t result = kFALSE;
2677 if ((detectors.CompareTo(detName) == 0) ||
2678 detectors.BeginsWith(detName+" ") ||
2679 detectors.EndsWith(" "+detName) ||
2680 detectors.Contains(" "+detName+" ")) {
2681 detectors.ReplaceAll(detName, "");
2685 // clean up the detectors string
2686 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2687 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2688 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2693 //_____________________________________________________________________________
2694 Bool_t AliReconstruction::InitRunLoader()
2696 // get or create the run loader
2698 if (gAlice) delete gAlice;
2701 TFile *gafile = TFile::Open(fGAliceFileName.Data());
2702 // if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2703 if (gafile) { // galice.root exists
2707 // load all base libraries to get the loader classes
2708 TString libs = gSystem->GetLibraries();
2709 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2710 TString detName = fgkDetectorName[iDet];
2711 if (detName == "HLT") continue;
2712 if (libs.Contains("lib" + detName + "base.so")) continue;
2713 gSystem->Load("lib" + detName + "base.so");
2715 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2717 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2722 fRunLoader->CdGAFile();
2723 fRunLoader->LoadgAlice();
2725 //PH This is a temporary fix to give access to the kinematics
2726 //PH that is needed for the labels of ITS clusters
2727 fRunLoader->LoadHeader();
2728 fRunLoader->LoadKinematics();
2730 } else { // galice.root does not exist
2732 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2734 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2735 AliConfig::GetDefaultEventFolderName(),
2738 AliError(Form("could not create run loader in file %s",
2739 fGAliceFileName.Data()));
2743 fIsNewRunLoader = kTRUE;
2744 fRunLoader->MakeTree("E");
2746 if (fNumberOfEventsPerFile > 0)
2747 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2749 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2755 //_____________________________________________________________________________
2756 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2758 // get the reconstructor object and the loader for a detector
2760 if (fReconstructor[iDet]) {
2761 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2762 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2763 fReconstructor[iDet]->SetRecoParam(par);
2764 fReconstructor[iDet]->SetRunInfo(fRunInfo);
2766 return fReconstructor[iDet];
2769 // load the reconstructor object
2770 TPluginManager* pluginManager = gROOT->GetPluginManager();
2771 TString detName = fgkDetectorName[iDet];
2772 TString recName = "Ali" + detName + "Reconstructor";
2774 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2776 AliReconstructor* reconstructor = NULL;
2777 // first check if a plugin is defined for the reconstructor
2778 TPluginHandler* pluginHandler =
2779 pluginManager->FindHandler("AliReconstructor", detName);
2780 // if not, add a plugin for it
2781 if (!pluginHandler) {
2782 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2783 TString libs = gSystem->GetLibraries();
2784 if (libs.Contains("lib" + detName + "base.so") ||
2785 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2786 pluginManager->AddHandler("AliReconstructor", detName,
2787 recName, detName + "rec", recName + "()");
2789 pluginManager->AddHandler("AliReconstructor", detName,
2790 recName, detName, recName + "()");
2792 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2794 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2795 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2797 if (reconstructor) {
2798 TObject* obj = fOptions.FindObject(detName.Data());
2799 if (obj) reconstructor->SetOption(obj->GetTitle());
2800 reconstructor->SetRunInfo(fRunInfo);
2801 reconstructor->Init();
2802 fReconstructor[iDet] = reconstructor;
2805 // get or create the loader
2806 if (detName != "HLT") {
2807 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2808 if (!fLoader[iDet]) {
2809 AliConfig::Instance()
2810 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2812 // first check if a plugin is defined for the loader
2814 pluginManager->FindHandler("AliLoader", detName);
2815 // if not, add a plugin for it
2816 if (!pluginHandler) {
2817 TString loaderName = "Ali" + detName + "Loader";
2818 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2819 pluginManager->AddHandler("AliLoader", detName,
2820 loaderName, detName + "base",
2821 loaderName + "(const char*, TFolder*)");
2822 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2824 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2826 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2827 fRunLoader->GetEventFolder());
2829 if (!fLoader[iDet]) { // use default loader
2830 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2832 if (!fLoader[iDet]) {
2833 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2834 if (fStopOnError) return NULL;
2836 fRunLoader->AddLoader(fLoader[iDet]);
2837 fRunLoader->CdGAFile();
2838 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2839 fRunLoader->Write(0, TObject::kOverwrite);
2844 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2845 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2846 reconstructor->SetRecoParam(par);
2847 reconstructor->SetRunInfo(fRunInfo);
2849 return reconstructor;
2852 //_____________________________________________________________________________
2853 AliVertexer* AliReconstruction::CreateVertexer()
2855 // create the vertexer
2856 // Please note that the caller is the owner of the
2859 AliVertexer* vertexer = NULL;
2860 AliReconstructor* itsReconstructor = GetReconstructor(0);
2861 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) || fRunTracking.Contains("ITS"))) {
2862 vertexer = itsReconstructor->CreateVertexer();
2865 AliWarning("couldn't create a vertexer for ITS");
2871 //_____________________________________________________________________________
2872 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2874 // create the trackers
2875 AliInfo("Creating trackers");
2877 TString detStr = detectors;
2878 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2879 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2880 AliReconstructor* reconstructor = GetReconstructor(iDet);
2881 if (!reconstructor) continue;
2882 TString detName = fgkDetectorName[iDet];
2883 if (detName == "HLT") {
2884 fRunHLTTracking = kTRUE;
2887 if (detName == "MUON") {
2888 fRunMuonTracking = kTRUE;
2893 fTracker[iDet] = reconstructor->CreateTracker();
2894 if (!fTracker[iDet] && (iDet < 7)) {
2895 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2896 if (fStopOnError) return kFALSE;
2898 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2904 //_____________________________________________________________________________
2905 void AliReconstruction::CleanUp()
2907 // delete trackers and the run loader and close and delete the file
2909 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2910 delete fReconstructor[iDet];
2911 fReconstructor[iDet] = NULL;
2912 fLoader[iDet] = NULL;
2913 delete fTracker[iDet];
2914 fTracker[iDet] = NULL;
2919 delete fSPDTrackleter;
2920 fSPDTrackleter = NULL;
2929 delete fParentRawReader;
2930 fParentRawReader=NULL;
2938 if (AliQAManager::QAManager())
2939 AliQAManager::QAManager()->ShowQA() ;
2940 AliQAManager::Destroy() ;
2942 TGeoGlobalMagField::Instance()->SetField(NULL);
2945 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2947 // Write space-points which are then used in the alignment procedures
2948 // For the moment only ITS, TPC, TRD and TOF
2950 Int_t ntracks = esd->GetNumberOfTracks();
2951 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2953 AliESDtrack *track = esd->GetTrack(itrack);
2956 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
2957 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2958 nsp += track->GetNcls(iDet);
2960 if (iDet==0) { // ITS "extra" clusters
2961 track->GetClusters(iDet,idx);
2962 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2967 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2968 track->SetTrackPointArray(sp);
2970 for (Int_t iDet = 5; iDet >= 0; iDet--) {
2971 AliTracker *tracker = fTracker[iDet];
2972 if (!tracker) continue;
2973 Int_t nspdet = track->GetClusters(iDet,idx);
2975 if (iDet==0) // ITS "extra" clusters
2976 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2978 if (nspdet <= 0) continue;
2982 while (isp2 < nspdet) {
2983 Bool_t isvalid=kTRUE;
2985 Int_t index=idx[isp++];
2986 if (index < 0) continue;
2988 TString dets = fgkDetectorName[iDet];
2989 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2990 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2991 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2992 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2993 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2995 isvalid = tracker->GetTrackPoint(index,p);
2998 if (!isvalid) continue;
2999 if (iDet==0 && (isp-1)>=6) p.SetExtra();
3000 sp->AddPoint(isptrack,&p); isptrack++;
3007 //_____________________________________________________________________________
3008 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
3010 // The method reads the raw-data error log
3011 // accumulated within the rawReader.
3012 // It extracts the raw-data errors related to
3013 // the current event and stores them into
3014 // a TClonesArray inside the esd object.
3016 if (!fRawReader) return;
3018 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
3020 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
3022 if (iEvent != log->GetEventNumber()) continue;
3024 esd->AddRawDataErrorLog(log);
3029 //_____________________________________________________________________________
3030 void AliReconstruction::CheckQA()
3032 // check the QA of SIM for this run and remove the detectors
3033 // with status Fatal
3035 // TString newRunLocalReconstruction ;
3036 // TString newRunTracking ;
3037 // TString newFillESD ;
3039 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
3040 // TString detName(AliQAv1::GetDetName(iDet)) ;
3041 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
3042 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
3043 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
3044 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
3046 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
3047 // fRunLocalReconstruction.Contains("ALL") ) {
3048 // newRunLocalReconstruction += detName ;
3049 // newRunLocalReconstruction += " " ;
3051 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
3052 // fRunTracking.Contains("ALL") ) {
3053 // newRunTracking += detName ;
3054 // newRunTracking += " " ;
3056 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
3057 // fFillESD.Contains("ALL") ) {
3058 // newFillESD += detName ;
3059 // newFillESD += " " ;
3063 // fRunLocalReconstruction = newRunLocalReconstruction ;
3064 // fRunTracking = newRunTracking ;
3065 // fFillESD = newFillESD ;
3068 //_____________________________________________________________________________
3069 Int_t AliReconstruction::GetDetIndex(const char* detector)
3071 // return the detector index corresponding to detector
3073 for (index = 0; index < kNDetectors ; index++) {
3074 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3079 //_____________________________________________________________________________
3080 Bool_t AliReconstruction::FinishPlaneEff() {
3082 // Here execute all the necessary operationis, at the end of the tracking phase,
3083 // in case that evaluation of PlaneEfficiencies was required for some detector.
3084 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3086 // This Preliminary version works only FOR ITS !!!!!
3087 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3090 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3093 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3094 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3095 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3096 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
3097 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3098 TString name=planeeff->GetName();
3100 TFile* pefile = TFile::Open(name, "RECREATE");
3101 ret=(Bool_t)planeeff->Write();
3103 if(planeeff->GetCreateHistos()) {
3104 TString hname=planeeff->GetName();
3105 hname+="Histo.root";
3106 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3109 if(fSPDTrackleter) {
3110 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3111 TString name="AliITSPlaneEffSPDtracklet.root";
3112 TFile* pefile = TFile::Open(name, "RECREATE");
3113 ret=(Bool_t)planeeff->Write();
3115 AliESDEvent *dummy=NULL;
3116 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3121 //_____________________________________________________________________________
3122 Bool_t AliReconstruction::InitPlaneEff() {
3124 // Here execute all the necessary operations, before of the tracking phase,
3125 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3126 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3127 // which should be updated/recalculated.
3129 // This Preliminary version will work only FOR ITS !!!!!
3130 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3133 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3135 AliWarning(Form("Implementation of this method not yet completed !! Method return kTRUE"));
3137 fSPDTrackleter = NULL;
3138 AliReconstructor* itsReconstructor = GetReconstructor(0);
3139 if (itsReconstructor) {
3140 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3142 if (fSPDTrackleter) {
3143 AliInfo("Trackleter for SPD has been created");
3149 //_____________________________________________________________________________
3150 Bool_t AliReconstruction::InitAliEVE()
3152 // This method should be called only in case
3153 // AliReconstruction is run
3154 // within the alieve environment.
3155 // It will initialize AliEVE in a way
3156 // so that it can visualize event processed
3157 // by AliReconstruction.
3158 // The return flag shows whenever the
3159 // AliEVE initialization was successful or not.
3162 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3163 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3164 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3166 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3167 gROOT->ProcessLine("alieve_online_init()");
3172 //_____________________________________________________________________________
3173 void AliReconstruction::RunAliEVE()
3175 // Runs AliEVE visualisation of
3176 // the current event.
3177 // Should be executed only after
3178 // successful initialization of AliEVE.
3180 AliInfo("Running AliEVE...");
3181 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3185 //_____________________________________________________________________________
3186 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3188 // Allows to run QA for a selected set of detectors
3189 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3190 // all selected detectors run the same selected tasks
3192 if (!detAndAction.Contains(":")) {
3193 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3197 Int_t colon = detAndAction.Index(":") ;
3198 fQADetectors = detAndAction(0, colon) ;
3199 if (fQADetectors.Contains("ALL") )
3200 fQADetectors = fFillESD ;
3201 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3202 if (fQATasks.Contains("ALL") ) {
3203 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3205 fQATasks.ToUpper() ;
3207 if ( fQATasks.Contains("RAW") )
3208 tempo = Form("%d ", AliQAv1::kRAWS) ;
3209 if ( fQATasks.Contains("DIGIT") )
3210 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3211 if ( fQATasks.Contains("RECPOINT") )
3212 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3213 if ( fQATasks.Contains("ESD") )
3214 tempo += Form("%d ", AliQAv1::kESDS) ;
3216 if (fQATasks.IsNull()) {
3217 AliInfo("No QA requested\n") ;
3222 TString tempo(fQATasks) ;
3223 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3224 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3225 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3226 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3227 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3232 //_____________________________________________________________________________
3233 Bool_t AliReconstruction::InitRecoParams()
3235 // The method accesses OCDB and retrieves all
3236 // the available reco-param objects from there.
3238 Bool_t isOK = kTRUE;
3240 if (fRecoParam.GetDetRecoParamArray(kNDetectors)) {
3241 AliInfo("Using custom GRP reconstruction parameters");
3244 AliInfo("Loading GRP reconstruction parameter objects");
3246 AliCDBPath path("GRP","Calib","RecoParam");
3247 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3249 AliWarning("Couldn't find GRP RecoParam entry in OCDB");
3253 TObject *recoParamObj = entry->GetObject();
3254 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3255 // GRP has a normal TobjArray of AliDetectorRecoParam objects
3256 // Registering them in AliRecoParam
3257 fRecoParam.AddDetRecoParamArray(kNDetectors,dynamic_cast<TObjArray*>(recoParamObj));
3259 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3260 // GRP has only onse set of reco parameters
3261 // Registering it in AliRecoParam
3262 AliInfo("Single set of GRP reconstruction parameters found");
3263 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3264 fRecoParam.AddDetRecoParam(kNDetectors,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3267 AliError("No valid GRP RecoParam object found in the OCDB");
3274 TString detStr = fLoadCDB;
3275 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3277 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3279 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3280 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3284 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3286 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3287 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3289 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3293 TObject *recoParamObj = entry->GetObject();
3294 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3295 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3296 // Registering them in AliRecoParam
3297 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3299 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3300 // The detector has only onse set of reco parameters
3301 // Registering it in AliRecoParam
3302 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3303 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3304 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3307 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3311 // FIX ME: We have to disable the unloading of reco-param CDB
3312 // entries because QA framework is using them. Has to be fix in
3313 // a way that the QA takes the objects already constructed in
3315 // AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3319 if (AliDebugLevel() > 0) fRecoParam.Print();
3324 //_____________________________________________________________________________
3325 Bool_t AliReconstruction::GetEventInfo()
3327 // Fill the event info object
3329 AliCodeTimerAuto("")
3331 AliCentralTrigger *aCTP = NULL;
3333 fEventInfo.SetEventType(fRawReader->GetType());
3335 ULong64_t mask = fRawReader->GetClassMask();
3336 fEventInfo.SetTriggerMask(mask);
3337 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3338 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3340 aCTP = new AliCentralTrigger();
3341 TString configstr("");
3342 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3343 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3347 aCTP->SetClassMask(mask);
3348 aCTP->SetClusterMask(clmask);
3351 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3353 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3354 aCTP = fRunLoader->GetTrigger();
3355 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3356 // get inputs from actp - just get
3357 AliESDHeader* esdheader = fesd->GetHeader();
3358 esdheader->SetL0TriggerInputs(aCTP->GetL0TriggerInputs());
3359 esdheader->SetL1TriggerInputs(aCTP->GetL1TriggerInputs());
3360 esdheader->SetL2TriggerInputs(aCTP->GetL2TriggerInputs());
3361 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3364 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3369 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3371 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3372 if (fRawReader) delete aCTP;
3376 UChar_t clustmask = 0;
3378 ULong64_t trmask = fEventInfo.GetTriggerMask();
3379 const TObjArray& classesArray = config->GetClasses();
3380 Int_t nclasses = classesArray.GetEntriesFast();
3381 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3382 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3384 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3385 fesd->SetTriggerClass(trclass->GetName(),trindex);
3386 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3387 if (trmask & (1ull << trindex)) {
3389 trclasses += trclass->GetName();
3391 clustmask |= trclass->GetCluster()->GetClusterMask();
3395 fEventInfo.SetTriggerClasses(trclasses);
3397 // Set the information in ESD
3398 fesd->SetTriggerMask(trmask);
3399 fesd->SetTriggerCluster(clustmask);
3401 if (!aCTP->CheckTriggeredDetectors()) {
3402 if (fRawReader) delete aCTP;
3406 if (fRawReader) delete aCTP;
3408 // We have to fill also the HLT decision here!!
3414 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3416 // Match the detector list found in the rec.C or the default 'ALL'
3417 // to the list found in the GRP (stored there by the shuttle PP which
3418 // gets the information from ECS)
3419 static TString resultList;
3420 TString detList = detectorList;
3424 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3425 if ((detectorMask >> iDet) & 0x1) {
3426 TString det = AliDAQ::OfflineModuleName(iDet);
3427 if ((detList.CompareTo("ALL") == 0) ||
3428 ((detList.BeginsWith("ALL ") ||
3429 detList.EndsWith(" ALL") ||
3430 detList.Contains(" ALL ")) &&
3431 !(detList.BeginsWith("-"+det+" ") ||
3432 detList.EndsWith(" -"+det) ||
3433 detList.Contains(" -"+det+" "))) ||
3434 (detList.CompareTo(det) == 0) ||
3435 detList.BeginsWith(det+" ") ||
3436 detList.EndsWith(" "+det) ||
3437 detList.Contains( " "+det+" " )) {
3438 if (!resultList.EndsWith(det + " ")) {
3447 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3448 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3449 if ((detList.CompareTo("ALL") == 0) ||
3450 ((detList.BeginsWith("ALL ") ||
3451 detList.EndsWith(" ALL") ||
3452 detList.Contains(" ALL ")) &&
3453 !(detList.BeginsWith("-"+hltDet+" ") ||
3454 detList.EndsWith(" -"+hltDet) ||
3455 detList.Contains(" -"+hltDet+" "))) ||
3456 (detList.CompareTo(hltDet) == 0) ||
3457 detList.BeginsWith(hltDet+" ") ||
3458 detList.EndsWith(" "+hltDet) ||
3459 detList.Contains( " "+hltDet+" " )) {
3460 resultList += hltDet;
3464 return resultList.Data();
3468 //______________________________________________________________________________
3469 void AliReconstruction::Abort(const char *method, EAbort what)
3471 // Abort processing. If what = kAbortProcess, the Process() loop will be
3472 // aborted. If what = kAbortFile, the current file in a chain will be
3473 // aborted and the processing will continue with the next file, if there
3474 // is no next file then Process() will be aborted. Abort() can also be
3475 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3476 // the SlaveTerminate() and Terminate() are always called. The abort flag
3477 // can be checked in these methods using GetAbort().
3479 // The method is overwritten in AliReconstruction for better handling of
3480 // reco specific errors
3482 if (!fStopOnError) return;
3486 TString whyMess = method;
3487 whyMess += " failed! Aborting...";
3489 AliError(whyMess.Data());
3492 TString mess = "Abort";
3493 if (fAbort == kAbortProcess)
3494 mess = "AbortProcess";
3495 else if (fAbort == kAbortFile)
3498 Info(mess, whyMess.Data());
3501 //______________________________________________________________________________
3502 Bool_t AliReconstruction::ProcessEvent(void* event)
3504 // Method that is used in case the event loop
3505 // is steered from outside, for example by AMORE
3506 // 'event' is a pointer to the DATE event in the memory
3508 if (fRawReader) delete fRawReader;
3509 fRawReader = new AliRawReaderDate(event);
3510 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
3517 //______________________________________________________________________________
3518 Bool_t AliReconstruction::ParseOutput()
3520 // The method parses the output file
3521 // location string in order to steer
3522 // properly the selector
3524 TPMERegexp re1("(\\w+\\.zip#\\w+\\.root):([,*\\w+\\.root,*]+)@dataset://(\\w++)");
3525 TPMERegexp re2("(\\w+\\.root)?@?dataset://(\\w++)");
3527 if (re1.Match(fESDOutput) == 4) {
3528 // root archive with output files stored and regustered
3530 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",re1[1].Data()));
3531 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re1[3].Data()));
3532 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3533 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_ARCHIVE",re1[2].Data()));
3534 AliInfo(Form("%s files will be stored within %s in dataset %s",
3539 else if (re2.Match(fESDOutput) == 3) {
3540 // output file stored and registered
3542 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",(re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data()));
3543 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re2[2].Data()));
3544 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3545 AliInfo(Form("%s will be stored in dataset %s",
3546 (re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data(),
3550 if (fESDOutput.IsNull()) {
3551 // Output location not given.
3552 // Assuming xrootd has been already started and
3553 // the output file has to be sent back
3554 // to the client machine
3555 TString esdUrl(Form("root://%s/%s/",
3556 TUrl(gSystem->HostName()).GetHostFQDN(),
3558 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE","AliESDs.root"));
3559 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",esdUrl.Data()));
3560 AliInfo(Form("AliESDs.root will be stored in %s",
3564 // User specified an output location.
3565 // Ones has just to parse it here
3566 TUrl outputUrl(fESDOutput.Data());
3567 TString outputFile(gSystem->BaseName(outputUrl.GetFile()));
3568 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.IsNull() ? "AliESDs.root" : outputFile.Data()));
3569 TString outputLocation(outputUrl.GetUrl());
3570 outputLocation.ReplaceAll(outputFile.Data(),"");
3571 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",outputLocation.Data()));
3572 AliInfo(Form("%s will be stored in %s",
3573 outputFile.IsNull() ? "AliESDs.root" : outputFile.Data(),
3574 outputLocation.Data()));