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"
192 #include "AliTriggerRunScalers.h"
193 #include "AliCTPTimeParams.h"
194 #include "AliESDHLTDecision.h"
196 ClassImp(AliReconstruction)
198 //_____________________________________________________________________________
199 const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
201 //_____________________________________________________________________________
202 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
204 fRunVertexFinder(kTRUE),
205 fRunVertexFinderTracks(kTRUE),
206 fRunHLTTracking(kFALSE),
207 fRunMuonTracking(kFALSE),
209 fRunCascadeFinder(kTRUE),
211 fWriteAlignmentData(kFALSE),
212 fWriteESDfriend(kFALSE),
213 fFillTriggerESD(kTRUE),
221 fRunLocalReconstruction("ALL"),
225 fUseTrackingErrorsForAlignment(""),
226 fGAliceFileName(gAliceFilename),
229 fProofOutputFileName(""),
230 fProofOutputLocation(""),
231 fProofOutputDataset(kFALSE),
232 fProofOutputArchive(""),
236 fNumberOfEventsPerFile((UInt_t)-1),
238 fLoadAlignFromCDB(kTRUE),
239 fLoadAlignData("ALL"),
244 fCTPTimeParams(NULL),
248 fParentRawReader(NULL),
252 fSPDTrackleter(NULL),
254 fDiamondProfileSPD(NULL),
255 fDiamondProfile(NULL),
256 fDiamondProfileTPC(NULL),
257 fListOfCosmicTriggers(NULL),
261 fAlignObjArray(NULL),
265 fInitCDBCalled(kFALSE),
266 fSetRunNumberFromDataCalled(kFALSE),
271 fSameQACycle(kFALSE),
272 fInitQACalled(kFALSE),
273 fWriteQAExpertData(kTRUE),
274 fRunPlaneEff(kFALSE),
283 fIsNewRunLoader(kFALSE),
287 // create reconstruction object with default parameters
290 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
291 fReconstructor[iDet] = NULL;
292 fLoader[iDet] = NULL;
293 fTracker[iDet] = NULL;
295 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
296 fQACycles[iDet] = 999999 ;
297 fQAWriteExpert[iDet] = kFALSE ;
303 //_____________________________________________________________________________
304 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
306 fRunVertexFinder(rec.fRunVertexFinder),
307 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
308 fRunHLTTracking(rec.fRunHLTTracking),
309 fRunMuonTracking(rec.fRunMuonTracking),
310 fRunV0Finder(rec.fRunV0Finder),
311 fRunCascadeFinder(rec.fRunCascadeFinder),
312 fStopOnError(rec.fStopOnError),
313 fWriteAlignmentData(rec.fWriteAlignmentData),
314 fWriteESDfriend(rec.fWriteESDfriend),
315 fFillTriggerESD(rec.fFillTriggerESD),
317 fCleanESD(rec.fCleanESD),
318 fV0DCAmax(rec.fV0DCAmax),
319 fV0CsPmin(rec.fV0CsPmin),
323 fRunLocalReconstruction(rec.fRunLocalReconstruction),
324 fRunTracking(rec.fRunTracking),
325 fFillESD(rec.fFillESD),
326 fLoadCDB(rec.fLoadCDB),
327 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
328 fGAliceFileName(rec.fGAliceFileName),
329 fRawInput(rec.fRawInput),
330 fESDOutput(rec.fESDOutput),
331 fProofOutputFileName(rec.fProofOutputFileName),
332 fProofOutputLocation(rec.fProofOutputLocation),
333 fProofOutputDataset(rec.fProofOutputDataset),
334 fProofOutputArchive(rec.fProofOutputArchive),
335 fEquipIdMap(rec.fEquipIdMap),
336 fFirstEvent(rec.fFirstEvent),
337 fLastEvent(rec.fLastEvent),
338 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
340 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
341 fLoadAlignData(rec.fLoadAlignData),
342 fUseHLTData(rec.fUseHLTData),
346 fCTPTimeParams(NULL),
350 fParentRawReader(NULL),
352 fRecoParam(rec.fRecoParam),
354 fSPDTrackleter(NULL),
356 fDiamondProfileSPD(rec.fDiamondProfileSPD),
357 fDiamondProfile(rec.fDiamondProfile),
358 fDiamondProfileTPC(rec.fDiamondProfileTPC),
359 fListOfCosmicTriggers(NULL),
363 fAlignObjArray(rec.fAlignObjArray),
364 fCDBUri(rec.fCDBUri),
365 fQARefUri(rec.fQARefUri),
367 fInitCDBCalled(rec.fInitCDBCalled),
368 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
369 fQADetectors(rec.fQADetectors),
370 fQATasks(rec.fQATasks),
372 fRunGlobalQA(rec.fRunGlobalQA),
373 fSameQACycle(rec.fSameQACycle),
374 fInitQACalled(rec.fInitQACalled),
375 fWriteQAExpertData(rec.fWriteQAExpertData),
376 fRunPlaneEff(rec.fRunPlaneEff),
385 fIsNewRunLoader(rec.fIsNewRunLoader),
391 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
392 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
394 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
395 fReconstructor[iDet] = NULL;
396 fLoader[iDet] = NULL;
397 fTracker[iDet] = NULL;
400 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
401 fQACycles[iDet] = rec.fQACycles[iDet];
402 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
405 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
406 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
410 //_____________________________________________________________________________
411 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
413 // assignment operator
414 // Used in PROOF mode
415 // Be very careful while modifing it!
416 // Simple rules to follow:
417 // for persistent data members - use their assignment operators
418 // for non-persistent ones - do nothing or take the default values from constructor
419 // TSelector members should not be touched
420 if(&rec == this) return *this;
422 fRunVertexFinder = rec.fRunVertexFinder;
423 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
424 fRunHLTTracking = rec.fRunHLTTracking;
425 fRunMuonTracking = rec.fRunMuonTracking;
426 fRunV0Finder = rec.fRunV0Finder;
427 fRunCascadeFinder = rec.fRunCascadeFinder;
428 fStopOnError = rec.fStopOnError;
429 fWriteAlignmentData = rec.fWriteAlignmentData;
430 fWriteESDfriend = rec.fWriteESDfriend;
431 fFillTriggerESD = rec.fFillTriggerESD;
433 fCleanESD = rec.fCleanESD;
434 fV0DCAmax = rec.fV0DCAmax;
435 fV0CsPmin = rec.fV0CsPmin;
439 fRunLocalReconstruction = rec.fRunLocalReconstruction;
440 fRunTracking = rec.fRunTracking;
441 fFillESD = rec.fFillESD;
442 fLoadCDB = rec.fLoadCDB;
443 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
444 fGAliceFileName = rec.fGAliceFileName;
445 fRawInput = rec.fRawInput;
446 fESDOutput = rec.fESDOutput;
447 fProofOutputFileName = rec.fProofOutputFileName;
448 fProofOutputLocation = rec.fProofOutputLocation;
449 fProofOutputDataset = rec.fProofOutputDataset;
450 fProofOutputArchive = rec.fProofOutputArchive;
451 fEquipIdMap = rec.fEquipIdMap;
452 fFirstEvent = rec.fFirstEvent;
453 fLastEvent = rec.fLastEvent;
454 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
456 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
457 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
460 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
461 fLoadAlignData = rec.fLoadAlignData;
462 fUseHLTData = rec.fUseHLTData;
464 delete fRunInfo; fRunInfo = NULL;
465 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
467 fEventInfo = rec.fEventInfo;
469 delete fRunScalers; fRunScalers = NULL;
470 if (rec.fRunScalers) fRunScalers = new AliTriggerRunScalers(*rec.fRunScalers);
472 delete fCTPTimeParams; fCTPTimeParams = NULL;
473 if (rec.fCTPTimeParams) fCTPTimeParams = new AliCTPTimeParams(*rec.fCTPTimeParams);
477 fParentRawReader = NULL;
479 fRecoParam = rec.fRecoParam;
481 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
482 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
483 delete fLoader[iDet]; fLoader[iDet] = NULL;
484 delete fTracker[iDet]; fTracker[iDet] = NULL;
487 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
488 fQACycles[iDet] = rec.fQACycles[iDet];
489 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
492 delete fSPDTrackleter; fSPDTrackleter = NULL;
494 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
495 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
496 delete fDiamondProfile; fDiamondProfile = NULL;
497 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
498 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
499 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
501 delete fListOfCosmicTriggers; fListOfCosmicTriggers = NULL;
502 if (rec.fListOfCosmicTriggers) fListOfCosmicTriggers = (THashTable*)((rec.fListOfCosmicTriggers)->Clone());
504 delete fGRPData; fGRPData = NULL;
505 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
506 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
508 delete fAlignObjArray; fAlignObjArray = NULL;
511 fQARefUri = rec.fQARefUri;
512 fSpecCDBUri.Delete();
513 fInitCDBCalled = rec.fInitCDBCalled;
514 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
515 fQADetectors = rec.fQADetectors;
516 fQATasks = rec.fQATasks;
518 fRunGlobalQA = rec.fRunGlobalQA;
519 fSameQACycle = rec.fSameQACycle;
520 fInitQACalled = rec.fInitQACalled;
521 fWriteQAExpertData = rec.fWriteQAExpertData;
522 fRunPlaneEff = rec.fRunPlaneEff;
531 fIsNewRunLoader = rec.fIsNewRunLoader;
538 //_____________________________________________________________________________
539 AliReconstruction::~AliReconstruction()
544 if (fListOfCosmicTriggers) {
545 fListOfCosmicTriggers->Delete();
546 delete fListOfCosmicTriggers;
550 delete fCTPTimeParams;
552 if (fAlignObjArray) {
553 fAlignObjArray->Delete();
554 delete fAlignObjArray;
556 fSpecCDBUri.Delete();
558 AliCodeTimer::Instance()->Print();
561 //_____________________________________________________________________________
562 void AliReconstruction::InitQA()
564 //Initialize the QA and start of cycle
565 AliCodeTimerAuto("",0);
567 if (fInitQACalled) return;
568 fInitQACalled = kTRUE;
570 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kRECMODE) ;
571 if (fWriteQAExpertData)
572 qam->SetWriteExpert() ;
574 if (qam->IsDefaultStorageSet()) {
575 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
576 AliWarning("Default QA reference storage has been already set !");
577 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fQARefUri.Data()));
578 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
579 fQARefUri = qam->GetDefaultStorage()->GetURI();
581 if (fQARefUri.Length() > 0) {
582 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
583 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
584 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
586 fQARefUri="local://$ALICE_ROOT/QAref";
587 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
588 AliWarning("Default QA refeference storage not yet set !!!!");
589 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
590 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
593 qam->SetDefaultStorage(fQARefUri);
597 qam->SetActiveDetectors(fQADetectors) ;
598 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) {
599 qam->SetCycleLength(AliQAv1::DETECTORINDEX_t(det), fQACycles[det]) ;
600 qam->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
602 if (!fRawReader && !fInput && IsInTasks(AliQAv1::kRAWS))
603 fQATasks.ReplaceAll(Form("%d",AliQAv1::kRAWS), "") ;
604 qam->SetTasks(fQATasks) ;
605 qam->InitQADataMaker(AliCDBManager::Instance()->GetRun()) ;
608 Bool_t sameCycle = kFALSE ;
609 AliQADataMaker *qadm = qam->GetQADataMaker(AliQAv1::kGLOBAL);
610 AliInfo(Form("Initializing the global QA data maker"));
611 if (IsInTasks(AliQAv1::kRECPOINTS)) {
612 qadm->StartOfCycle(AliQAv1::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
613 TObjArray **arr=qadm->Init(AliQAv1::kRECPOINTS);
614 AliTracker::SetResidualsArray(arr);
617 if (IsInTasks(AliQAv1::kESDS)) {
618 qadm->StartOfCycle(AliQAv1::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
619 qadm->Init(AliQAv1::kESDS);
622 AliSysInfo::AddStamp("InitQA") ;
625 //_____________________________________________________________________________
626 void AliReconstruction::MergeQA(const char *fileName)
628 //Initialize the QA and start of cycle
629 AliCodeTimerAuto("",0) ;
630 AliQAManager::QAManager()->Merge(AliCDBManager::Instance()->GetRun(),fileName) ;
631 AliSysInfo::AddStamp("MergeQA") ;
634 //_____________________________________________________________________________
635 void AliReconstruction::InitCDB()
637 // activate a default CDB storage
638 // First check if we have any CDB storage set, because it is used
639 // to retrieve the calibration and alignment constants
640 AliCodeTimerAuto("",0);
642 if (fInitCDBCalled) return;
643 fInitCDBCalled = kTRUE;
645 AliCDBManager* man = AliCDBManager::Instance();
646 if (man->IsDefaultStorageSet())
648 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
649 AliWarning("Default CDB storage has been already set !");
650 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
651 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
652 fCDBUri = man->GetDefaultStorage()->GetURI();
655 if (fCDBUri.Length() > 0)
657 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
658 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
659 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
660 man->SetDefaultStorage(fCDBUri);
662 else if (!man->GetRaw()){
663 fCDBUri="local://$ALICE_ROOT/OCDB";
664 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
665 AliWarning("Default CDB storage not yet set !!!!");
666 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
667 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
668 man->SetDefaultStorage(fCDBUri);
671 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
672 AliWarning("Default storage will be set after setting the Run Number!!!");
673 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
677 // Now activate the detector specific CDB storage locations
678 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
679 TObject* obj = fSpecCDBUri[i];
681 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
682 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
683 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
684 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
686 AliSysInfo::AddStamp("InitCDB");
689 //_____________________________________________________________________________
690 void AliReconstruction::SetDefaultStorage(const char* uri) {
691 // Store the desired default CDB storage location
692 // Activate it later within the Run() method
698 //_____________________________________________________________________________
699 void AliReconstruction::SetQARefDefaultStorage(const char* uri) {
700 // Store the desired default CDB storage location
701 // Activate it later within the Run() method
704 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
707 //_____________________________________________________________________________
708 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
709 // Store a detector-specific CDB storage location
710 // Activate it later within the Run() method
712 AliCDBPath aPath(calibType);
713 if(!aPath.IsValid()){
714 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
715 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
716 if(!strcmp(calibType, fgkDetectorName[iDet])) {
717 aPath.SetPath(Form("%s/*", calibType));
718 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
722 if(!aPath.IsValid()){
723 AliError(Form("Not a valid path or detector: %s", calibType));
728 // // check that calibType refers to a "valid" detector name
729 // Bool_t isDetector = kFALSE;
730 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
731 // TString detName = fgkDetectorName[iDet];
732 // if(aPath.GetLevel0() == detName) {
733 // isDetector = kTRUE;
739 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
743 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
744 if (obj) fSpecCDBUri.Remove(obj);
745 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
749 //_____________________________________________________________________________
750 Bool_t AliReconstruction::SetRunNumberFromData()
752 // The method is called in Run() in order
753 // to set a correct run number.
754 // In case of raw data reconstruction the
755 // run number is taken from the raw data header
757 if (fSetRunNumberFromDataCalled) return kTRUE;
758 fSetRunNumberFromDataCalled = kTRUE;
760 AliCDBManager* man = AliCDBManager::Instance();
763 if(fRawReader->NextEvent()) {
764 if(man->GetRun() > 0) {
765 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
767 man->SetRun(fRawReader->GetRunNumber());
768 fRawReader->RewindEvents();
771 if(man->GetRun() > 0) {
772 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
775 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
781 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
783 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
788 // read run number from gAlice
789 if(rl->GetHeader()) {
790 man->SetRun(rl->GetHeader()->GetRun());
795 AliError("Neither run-loader header nor RawReader objects are found !");
807 //_____________________________________________________________________________
808 void AliReconstruction::SetCDBLock() {
809 // Set CDB lock: from now on it is forbidden to reset the run number
810 // or the default storage or to activate any further storage!
812 AliCDBManager::Instance()->SetLock(1);
815 //_____________________________________________________________________________
816 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
818 // Read the alignment objects from CDB.
819 // Each detector is supposed to have the
820 // alignment objects in DET/Align/Data CDB path.
821 // All the detector objects are then collected,
822 // sorted by geometry level (starting from ALIC) and
823 // then applied to the TGeo geometry.
824 // Finally an overlaps check is performed.
826 // Load alignment data from CDB and fill fAlignObjArray
827 if(fLoadAlignFromCDB){
829 TString detStr = detectors;
830 TString loadAlObjsListOfDets = "";
832 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
833 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
834 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
836 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
838 loadAlObjsListOfDets += fgkDetectorName[iDet];
839 loadAlObjsListOfDets += " ";
841 } // end loop over detectors
843 if(AliGeomManager::GetNalignable("GRP") != 0)
844 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
845 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
846 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
848 // Check if the array with alignment objects was
849 // provided by the user. If yes, apply the objects
850 // to the present TGeo geometry
851 if (fAlignObjArray) {
852 if (gGeoManager && gGeoManager->IsClosed()) {
853 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
854 AliError("The misalignment of one or more volumes failed!"
855 "Compare the list of simulated detectors and the list of detector alignment data!");
860 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
866 if (fAlignObjArray) {
867 fAlignObjArray->Delete();
868 delete fAlignObjArray; fAlignObjArray=NULL;
874 //_____________________________________________________________________________
875 void AliReconstruction::SetGAliceFile(const char* fileName)
877 // set the name of the galice file
879 fGAliceFileName = fileName;
882 //_____________________________________________________________________________
883 void AliReconstruction::SetInput(const char* input)
885 // In case the input string starts with 'mem://', we run in an online mode
886 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
887 // file is assumed. One can give as an input:
888 // mem://: - events taken from DAQ monitoring libs online
890 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
891 if (input) fRawInput = input;
894 //_____________________________________________________________________________
895 void AliReconstruction::SetOutput(const char* output)
897 // Set the output ESD filename
898 // 'output' is a normalt ROOT url
899 // The method is used in case of raw-data reco with PROOF
900 if (output) fESDOutput = output;
903 //_____________________________________________________________________________
904 void AliReconstruction::SetOption(const char* detector, const char* option)
906 // set options for the reconstruction of a detector
908 TObject* obj = fOptions.FindObject(detector);
909 if (obj) fOptions.Remove(obj);
910 fOptions.Add(new TNamed(detector, option));
913 //_____________________________________________________________________________
914 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
916 // Set custom reconstruction parameters for a given detector
917 // Single set of parameters for all the events
919 // First check if the reco-params are global
920 if(!strcmp(detector, "GRP")) {
922 fRecoParam.AddDetRecoParam(kNDetectors,par);
926 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
927 if(!strcmp(detector, fgkDetectorName[iDet])) {
929 fRecoParam.AddDetRecoParam(iDet,par);
936 //_____________________________________________________________________________
937 Bool_t AliReconstruction::InitGRP() {
938 //------------------------------------
939 // Initialization of the GRP entry
940 //------------------------------------
941 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
945 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
948 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
950 fGRPData = new AliGRPObject();
951 fGRPData->ReadValuesFromMap(m);
955 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
956 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
960 // FIX ME: The unloading of GRP entry is temporarily disabled
961 // because ZDC and VZERO are using it in order to initialize
962 // their reconstructor objects. In the future one has to think
963 // of propagating AliRunInfo to the reconstructors.
964 // AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
968 AliError("No GRP entry found in OCDB!");
972 TString lhcState = fGRPData->GetLHCState();
973 if (lhcState==AliGRPObject::GetInvalidString()) {
974 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
975 lhcState = "UNKNOWN";
978 TString beamType = fGRPData->GetBeamType();
979 if (beamType==AliGRPObject::GetInvalidString()) {
980 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
981 beamType = "UNKNOWN";
984 Float_t beamEnergy = fGRPData->GetBeamEnergy();
985 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
986 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
990 TString runType = fGRPData->GetRunType();
991 if (runType==AliGRPObject::GetInvalidString()) {
992 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
996 Int_t activeDetectors = fGRPData->GetDetectorMask();
997 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
998 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
999 activeDetectors = 1074790399;
1002 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1006 // Process the list of active detectors
1007 if (activeDetectors) {
1008 UInt_t detMask = activeDetectors;
1009 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1010 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1011 fFillESD = MatchDetectorList(fFillESD,detMask);
1012 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1013 fLoadCDB.Form("%s %s %s %s",
1014 fRunLocalReconstruction.Data(),
1015 fRunTracking.Data(),
1017 fQADetectors.Data());
1018 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1019 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1) &&
1020 !((detMask >> AliDAQ::DetectorID("ITSSDD")) & 0x1) &&
1021 !((detMask >> AliDAQ::DetectorID("ITSSSD")) & 0x1) ) {
1022 // switch off the vertexer
1023 AliInfo("SPD,SDD,SSD is not in the list of active detectors. Vertexer switched off.");
1024 fRunVertexFinder = kFALSE;
1026 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1027 // switch off the reading of CTP raw-data payload
1028 if (fFillTriggerESD) {
1029 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1030 fFillTriggerESD = kFALSE;
1035 AliInfo("===================================================================================");
1036 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1037 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1038 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1039 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1040 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1041 AliInfo("===================================================================================");
1043 //*** Dealing with the magnetic field map
1044 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {
1045 if (TGeoGlobalMagField::Instance()->GetField()->TestBit(AliMagF::kOverrideGRP)) {
1046 AliInfo("ExpertMode!!! GRP information will be ignored !");
1047 AliInfo("ExpertMode!!! Running with the externally locked B field !");
1050 AliInfo("Destroying existing B field instance!");
1051 delete TGeoGlobalMagField::Instance();
1054 if ( !TGeoGlobalMagField::Instance()->IsLocked() ) {
1055 // Construct the field map out of the information retrieved from GRP.
1058 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1059 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1060 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1064 Char_t l3Polarity = fGRPData->GetL3Polarity();
1065 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1066 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1071 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1072 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1073 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1077 Char_t diPolarity = fGRPData->GetDipolePolarity();
1078 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1079 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1083 // read special bits for the polarity convention and map type
1084 Int_t polConvention = fGRPData->IsPolarityConventionLHC() ? AliMagF::kConvLHC : AliMagF::kConvDCS2008;
1085 Bool_t uniformB = fGRPData->IsUniformBMap();
1088 AliMagF* fld = AliMagF::CreateFieldMap(TMath::Abs(l3Current) * (l3Polarity ? -1:1),
1089 TMath::Abs(diCurrent) * (diPolarity ? -1:1),
1090 polConvention,uniformB,beamEnergy, beamType.Data());
1092 TGeoGlobalMagField::Instance()->SetField( fld );
1093 TGeoGlobalMagField::Instance()->Lock();
1094 AliInfo("Running with the B field constructed out of GRP !");
1096 else AliFatal("Failed to create a B field map !");
1098 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1101 //*** Get the diamond profiles from OCDB
1102 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1104 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1106 AliError("No SPD diamond profile found in OCDB!");
1109 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1111 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1113 AliError("No diamond profile found in OCDB!");
1116 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1118 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1120 AliError("No TPC diamond profile found in OCDB!");
1123 entry = AliCDBManager::Instance()->Get("GRP/Calib/CosmicTriggers");
1125 fListOfCosmicTriggers = dynamic_cast<THashTable*>(entry->GetObject());
1127 AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
1130 if (!fListOfCosmicTriggers) {
1131 AliWarning("Can not get list of cosmic triggers from OCDB! Cosmic event specie will be effectively disabled!");
1137 //_____________________________________________________________________________
1138 Bool_t AliReconstruction::LoadCDB()
1140 AliCodeTimerAuto("",0);
1142 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1144 TString detStr = fLoadCDB;
1145 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1146 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1147 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1150 // Temporary fix - one has to define the correct policy in order
1151 // to load the trigger OCDB entries only for the detectors that
1152 // in the trigger or that are needed in order to put correct
1153 // information in ESD
1154 AliCDBManager::Instance()->GetAll("TRIGGER/*/*");
1158 //_____________________________________________________________________________
1159 Bool_t AliReconstruction::LoadTriggerScalersCDB()
1161 AliCodeTimerAuto("",0);
1163 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
1167 AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
1168 fRunScalers = dynamic_cast<AliTriggerRunScalers*> (entry->GetObject());
1170 if (fRunScalers->CorrectScalersOverflow() == 0) AliInfo("32bit Trigger counters corrected for overflow");
1175 //_____________________________________________________________________________
1176 Bool_t AliReconstruction::LoadCTPTimeParamsCDB()
1178 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
1182 AliInfo("Found an AliCTPTimeParams in GRP/CTP/CTPtiming, reading it");
1183 fCTPTimeParams = dynamic_cast<AliCTPTimeParams*> (entry->GetObject());
1190 //_____________________________________________________________________________
1191 Bool_t AliReconstruction::Run(const char* input)
1194 AliCodeTimerAuto("",0);
1197 if (GetAbort() != TSelector::kContinue) return kFALSE;
1199 TChain *chain = NULL;
1200 if (fRawReader && (chain = fRawReader->GetChain())) {
1201 Long64_t nEntries = (fLastEvent < 0) ? (TChain::kBigNumber) : (fLastEvent - fFirstEvent + 1);
1204 // Temporary fix for long raw-data runs (until socket timeout handling in PROOF is revised)
1205 gProof->Exec("gEnv->SetValue(\"Proof.SocketActivityTimeout\",-1)", kTRUE);
1208 gProof->Exec("TGrid::Connect(\"alien://\")",kTRUE);
1210 TMessage::EnableSchemaEvolutionForAll(kTRUE);
1211 gProof->Exec("TMessage::EnableSchemaEvolutionForAll(kTRUE)",kTRUE);
1213 gProof->AddInput(this);
1215 if (!ParseOutput()) return kFALSE;
1217 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1219 chain->Process("AliReconstruction","",nEntries,fFirstEvent);
1222 chain->Process(this,"",nEntries,fFirstEvent);
1227 if (GetAbort() != TSelector::kContinue) return kFALSE;
1229 if (GetAbort() != TSelector::kContinue) return kFALSE;
1230 //******* The loop over events
1231 AliInfo("Starting looping over events");
1233 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1234 (fRawReader && fRawReader->NextEvent())) {
1235 if (!ProcessEvent(iEvent)) {
1236 Abort("ProcessEvent",TSelector::kAbortFile);
1242 if (GetAbort() != TSelector::kContinue) return kFALSE;
1244 if (GetAbort() != TSelector::kContinue) return kFALSE;
1250 //_____________________________________________________________________________
1251 void AliReconstruction::InitRawReader(const char* input)
1253 AliCodeTimerAuto("",0);
1255 // Init raw-reader and
1256 // set the input in case of raw data
1257 if (input) fRawInput = input;
1258 fRawReader = AliRawReader::Create(fRawInput.Data());
1260 if (fRawInput.IsNull()) {
1261 AliInfo("Reconstruction will run over digits");
1264 AliFatal("Can not create raw-data reader ! Exiting...");
1268 if (!fEquipIdMap.IsNull() && fRawReader)
1269 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1271 if (!fUseHLTData.IsNull()) {
1272 // create the RawReaderHLT which performs redirection of HLT input data for
1273 // the specified detectors
1274 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1276 fParentRawReader=fRawReader;
1277 fRawReader=pRawReader;
1279 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1282 AliSysInfo::AddStamp("CreateRawReader");
1285 //_____________________________________________________________________________
1286 void AliReconstruction::InitRun(const char* input)
1288 // Initialization of raw-reader,
1289 // run number, CDB etc.
1290 AliCodeTimerAuto("",0);
1291 AliSysInfo::AddStamp("Start");
1293 // Initialize raw-reader if any
1294 InitRawReader(input);
1296 // Initialize the CDB storage
1299 // Set run number in CDBManager (if it is not already set by the user)
1300 if (!SetRunNumberFromData()) {
1301 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1305 // Set CDB lock: from now on it is forbidden to reset the run number
1306 // or the default storage or to activate any further storage!
1311 //_____________________________________________________________________________
1312 void AliReconstruction::Begin(TTree *)
1314 // Initialize AlReconstruction before
1315 // going into the event loop
1316 // Should follow the TSelector convention
1317 // i.e. initialize only the object on the client side
1318 AliCodeTimerAuto("",0);
1320 AliReconstruction *reco = NULL;
1322 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1325 AliSysInfo::AddStamp("ReadInputInBegin");
1328 // Import ideal TGeo geometry and apply misalignment
1330 TString geom(gSystem->DirName(fGAliceFileName));
1331 geom += "/geometry.root";
1332 AliGeomManager::LoadGeometry(geom.Data());
1334 Abort("LoadGeometry", TSelector::kAbortProcess);
1337 AliSysInfo::AddStamp("LoadGeom");
1338 TString detsToCheck=fRunLocalReconstruction;
1339 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1340 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1343 AliSysInfo::AddStamp("CheckGeom");
1346 if (!MisalignGeometry(fLoadAlignData)) {
1347 Abort("MisalignGeometry", TSelector::kAbortProcess);
1350 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1351 AliSysInfo::AddStamp("MisalignGeom");
1354 Abort("InitGRP", TSelector::kAbortProcess);
1357 AliSysInfo::AddStamp("InitGRP");
1360 Abort("LoadCDB", TSelector::kAbortProcess);
1363 AliSysInfo::AddStamp("LoadCDB");
1365 if (!LoadTriggerScalersCDB()) {
1366 Abort("LoadTriggerScalersCDB", TSelector::kAbortProcess);
1369 AliSysInfo::AddStamp("LoadTriggerScalersCDB");
1371 if (!LoadCTPTimeParamsCDB()) {
1372 Abort("LoadCTPTimeParamsCDB", TSelector::kAbortProcess);
1375 AliSysInfo::AddStamp("LoadCTPTimeParamsCDB");
1377 // Read the reconstruction parameters from OCDB
1378 if (!InitRecoParams()) {
1379 AliWarning("Not all detectors have correct RecoParam objects initialized");
1381 AliSysInfo::AddStamp("InitRecoParams");
1383 if (fInput && gProof) {
1384 if (reco) *reco = *this;
1386 gGeoManager->SetName("Geometry");
1387 gProof->AddInputData(gGeoManager,kTRUE);
1389 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1390 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1391 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1392 magFieldMap->SetName("MagneticFieldMap");
1393 gProof->AddInputData(magFieldMap,kTRUE);
1398 //_____________________________________________________________________________
1399 void AliReconstruction::SlaveBegin(TTree*)
1401 // Initialization related to run-loader,
1402 // vertexer, trackers, recontructors
1403 // In proof mode it is executed on the slave
1404 AliCodeTimerAuto("",0);
1406 TProofOutputFile *outProofFile = NULL;
1408 if (AliDebugLevel() > 0) fInput->Print();
1409 if (AliDebugLevel() > 10) fInput->Dump();
1410 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1413 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1415 AliGeomManager::SetGeometry(tgeo);
1417 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1418 Int_t runNumber = -1;
1419 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1420 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1421 man->SetCacheFlag(kTRUE);
1422 man->SetLock(kTRUE);
1426 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1427 AliMagF *newMap = new AliMagF(*map);
1428 if (!newMap->LoadParameterization()) {
1429 Abort("AliMagF::LoadParameterization", TSelector::kAbortProcess);
1432 TGeoGlobalMagField::Instance()->SetField(newMap);
1433 TGeoGlobalMagField::Instance()->Lock();
1435 if (TNamed *outputFileName = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE"))
1436 fProofOutputFileName = outputFileName->GetTitle();
1437 if (TNamed *outputLocation = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_LOCATION"))
1438 fProofOutputLocation = outputLocation->GetTitle();
1439 if (fInput->FindObject("PROOF_OUTPUTFILE_DATASET"))
1440 fProofOutputDataset = kTRUE;
1441 if (TNamed *archiveList = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_ARCHIVE"))
1442 fProofOutputArchive = archiveList->GetTitle();
1443 if (!fProofOutputFileName.IsNull() &&
1444 !fProofOutputLocation.IsNull() &&
1445 fProofOutputArchive.IsNull()) {
1446 if (!fProofOutputDataset) {
1447 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"M");
1448 outProofFile->SetOutputFileName(Form("%s%s",fProofOutputLocation.Data(),fProofOutputFileName.Data()));
1451 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"DROV",fProofOutputLocation.Data());
1453 if (AliDebugLevel() > 0) outProofFile->Dump();
1454 fOutput->Add(outProofFile);
1456 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1459 // get the run loader
1460 if (!InitRunLoader()) {
1461 Abort("InitRunLoader", TSelector::kAbortProcess);
1464 AliSysInfo::AddStamp("LoadLoader");
1466 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1469 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1470 Abort("CreateTrackers", TSelector::kAbortProcess);
1473 AliSysInfo::AddStamp("CreateTrackers");
1475 // create the ESD output file and tree
1476 if (!outProofFile) {
1477 ffile = TFile::Open("AliESDs.root", "RECREATE");
1478 ffile->SetCompressionLevel(2);
1479 if (!ffile->IsOpen()) {
1480 Abort("OpenESDFile", TSelector::kAbortProcess);
1485 AliInfo(Form("Opening output PROOF file: %s/%s",
1486 outProofFile->GetDir(), outProofFile->GetFileName()));
1487 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1488 Abort(Form("Problems opening output PROOF file: %s/%s",
1489 outProofFile->GetDir(), outProofFile->GetFileName()),
1490 TSelector::kAbortProcess);
1495 ftree = new TTree("esdTree", "Tree with ESD objects");
1496 fesd = new AliESDEvent();
1497 fesd->CreateStdContent();
1498 if (fesd->GetESDRun()) ((AliESDRun*)fesd->GetESDRun())->SetBeamEnergyIsSqrtSHalfGeV();
1499 // add a so far non-std object to the ESD, this will
1500 // become part of the std content
1501 fesd->AddObject(new AliESDHLTDecision);
1503 fesd->WriteToTree(ftree);
1504 if (fWriteESDfriend) {
1506 // Since we add the branch manually we must
1507 // book and add it after WriteToTree
1508 // otherwise it is created twice,
1509 // once via writetotree and once here.
1510 // The case for AliESDfriend is now
1511 // caught also in AlIESDEvent::WriteToTree but
1512 // be careful when changing the name (AliESDfriend is not
1513 // a TNamed so we had to hardwire it)
1514 fesdf = new AliESDfriend();
1515 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1516 br->SetFile("AliESDfriends.root");
1517 fesd->AddObject(fesdf);
1519 ftree->GetUserInfo()->Add(fesd);
1521 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1522 fhltesd = new AliESDEvent();
1523 fhltesd->CreateStdContent();
1524 if (fhltesd->GetESDRun()) ((AliESDRun*)fhltesd->GetESDRun())->SetBeamEnergyIsSqrtSHalfGeV();
1525 // read the ESD template from CDB
1526 // HLT is allowed to put non-std content to its ESD, the non-std
1527 // objects need to be created before invocation of WriteToTree in
1528 // order to create all branches. Initialization is done from an
1529 // ESD layout template in CDB
1530 AliCDBManager* man = AliCDBManager::Instance();
1531 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1532 AliCDBEntry* hltESDConfig=NULL;
1533 if (man->GetId(hltESDConfigPath)!=NULL &&
1534 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1535 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1537 // init all internal variables from the list of objects
1538 pESDLayout->GetStdContent();
1540 // copy content and create non-std objects
1541 *fhltesd=*pESDLayout;
1544 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1545 hltESDConfigPath.GetPath().Data()));
1549 fhltesd->WriteToTree(fhlttree);
1550 fhlttree->GetUserInfo()->Add(fhltesd);
1552 ProcInfo_t procInfo;
1553 gSystem->GetProcInfo(&procInfo);
1554 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1557 //Initialize the QA and start of cycle
1558 if (fRunQA || fRunGlobalQA)
1561 //Initialize the Plane Efficiency framework
1562 if (fRunPlaneEff && !InitPlaneEff()) {
1563 Abort("InitPlaneEff", TSelector::kAbortProcess);
1567 if (strcmp(gProgName,"alieve") == 0)
1568 fRunAliEVE = InitAliEVE();
1573 //_____________________________________________________________________________
1574 Bool_t AliReconstruction::Process(Long64_t entry)
1576 // run the reconstruction over a single entry
1577 // from the chain with raw data
1578 AliCodeTimerAuto("",0);
1580 TTree *currTree = fChain->GetTree();
1581 AliRawVEvent *event = NULL;
1582 currTree->SetBranchAddress("rawevent",&event);
1583 currTree->GetEntry(entry);
1584 fRawReader = new AliRawReaderRoot(event);
1585 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1593 //_____________________________________________________________________________
1594 void AliReconstruction::Init(TTree *tree)
1597 AliError("The input tree is not found!");
1603 //_____________________________________________________________________________
1604 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1606 // run the reconstruction over a single event
1607 // The event loop is steered in Run method
1609 AliCodeTimerAuto("",0);
1611 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1612 fRunLoader->SetEventNumber(iEvent);
1613 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1615 fRunLoader->TreeE()->Fill();
1616 if (fRawReader && fRawReader->UseAutoSaveESD())
1617 fRunLoader->TreeE()->AutoSave("SaveSelf");
1620 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1624 AliInfo(Form("processing event %d", iEvent));
1626 fRunLoader->GetEvent(iEvent);
1628 // Fill Event-info object
1630 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo,fListOfCosmicTriggers);
1631 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1633 // Set the reco-params
1635 TString detStr = fLoadCDB;
1636 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1637 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1638 AliReconstructor *reconstructor = GetReconstructor(iDet);
1639 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1640 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1641 reconstructor->SetRecoParam(par);
1642 reconstructor->SetEventInfo(&fEventInfo);
1644 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1645 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(par->GetEventSpecie())) ;
1649 const AliDetectorRecoParam *grppar = fRecoParam.GetDetRecoParam(kNDetectors);
1650 AliQAManager::QAManager()->SetRecoParam(AliQAv1::kGLOBAL, grppar) ;
1651 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(grppar->GetEventSpecie())) ;
1655 if (fRunQA && IsInTasks(AliQAv1::kRAWS)) {
1656 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1657 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1659 // local single event reconstruction
1660 if (!fRunLocalReconstruction.IsNull()) {
1661 TString detectors=fRunLocalReconstruction;
1662 // run HLT event reconstruction first
1663 // ;-( IsSelected changes the string
1664 if (IsSelected("HLT", detectors) &&
1665 !RunLocalEventReconstruction("HLT")) {
1666 if (fStopOnError) {CleanUp(); return kFALSE;}
1668 detectors=fRunLocalReconstruction;
1669 detectors.ReplaceAll("HLT", "");
1670 if (!RunLocalEventReconstruction(detectors)) {
1679 // fill Event header information from the RawEventHeader
1680 if (fRawReader){FillRawEventHeaderESD(fesd);}
1681 if (fRawReader){FillRawEventHeaderESD(fhltesd);}
1683 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1684 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1685 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1686 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1688 fesd->SetEventSpecie(fRecoParam.GetEventSpecie());
1689 fhltesd->SetEventSpecie(fRecoParam.GetEventSpecie());
1691 // Set magnetic field from the tracker
1692 fesd->SetMagneticField(AliTracker::GetBz());
1693 fhltesd->SetMagneticField(AliTracker::GetBz());
1695 AliMagF* fld = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1696 if (fld) { // set info needed for field initialization
1697 fesd->SetCurrentL3(fld->GetCurrentSol());
1698 fesd->SetCurrentDip(fld->GetCurrentDip());
1699 fesd->SetBeamEnergy(fld->GetBeamEnergy());
1700 fesd->SetBeamType(fld->GetBeamTypeText());
1701 fesd->SetUniformBMap(fld->IsUniform());
1702 fesd->SetBInfoStored();
1704 fhltesd->SetCurrentL3(fld->GetCurrentSol());
1705 fhltesd->SetCurrentDip(fld->GetCurrentDip());
1706 fhltesd->SetBeamEnergy(fld->GetBeamEnergy());
1707 fhltesd->SetBeamType(fld->GetBeamTypeText());
1708 fhltesd->SetUniformBMap(fld->IsUniform());
1709 fhltesd->SetBInfoStored();
1712 // Set most probable pt, for B=0 tracking
1713 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1714 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1715 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1717 // Fill raw-data error log into the ESD
1718 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1721 if (fRunVertexFinder) {
1722 if (!RunVertexFinder(fesd)) {
1723 if (fStopOnError) {CleanUp(); return kFALSE;}
1727 // For Plane Efficiency: run the SPD trackleter
1728 if (fRunPlaneEff && fSPDTrackleter) {
1729 if (!RunSPDTrackleting(fesd)) {
1730 if (fStopOnError) {CleanUp(); return kFALSE;}
1735 if (!fRunTracking.IsNull()) {
1736 if (fRunMuonTracking) {
1737 if (!RunMuonTracking(fesd)) {
1738 if (fStopOnError) {CleanUp(); return kFALSE;}
1744 if (!fRunTracking.IsNull()) {
1745 if (!RunTracking(fesd)) {
1746 if (fStopOnError) {CleanUp(); return kFALSE;}
1751 if (!fFillESD.IsNull()) {
1752 TString detectors=fFillESD;
1753 // run HLT first and on hltesd
1754 // ;-( IsSelected changes the string
1755 if (IsSelected("HLT", detectors) &&
1756 !FillESD(fhltesd, "HLT")) {
1757 if (fStopOnError) {CleanUp(); return kFALSE;}
1760 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1761 if (detectors.Contains("ALL")) {
1763 for (Int_t idet=0; idet<kNDetectors; ++idet){
1764 detectors += fgkDetectorName[idet];
1768 detectors.ReplaceAll("HLT", "");
1769 if (!FillESD(fesd, detectors)) {
1770 if (fStopOnError) {CleanUp(); return kFALSE;}
1775 AliESDpid::MakePID(fesd);
1777 if (fFillTriggerESD) {
1778 if (!FillTriggerESD(fesd)) {
1779 if (fStopOnError) {CleanUp(); return kFALSE;}
1782 // Always fill scalers
1783 if (!FillTriggerScalers(fesd)) {
1784 if (fStopOnError) {CleanUp(); return kFALSE;}
1791 // Propagate track to the beam pipe (if not already done by ITS)
1793 const Int_t ntracks = fesd->GetNumberOfTracks();
1794 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1797 UShort_t *selectedIdx=new UShort_t[ntracks];
1799 for (Int_t itrack=0; itrack<ntracks; itrack++){
1800 const Double_t kMaxStep = 1; //max step over the material
1803 AliESDtrack *track = fesd->GetTrack(itrack);
1804 if (!track) continue;
1806 AliExternalTrackParam *tpcTrack =
1807 (AliExternalTrackParam *)track->GetTPCInnerParam();
1811 PropagateTrackToBxByBz(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1814 Int_t n=trkArray.GetEntriesFast();
1815 selectedIdx[n]=track->GetID();
1816 trkArray.AddLast(tpcTrack);
1819 //Tracks refitted by ITS should already be at the SPD vertex
1820 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1823 PropagateTrackToBxByBz(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1824 Double_t x[3]; track->GetXYZ(x);
1825 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1826 track->RelateToVertexBxByBz(fesd->GetPrimaryVertexSPD(), b, kVeryBig);
1831 // Improve the reconstructed primary vertex position using the tracks
1833 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1834 if(fesd->GetPrimaryVertexSPD()) {
1835 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1836 if(vtitle.Contains("cosmics")) {
1837 runVertexFinderTracks=kFALSE;
1841 if (runVertexFinderTracks) {
1842 // TPC + ITS primary vertex
1843 ftVertexer->SetITSMode();
1844 ftVertexer->SetConstraintOff();
1845 // get cuts for vertexer from AliGRPRecoParam
1847 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1848 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1849 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1850 ftVertexer->SetCuts(cutsVertexer);
1851 delete [] cutsVertexer; cutsVertexer = NULL;
1852 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1853 ftVertexer->SetVtxStart(fDiamondProfile);
1855 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1857 if (pvtx->GetStatus()) {
1858 fesd->SetPrimaryVertexTracks(pvtx);
1859 for (Int_t i=0; i<ntracks; i++) {
1860 AliESDtrack *t = fesd->GetTrack(i);
1861 Double_t x[3]; t->GetXYZ(x);
1862 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1863 t->RelateToVertexBxByBz(pvtx, b, kVeryBig);
1866 delete pvtx; pvtx=NULL;
1869 // TPC-only primary vertex
1870 ftVertexer->SetTPCMode();
1871 ftVertexer->SetConstraintOff();
1872 // get cuts for vertexer from AliGRPRecoParam
1874 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1875 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1876 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1877 ftVertexer->SetCuts(cutsVertexer);
1878 delete [] cutsVertexer; cutsVertexer = NULL;
1879 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1880 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1882 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1884 if (pvtx->GetStatus()) {
1885 fesd->SetPrimaryVertexTPC(pvtx);
1886 for (Int_t i=0; i<ntracks; i++) {
1887 AliESDtrack *t = fesd->GetTrack(i);
1888 Double_t x[3]; t->GetXYZ(x);
1889 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1890 t->RelateToVertexTPCBxByBz(pvtx, b, kVeryBig);
1893 delete pvtx; pvtx=NULL;
1897 delete[] selectedIdx;
1899 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1904 AliV0vertexer vtxer;
1905 vtxer.Tracks2V0vertices(fesd);
1907 if (fRunCascadeFinder) {
1909 AliCascadeVertexer cvtxer;
1910 cvtxer.V0sTracks2CascadeVertices(fesd);
1915 if (fCleanESD) CleanESD(fesd);
1917 if (fRunQA && IsInTasks(AliQAv1::kESDS)) {
1918 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1919 AliQAManager::QAManager()->RunOneEvent(fesd, fhltesd) ;
1922 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1923 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1924 if (qadm && IsInTasks(AliQAv1::kESDS))
1925 qadm->Exec(AliQAv1::kESDS, fesd);
1928 // copy HLT decision from HLTesd to esd
1929 // the most relevant information is stored in a reduced container in the esd,
1930 // while the full information can be found in the HLTesd
1931 TObject* pHLTSrc=fhltesd->FindListObject(AliESDHLTDecision::Name());
1932 TObject* pHLTTgt=fesd->FindListObject(AliESDHLTDecision::Name());
1933 if (pHLTSrc && pHLTTgt) {
1934 pHLTSrc->Copy(*pHLTTgt);
1937 if (fWriteESDfriend) {
1938 // fesdf->~AliESDfriend();
1939 // new (fesdf) AliESDfriend(); // Reset...
1940 fesd->GetESDfriend(fesdf);
1944 // Auto-save the ESD tree in case of prompt reco @P2
1945 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1946 ftree->AutoSave("SaveSelf");
1947 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1948 if (friendfile) friendfile->Save();
1955 if (fRunAliEVE) RunAliEVE();
1959 if (fWriteESDfriend) {
1960 fesdf->~AliESDfriend();
1961 new (fesdf) AliESDfriend(); // Reset...
1964 ProcInfo_t procInfo;
1965 gSystem->GetProcInfo(&procInfo);
1966 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1969 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1970 if (fReconstructor[iDet]) {
1971 fReconstructor[iDet]->SetRecoParam(NULL);
1972 fReconstructor[iDet]->SetEventInfo(NULL);
1974 if (fTracker[iDet]) fTracker[iDet]->SetEventInfo(NULL);
1977 if (fRunQA || fRunGlobalQA)
1978 AliQAManager::QAManager()->Increment() ;
1983 //_____________________________________________________________________________
1984 void AliReconstruction::SlaveTerminate()
1986 // Finalize the run on the slave side
1987 // Called after the exit
1988 // from the event loop
1989 AliCodeTimerAuto("",0);
1991 if (fIsNewRunLoader) { // galice.root didn't exist
1992 fRunLoader->WriteHeader("OVERWRITE");
1993 fRunLoader->CdGAFile();
1994 fRunLoader->Write(0, TObject::kOverwrite);
1997 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1998 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
2000 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
2001 cdbMapCopy->SetOwner(1);
2002 cdbMapCopy->SetName("cdbMap");
2003 TIter iter(cdbMap->GetTable());
2006 while((pair = dynamic_cast<TPair*> (iter.Next()))){
2007 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
2008 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
2009 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
2012 TList *cdbListCopy = new TList();
2013 cdbListCopy->SetOwner(1);
2014 cdbListCopy->SetName("cdbList");
2016 TIter iter2(cdbList);
2019 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
2020 cdbListCopy->Add(new TObjString(id->ToString().Data()));
2023 ftree->GetUserInfo()->Add(cdbMapCopy);
2024 ftree->GetUserInfo()->Add(cdbListCopy);
2029 if (fWriteESDfriend)
2030 ftree->SetBranchStatus("ESDfriend*",0);
2031 // we want to have only one tree version number
2032 ftree->Write(ftree->GetName(),TObject::kOverwrite);
2033 fhlttree->Write(fhlttree->GetName(),TObject::kOverwrite);
2035 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
2036 if (fRunPlaneEff && !FinishPlaneEff()) {
2037 AliWarning("Finish PlaneEff evaluation failed");
2040 // End of cycle for the in-loop
2042 if (fRunQA || fRunGlobalQA) {
2043 AliQAManager::QAManager()->EndOfCycle() ;
2045 !fProofOutputLocation.IsNull() &&
2046 fProofOutputArchive.IsNull() &&
2047 !fProofOutputDataset) {
2048 TString qaOutputFile(Form("%sMerged.%s.Data.root",
2049 fProofOutputLocation.Data(),
2050 AliQAv1::GetQADataFileName()));
2051 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",
2052 AliQAv1::GetQADataFileName()));
2053 qaProofFile->SetOutputFileName(qaOutputFile.Data());
2054 if (AliDebugLevel() > 0) qaProofFile->Dump();
2055 fOutput->Add(qaProofFile);
2056 MergeQA(qaProofFile->GetFileName());
2067 if (!fProofOutputFileName.IsNull() &&
2068 !fProofOutputLocation.IsNull() &&
2069 fProofOutputDataset &&
2070 !fProofOutputArchive.IsNull()) {
2071 TProofOutputFile *zipProofFile = new TProofOutputFile(fProofOutputFileName.Data(),
2073 fProofOutputLocation.Data());
2074 if (AliDebugLevel() > 0) zipProofFile->Dump();
2075 fOutput->Add(zipProofFile);
2076 TString fileList(fProofOutputArchive.Data());
2077 fileList.ReplaceAll(","," ");
2079 #if ROOT_SVN_REVISION >= 30174
2080 command.Form("zip -n root %s/%s %s",zipProofFile->GetDir(kTRUE),zipProofFile->GetFileName(),fileList.Data());
2082 command.Form("zip -n root %s/%s %s",zipProofFile->GetDir(),zipProofFile->GetFileName(),fileList.Data());
2084 AliInfo(Form("Executing: %s",command.Data()));
2085 gSystem->Exec(command.Data());
2090 //_____________________________________________________________________________
2091 void AliReconstruction::Terminate()
2093 // Create tags for the events in the ESD tree (the ESD tree is always present)
2094 // In case of empty events the tags will contain dummy values
2095 AliCodeTimerAuto("",0);
2097 // Do not call the ESD tag creator in case of PROOF-based reconstruction
2099 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
2100 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
2101 delete esdtagCreator;
2104 // Cleanup of CDB manager: cache and active storages!
2105 AliCDBManager::Instance()->ClearCache();
2108 //_____________________________________________________________________________
2109 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
2111 // run the local reconstruction
2113 static Int_t eventNr=0;
2114 AliCodeTimerAuto("",0)
2116 TString detStr = detectors;
2117 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2118 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2119 AliReconstructor* reconstructor = GetReconstructor(iDet);
2120 if (!reconstructor) continue;
2121 AliLoader* loader = fLoader[iDet];
2122 // Matthias April 2008: temporary fix to run HLT reconstruction
2123 // although the HLT loader is missing
2124 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2126 reconstructor->Reconstruct(fRawReader, NULL);
2129 reconstructor->Reconstruct(dummy, NULL);
2134 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2137 // conversion of digits
2138 if (fRawReader && reconstructor->HasDigitConversion()) {
2139 AliInfo(Form("converting raw data digits into root objects for %s",
2140 fgkDetectorName[iDet]));
2141 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2142 // fgkDetectorName[iDet]),0);
2143 loader->LoadDigits("update");
2144 loader->CleanDigits();
2145 loader->MakeDigitsContainer();
2146 TTree* digitsTree = loader->TreeD();
2147 reconstructor->ConvertDigits(fRawReader, digitsTree);
2148 loader->WriteDigits("OVERWRITE");
2149 loader->UnloadDigits();
2151 // local reconstruction
2152 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2153 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]),0);
2154 loader->LoadRecPoints("update");
2155 loader->CleanRecPoints();
2156 loader->MakeRecPointsContainer();
2157 TTree* clustersTree = loader->TreeR();
2158 if (fRawReader && !reconstructor->HasDigitConversion()) {
2159 reconstructor->Reconstruct(fRawReader, clustersTree);
2161 loader->LoadDigits("read");
2162 TTree* digitsTree = loader->TreeD();
2164 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2168 reconstructor->Reconstruct(digitsTree, clustersTree);
2169 if (fRunQA && IsInTasks(AliQAv1::kDIGITSR)) {
2170 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2171 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2174 loader->UnloadDigits();
2176 if (fRunQA && IsInTasks(AliQAv1::kRECPOINTS)) {
2177 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2178 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2180 loader->WriteRecPoints("OVERWRITE");
2181 loader->UnloadRecPoints();
2182 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2184 IsSelected("CTP", detStr);
2185 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2186 AliError(Form("the following detectors were not found: %s",
2194 //_____________________________________________________________________________
2195 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2197 // run the SPD trackleting (for SPD efficiency purpouses)
2199 AliCodeTimerAuto("",0)
2201 Double_t vtxPos[3] = {0, 0, 0};
2202 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2204 TArrayF mcVertex(3);
2206 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2207 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2208 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2211 const AliESDVertex *vertex = esd->GetVertex();
2213 AliWarning("Vertex not found");
2216 vertex->GetXYZ(vtxPos);
2217 vertex->GetSigmaXYZ(vtxErr);
2218 if (fSPDTrackleter) {
2219 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2222 fLoader[0]->LoadRecPoints("read");
2223 TTree* tree = fLoader[0]->TreeR();
2225 AliError("Can't get the ITS cluster tree");
2228 fSPDTrackleter->LoadClusters(tree);
2229 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2231 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2232 AliError("AliITSTrackleterSPDEff Clusters2Tracks failed");
2233 // fLoader[0]->UnloadRecPoints();
2236 //fSPDTrackleter->UnloadRecPoints();
2238 AliWarning("SPDTrackleter not available");
2244 //_____________________________________________________________________________
2245 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2247 // run the barrel tracking
2249 AliCodeTimerAuto("",0)
2251 AliVertexer *vertexer = CreateVertexer();
2252 if (!vertexer) return kFALSE;
2254 AliInfo("running the ITS vertex finder");
2255 AliESDVertex* vertex = NULL;
2257 fLoader[0]->LoadRecPoints();
2258 TTree* cltree = fLoader[0]->TreeR();
2260 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2261 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2264 AliError("Can't get the ITS cluster tree");
2266 fLoader[0]->UnloadRecPoints();
2269 AliError("Can't get the ITS loader");
2272 AliWarning("Vertex not found");
2273 vertex = new AliESDVertex();
2274 vertex->SetName("default");
2277 vertex->SetName("reconstructed");
2282 vertex->GetXYZ(vtxPos);
2283 vertex->GetSigmaXYZ(vtxErr);
2285 esd->SetPrimaryVertexSPD(vertex);
2286 AliESDVertex *vpileup = NULL;
2287 Int_t novertices = 0;
2288 vpileup = vertexer->GetAllVertices(novertices);
2290 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2292 // if SPD multiplicity has been determined, it is stored in the ESD
2293 AliMultiplicity *mult = vertexer->GetMultiplicity();
2294 if(mult)esd->SetMultiplicity(mult);
2296 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2297 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2306 //_____________________________________________________________________________
2307 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2309 // run the HLT barrel tracking
2311 AliCodeTimerAuto("",0)
2314 AliError("Missing runLoader!");
2318 AliInfo("running HLT tracking");
2320 // Get a pointer to the HLT reconstructor
2321 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2322 if (!reconstructor) return kFALSE;
2325 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2326 TString detName = fgkDetectorName[iDet];
2327 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2328 reconstructor->SetOption(detName.Data());
2329 AliTracker *tracker = reconstructor->CreateTracker();
2331 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2332 if (fStopOnError) return kFALSE;
2336 Double_t vtxErr[3]={0.005,0.005,0.010};
2337 const AliESDVertex *vertex = esd->GetVertex();
2338 vertex->GetXYZ(vtxPos);
2339 tracker->SetVertex(vtxPos,vtxErr);
2341 fLoader[iDet]->LoadRecPoints("read");
2342 TTree* tree = fLoader[iDet]->TreeR();
2344 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2347 tracker->LoadClusters(tree);
2349 if (tracker->Clusters2Tracks(esd) != 0) {
2350 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2354 tracker->UnloadClusters();
2362 //_____________________________________________________________________________
2363 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2365 // run the muon spectrometer tracking
2367 AliCodeTimerAuto("",0)
2370 AliError("Missing runLoader!");
2373 Int_t iDet = 7; // for MUON
2375 AliInfo("is running...");
2377 // Get a pointer to the MUON reconstructor
2378 AliReconstructor *reconstructor = GetReconstructor(iDet);
2379 if (!reconstructor) return kFALSE;
2382 TString detName = fgkDetectorName[iDet];
2383 AliDebug(1, Form("%s tracking", detName.Data()));
2384 AliTracker *tracker = reconstructor->CreateTracker();
2386 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2391 fLoader[iDet]->LoadRecPoints("read");
2393 tracker->LoadClusters(fLoader[iDet]->TreeR());
2395 Int_t rv = tracker->Clusters2Tracks(esd);
2397 fLoader[iDet]->UnloadRecPoints();
2399 tracker->UnloadClusters();
2405 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2413 //_____________________________________________________________________________
2414 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2416 // run the barrel tracking
2417 static Int_t eventNr=0;
2418 AliCodeTimerAuto("",0)
2420 AliInfo("running tracking");
2422 // Set the event info which is used
2423 // by the trackers in order to obtain
2424 // information about read-out detectors,
2426 AliDebug(1, "Setting event info");
2427 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2428 if (!fTracker[iDet]) continue;
2429 fTracker[iDet]->SetEventInfo(&fEventInfo);
2432 //Fill the ESD with the T0 info (will be used by the TOF)
2433 if (fReconstructor[11] && fLoader[11]) {
2434 fLoader[11]->LoadRecPoints("READ");
2435 TTree *treeR = fLoader[11]->TreeR();
2437 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2441 // pass 1: TPC + ITS inwards
2442 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2443 if (!fTracker[iDet]) continue;
2444 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2447 fLoader[iDet]->LoadRecPoints("read");
2448 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2449 TTree* tree = fLoader[iDet]->TreeR();
2451 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2454 fTracker[iDet]->LoadClusters(tree);
2455 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2457 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2458 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2461 // preliminary PID in TPC needed by the ITS tracker
2463 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2464 AliESDpid::MakePID(esd);
2466 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2469 // pass 2: ALL backwards
2471 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2472 if (!fTracker[iDet]) continue;
2473 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2476 if (iDet > 1) { // all except ITS, TPC
2478 fLoader[iDet]->LoadRecPoints("read");
2479 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2480 tree = fLoader[iDet]->TreeR();
2482 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2485 fTracker[iDet]->LoadClusters(tree);
2486 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2490 if (iDet>1) // start filling residuals for the "outer" detectors
2492 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2493 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2495 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2496 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2497 if ( elem && (! elem->At(0)) ) {
2498 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2499 if (qadm) qadm->InitRecPointsForTracker() ;
2503 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2504 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2509 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2510 fTracker[iDet]->UnloadClusters();
2511 fLoader[iDet]->UnloadRecPoints();
2513 // updated PID in TPC needed by the ITS tracker -MI
2515 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2516 AliESDpid::MakePID(esd);
2518 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2520 //stop filling residuals for the "outer" detectors
2521 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2523 // pass 3: TRD + TPC + ITS refit inwards
2525 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2526 if (!fTracker[iDet]) continue;
2527 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2530 if (iDet<2) // start filling residuals for TPC and ITS
2532 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2533 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2535 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2536 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2537 if ( elem && (! elem->At(0)) ) {
2538 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2539 if (qadm) qadm->InitRecPointsForTracker() ;
2544 if (fTracker[iDet]->RefitInward(esd) != 0) {
2545 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2548 // run postprocessing
2549 if (fTracker[iDet]->PostProcess(esd) != 0) {
2550 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2553 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2556 // write space-points to the ESD in case alignment data output
2558 if (fWriteAlignmentData)
2559 WriteAlignmentData(esd);
2561 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2562 if (!fTracker[iDet]) continue;
2564 fTracker[iDet]->UnloadClusters();
2565 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2566 fLoader[iDet]->UnloadRecPoints();
2567 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2569 // stop filling residuals for TPC and ITS
2570 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2576 //_____________________________________________________________________________
2577 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2579 // Remove the data which are not needed for the physics analysis.
2582 Int_t nTracks=esd->GetNumberOfTracks();
2583 Int_t nV0s=esd->GetNumberOfV0s();
2585 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2587 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2588 Bool_t rc=esd->Clean(cleanPars);
2590 nTracks=esd->GetNumberOfTracks();
2591 nV0s=esd->GetNumberOfV0s();
2593 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2598 //_____________________________________________________________________________
2599 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2601 // fill the event summary data
2603 AliCodeTimerAuto("",0)
2604 static Int_t eventNr=0;
2605 TString detStr = detectors;
2607 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2608 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2609 AliReconstructor* reconstructor = GetReconstructor(iDet);
2610 if (!reconstructor) continue;
2611 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2612 TTree* clustersTree = NULL;
2613 if (fLoader[iDet]) {
2614 fLoader[iDet]->LoadRecPoints("read");
2615 clustersTree = fLoader[iDet]->TreeR();
2616 if (!clustersTree) {
2617 AliError(Form("Can't get the %s clusters tree",
2618 fgkDetectorName[iDet]));
2619 if (fStopOnError) return kFALSE;
2622 if (fRawReader && !reconstructor->HasDigitConversion()) {
2623 reconstructor->FillESD(fRawReader, clustersTree, esd);
2625 TTree* digitsTree = NULL;
2626 if (fLoader[iDet]) {
2627 fLoader[iDet]->LoadDigits("read");
2628 digitsTree = fLoader[iDet]->TreeD();
2630 AliError(Form("Can't get the %s digits tree",
2631 fgkDetectorName[iDet]));
2632 if (fStopOnError) return kFALSE;
2635 reconstructor->FillESD(digitsTree, clustersTree, esd);
2636 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2638 if (fLoader[iDet]) {
2639 fLoader[iDet]->UnloadRecPoints();
2643 IsSelected("CTP", detStr);
2644 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2645 AliError(Form("the following detectors were not found: %s",
2647 if (fStopOnError) return kFALSE;
2649 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2654 //_____________________________________________________________________________
2655 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2657 // Reads the trigger decision which is
2658 // stored in Trigger.root file and fills
2659 // the corresponding esd entries
2661 AliCodeTimerAuto("",0)
2663 AliInfo("Filling trigger information into the ESD");
2666 AliCTPRawStream input(fRawReader);
2667 if (!input.Next()) {
2668 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2671 if (esd->GetTriggerMask() != input.GetClassMask())
2672 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2673 input.GetClassMask(),esd->GetTriggerMask()));
2674 if (esd->GetOrbitNumber() != input.GetOrbitID())
2675 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2676 input.GetOrbitID(),esd->GetOrbitNumber()));
2677 if (esd->GetBunchCrossNumber() != input.GetBCID())
2678 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2679 input.GetBCID(),esd->GetBunchCrossNumber()));
2680 AliESDHeader* esdheader = esd->GetHeader();
2681 esdheader->SetL0TriggerInputs(input.GetL0Inputs());
2682 esdheader->SetL1TriggerInputs(input.GetL1Inputs());
2683 esdheader->SetL2TriggerInputs(input.GetL2Inputs());
2685 UInt_t orbit=input.GetOrbitID();
2686 for(Int_t i=0 ; i<input.GetNIRs() ; i++ )
2687 if(TMath::Abs(Int_t(orbit-(input.GetIR(i))->GetOrbit()))<=1){
2688 esdheader->AddTriggerIR(input.GetIR(i));
2694 //_____________________________________________________________________________
2695 Bool_t AliReconstruction::FillTriggerScalers(AliESDEvent*& esd)
2698 //fRunScalers->Print();
2699 if(fRunScalers && fRunScalers->CheckRunScalers()){
2700 AliTimeStamp* timestamp = new AliTimeStamp(esd->GetOrbitNumber(), esd->GetPeriodNumber(), esd->GetBunchCrossNumber());
2701 //AliTimeStamp* timestamp = new AliTimeStamp(10308000, 0, (ULong64_t)486238);
2702 AliESDHeader* esdheader = fesd->GetHeader();
2703 for(Int_t i=0;i<50;i++){
2704 if((1ull<<i) & esd->GetTriggerMask()){
2705 AliTriggerScalersESD* scalesd = fRunScalers->GetScalersForEventClass( timestamp, i+1);
2706 if(scalesd)esdheader->SetTriggerScalersRecord(scalesd);
2712 //_____________________________________________________________________________
2713 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2716 // Filling information from RawReader Header
2719 if (!fRawReader) return kFALSE;
2721 AliInfo("Filling information from RawReader Header");
2723 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2724 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2725 esd->SetPeriodNumber(fRawReader->GetPeriod());
2727 esd->SetTimeStamp(fRawReader->GetTimestamp());
2728 esd->SetEventType(fRawReader->GetType());
2734 //_____________________________________________________________________________
2735 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2737 // check whether detName is contained in detectors
2738 // if yes, it is removed from detectors
2740 // check if all detectors are selected
2741 if ((detectors.CompareTo("ALL") == 0) ||
2742 detectors.BeginsWith("ALL ") ||
2743 detectors.EndsWith(" ALL") ||
2744 detectors.Contains(" ALL ")) {
2749 // search for the given detector
2750 Bool_t result = kFALSE;
2751 if ((detectors.CompareTo(detName) == 0) ||
2752 detectors.BeginsWith(detName+" ") ||
2753 detectors.EndsWith(" "+detName) ||
2754 detectors.Contains(" "+detName+" ")) {
2755 detectors.ReplaceAll(detName, "");
2759 // clean up the detectors string
2760 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2761 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2762 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2767 //_____________________________________________________________________________
2768 Bool_t AliReconstruction::InitRunLoader()
2770 // get or create the run loader
2772 if (gAlice) delete gAlice;
2775 TFile *gafile = TFile::Open(fGAliceFileName.Data());
2776 // if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2777 if (gafile) { // galice.root exists
2781 // load all base libraries to get the loader classes
2782 TString libs = gSystem->GetLibraries();
2783 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2784 TString detName = fgkDetectorName[iDet];
2785 if (detName == "HLT") continue;
2786 if (libs.Contains("lib" + detName + "base.so")) continue;
2787 gSystem->Load("lib" + detName + "base.so");
2789 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2791 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2796 fRunLoader->CdGAFile();
2797 fRunLoader->LoadgAlice();
2799 //PH This is a temporary fix to give access to the kinematics
2800 //PH that is needed for the labels of ITS clusters
2801 fRunLoader->LoadHeader();
2802 fRunLoader->LoadKinematics();
2804 } else { // galice.root does not exist
2806 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2808 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2809 AliConfig::GetDefaultEventFolderName(),
2812 AliError(Form("could not create run loader in file %s",
2813 fGAliceFileName.Data()));
2817 fIsNewRunLoader = kTRUE;
2818 fRunLoader->MakeTree("E");
2820 if (fNumberOfEventsPerFile > 0)
2821 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2823 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2829 //_____________________________________________________________________________
2830 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2832 // get the reconstructor object and the loader for a detector
2834 if (fReconstructor[iDet]) {
2835 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2836 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2837 fReconstructor[iDet]->SetRecoParam(par);
2838 fReconstructor[iDet]->SetRunInfo(fRunInfo);
2840 return fReconstructor[iDet];
2843 // load the reconstructor object
2844 TPluginManager* pluginManager = gROOT->GetPluginManager();
2845 TString detName = fgkDetectorName[iDet];
2846 TString recName = "Ali" + detName + "Reconstructor";
2848 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2850 AliReconstructor* reconstructor = NULL;
2851 // first check if a plugin is defined for the reconstructor
2852 TPluginHandler* pluginHandler =
2853 pluginManager->FindHandler("AliReconstructor", detName);
2854 // if not, add a plugin for it
2855 if (!pluginHandler) {
2856 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2857 TString libs = gSystem->GetLibraries();
2858 if (libs.Contains("lib" + detName + "base.so") ||
2859 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2860 pluginManager->AddHandler("AliReconstructor", detName,
2861 recName, detName + "rec", recName + "()");
2863 pluginManager->AddHandler("AliReconstructor", detName,
2864 recName, detName, recName + "()");
2866 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2868 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2869 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2871 if (reconstructor) {
2872 TObject* obj = fOptions.FindObject(detName.Data());
2873 if (obj) reconstructor->SetOption(obj->GetTitle());
2874 reconstructor->SetRunInfo(fRunInfo);
2875 reconstructor->Init();
2876 fReconstructor[iDet] = reconstructor;
2879 // get or create the loader
2880 if (detName != "HLT") {
2881 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2882 if (!fLoader[iDet]) {
2883 AliConfig::Instance()
2884 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2886 // first check if a plugin is defined for the loader
2888 pluginManager->FindHandler("AliLoader", detName);
2889 // if not, add a plugin for it
2890 if (!pluginHandler) {
2891 TString loaderName = "Ali" + detName + "Loader";
2892 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2893 pluginManager->AddHandler("AliLoader", detName,
2894 loaderName, detName + "base",
2895 loaderName + "(const char*, TFolder*)");
2896 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2898 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2900 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2901 fRunLoader->GetEventFolder());
2903 if (!fLoader[iDet]) { // use default loader
2904 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2906 if (!fLoader[iDet]) {
2907 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2908 if (fStopOnError) return NULL;
2910 fRunLoader->AddLoader(fLoader[iDet]);
2911 fRunLoader->CdGAFile();
2912 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2913 fRunLoader->Write(0, TObject::kOverwrite);
2918 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2919 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2920 reconstructor->SetRecoParam(par);
2921 reconstructor->SetRunInfo(fRunInfo);
2923 return reconstructor;
2926 //_____________________________________________________________________________
2927 AliVertexer* AliReconstruction::CreateVertexer()
2929 // create the vertexer
2930 // Please note that the caller is the owner of the
2933 AliVertexer* vertexer = NULL;
2934 AliReconstructor* itsReconstructor = GetReconstructor(0);
2935 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) || fRunTracking.Contains("ITS"))) {
2936 vertexer = itsReconstructor->CreateVertexer();
2939 AliWarning("couldn't create a vertexer for ITS");
2945 //_____________________________________________________________________________
2946 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2948 // create the trackers
2949 AliInfo("Creating trackers");
2951 TString detStr = detectors;
2952 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2953 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2954 AliReconstructor* reconstructor = GetReconstructor(iDet);
2955 if (!reconstructor) continue;
2956 TString detName = fgkDetectorName[iDet];
2957 if (detName == "HLT") {
2958 fRunHLTTracking = kTRUE;
2961 if (detName == "MUON") {
2962 fRunMuonTracking = kTRUE;
2967 fTracker[iDet] = reconstructor->CreateTracker();
2968 if (!fTracker[iDet] && (iDet < 7)) {
2969 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2970 if (fStopOnError) return kFALSE;
2972 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2978 //_____________________________________________________________________________
2979 void AliReconstruction::CleanUp()
2981 // delete trackers and the run loader and close and delete the file
2983 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2984 delete fReconstructor[iDet];
2985 fReconstructor[iDet] = NULL;
2986 fLoader[iDet] = NULL;
2987 delete fTracker[iDet];
2988 fTracker[iDet] = NULL;
2993 delete fSPDTrackleter;
2994 fSPDTrackleter = NULL;
3003 delete fParentRawReader;
3004 fParentRawReader=NULL;
3012 if (AliQAManager::QAManager())
3013 AliQAManager::QAManager()->ShowQA() ;
3014 AliQAManager::Destroy() ;
3018 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
3020 // Write space-points which are then used in the alignment procedures
3021 // For the moment only ITS, TPC, TRD and TOF
3023 Int_t ntracks = esd->GetNumberOfTracks();
3024 for (Int_t itrack = 0; itrack < ntracks; itrack++)
3026 AliESDtrack *track = esd->GetTrack(itrack);
3029 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
3030 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
3031 nsp += (iDet==GetDetIndex("TRD")) ? track->GetTRDntracklets():track->GetNcls(iDet);
3033 if (iDet==GetDetIndex("ITS")) { // ITS "extra" clusters
3034 track->GetClusters(iDet,idx);
3035 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
3040 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
3041 track->SetTrackPointArray(sp);
3043 for (Int_t iDet = 5; iDet >= 0; iDet--) {
3044 AliTracker *tracker = fTracker[iDet];
3045 if (!tracker) continue;
3046 Int_t nspdet = (iDet==GetDetIndex("TRD")) ? track->GetTRDtracklets(idx):track->GetClusters(iDet,idx);
3048 if (iDet==GetDetIndex("ITS")) // ITS "extra" clusters
3049 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
3051 if (nspdet <= 0) continue;
3055 while (isp2 < nspdet) {
3056 Bool_t isvalid=kTRUE;
3058 Int_t index=idx[isp++];
3059 if (index < 0) continue;
3061 TString dets = fgkDetectorName[iDet];
3062 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
3063 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
3064 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
3065 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
3066 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
3068 isvalid = tracker->GetTrackPoint(index,p);
3071 if (!isvalid) continue;
3072 if (iDet==GetDetIndex("ITS") && (isp-1)>=6) p.SetExtra();
3073 sp->AddPoint(isptrack,&p); isptrack++;
3080 //_____________________________________________________________________________
3081 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
3083 // The method reads the raw-data error log
3084 // accumulated within the rawReader.
3085 // It extracts the raw-data errors related to
3086 // the current event and stores them into
3087 // a TClonesArray inside the esd object.
3089 if (!fRawReader) return;
3091 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
3093 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
3095 if (iEvent != log->GetEventNumber()) continue;
3097 esd->AddRawDataErrorLog(log);
3102 //_____________________________________________________________________________
3103 void AliReconstruction::CheckQA()
3105 // check the QA of SIM for this run and remove the detectors
3106 // with status Fatal
3108 // TString newRunLocalReconstruction ;
3109 // TString newRunTracking ;
3110 // TString newFillESD ;
3112 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
3113 // TString detName(AliQAv1::GetDetName(iDet)) ;
3114 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
3115 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
3116 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
3117 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
3119 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
3120 // fRunLocalReconstruction.Contains("ALL") ) {
3121 // newRunLocalReconstruction += detName ;
3122 // newRunLocalReconstruction += " " ;
3124 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
3125 // fRunTracking.Contains("ALL") ) {
3126 // newRunTracking += detName ;
3127 // newRunTracking += " " ;
3129 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
3130 // fFillESD.Contains("ALL") ) {
3131 // newFillESD += detName ;
3132 // newFillESD += " " ;
3136 // fRunLocalReconstruction = newRunLocalReconstruction ;
3137 // fRunTracking = newRunTracking ;
3138 // fFillESD = newFillESD ;
3141 //_____________________________________________________________________________
3142 Int_t AliReconstruction::GetDetIndex(const char* detector)
3144 // return the detector index corresponding to detector
3146 for (index = 0; index < kNDetectors ; index++) {
3147 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3152 //_____________________________________________________________________________
3153 Bool_t AliReconstruction::FinishPlaneEff() {
3155 // Here execute all the necessary operationis, at the end of the tracking phase,
3156 // in case that evaluation of PlaneEfficiencies was required for some detector.
3157 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3159 // This Preliminary version works only FOR ITS !!!!!
3160 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3163 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3166 TString detStr = fLoadCDB;
3167 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3168 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3169 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3170 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
3171 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3172 TString name=planeeff->GetName();
3174 TFile* pefile = TFile::Open(name, "RECREATE");
3175 ret=(Bool_t)planeeff->Write();
3177 if(planeeff->GetCreateHistos()) {
3178 TString hname=planeeff->GetName();
3179 hname+="Histo.root";
3180 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3183 if(fSPDTrackleter) {
3184 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3185 TString name="AliITSPlaneEffSPDtracklet.root";
3186 TFile* pefile = TFile::Open(name, "RECREATE");
3187 ret=(Bool_t)planeeff->Write();
3189 AliESDEvent *dummy=NULL;
3190 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3195 //_____________________________________________________________________________
3196 Bool_t AliReconstruction::InitPlaneEff() {
3198 // Here execute all the necessary operations, before of the tracking phase,
3199 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3200 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3201 // which should be updated/recalculated.
3203 // This Preliminary version will work only FOR ITS !!!!!
3204 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3207 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3209 AliWarning(Form("Implementation of this method not yet completed !! Method return kTRUE"));
3211 fSPDTrackleter = NULL;
3212 TString detStr = fLoadCDB;
3213 if (IsSelected(fgkDetectorName[0], detStr)) {
3214 AliReconstructor* itsReconstructor = GetReconstructor(0);
3215 if (itsReconstructor) {
3216 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3218 if (fSPDTrackleter) {
3219 AliInfo("Trackleter for SPD has been created");
3225 //_____________________________________________________________________________
3226 Bool_t AliReconstruction::InitAliEVE()
3228 // This method should be called only in case
3229 // AliReconstruction is run
3230 // within the alieve environment.
3231 // It will initialize AliEVE in a way
3232 // so that it can visualize event processed
3233 // by AliReconstruction.
3234 // The return flag shows whenever the
3235 // AliEVE initialization was successful or not.
3238 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3239 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3240 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3242 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3243 gROOT->ProcessLine("alieve_online_init()");
3248 //_____________________________________________________________________________
3249 void AliReconstruction::RunAliEVE()
3251 // Runs AliEVE visualisation of
3252 // the current event.
3253 // Should be executed only after
3254 // successful initialization of AliEVE.
3256 AliInfo("Running AliEVE...");
3257 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3261 //_____________________________________________________________________________
3262 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3264 // Allows to run QA for a selected set of detectors
3265 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3266 // all selected detectors run the same selected tasks
3268 if (!detAndAction.Contains(":")) {
3269 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3273 Int_t colon = detAndAction.Index(":") ;
3274 fQADetectors = detAndAction(0, colon) ;
3275 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3276 if (fQATasks.Contains("ALL") ) {
3277 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3279 fQATasks.ToUpper() ;
3281 if ( fQATasks.Contains("RAW") )
3282 tempo = Form("%d ", AliQAv1::kRAWS) ;
3283 if ( fQATasks.Contains("DIGIT") )
3284 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3285 if ( fQATasks.Contains("RECPOINT") )
3286 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3287 if ( fQATasks.Contains("ESD") )
3288 tempo += Form("%d ", AliQAv1::kESDS) ;
3290 if (fQATasks.IsNull()) {
3291 AliInfo("No QA requested\n") ;
3296 TString tempo(fQATasks) ;
3297 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3298 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3299 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3300 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3301 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3306 //_____________________________________________________________________________
3307 Bool_t AliReconstruction::InitRecoParams()
3309 // The method accesses OCDB and retrieves all
3310 // the available reco-param objects from there.
3312 Bool_t isOK = kTRUE;
3314 if (fRecoParam.GetDetRecoParamArray(kNDetectors)) {
3315 AliInfo("Using custom GRP reconstruction parameters");
3318 AliInfo("Loading GRP reconstruction parameter objects");
3320 AliCDBPath path("GRP","Calib","RecoParam");
3321 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3323 AliWarning("Couldn't find GRP RecoParam entry in OCDB");
3327 TObject *recoParamObj = entry->GetObject();
3328 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3329 // GRP has a normal TobjArray of AliDetectorRecoParam objects
3330 // Registering them in AliRecoParam
3331 fRecoParam.AddDetRecoParamArray(kNDetectors,dynamic_cast<TObjArray*>(recoParamObj));
3333 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3334 // GRP has only onse set of reco parameters
3335 // Registering it in AliRecoParam
3336 AliInfo("Single set of GRP reconstruction parameters found");
3337 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3338 fRecoParam.AddDetRecoParam(kNDetectors,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3341 AliError("No valid GRP RecoParam object found in the OCDB");
3348 TString detStr = fLoadCDB;
3349 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3351 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3353 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3354 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3358 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3360 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3361 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3363 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3367 TObject *recoParamObj = entry->GetObject();
3368 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3369 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3370 // Registering them in AliRecoParam
3371 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3373 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3374 // The detector has only onse set of reco parameters
3375 // Registering it in AliRecoParam
3376 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3377 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3378 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3381 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3385 // FIX ME: We have to disable the unloading of reco-param CDB
3386 // entries because QA framework is using them. Has to be fix in
3387 // a way that the QA takes the objects already constructed in
3389 // AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3393 if (AliDebugLevel() > 0) fRecoParam.Print();
3398 //_____________________________________________________________________________
3399 Bool_t AliReconstruction::GetEventInfo()
3401 // Fill the event info object
3403 AliCodeTimerAuto("",0)
3405 AliCentralTrigger *aCTP = NULL;
3407 fEventInfo.SetEventType(fRawReader->GetType());
3409 ULong64_t mask = fRawReader->GetClassMask();
3410 fEventInfo.SetTriggerMask(mask);
3411 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3412 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3414 aCTP = new AliCentralTrigger();
3415 TString configstr("");
3416 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3417 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3421 aCTP->SetClassMask(mask);
3422 aCTP->SetClusterMask(clmask);
3425 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3427 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3428 aCTP = fRunLoader->GetTrigger();
3429 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3430 // get inputs from actp - just get
3431 AliESDHeader* esdheader = fesd->GetHeader();
3432 esdheader->SetL0TriggerInputs(aCTP->GetL0TriggerInputs());
3433 esdheader->SetL1TriggerInputs(aCTP->GetL1TriggerInputs());
3434 esdheader->SetL2TriggerInputs(aCTP->GetL2TriggerInputs());
3435 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3438 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3443 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3445 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3446 if (fRawReader) delete aCTP;
3450 UChar_t clustmask = 0;
3452 ULong64_t trmask = fEventInfo.GetTriggerMask();
3453 const TObjArray& classesArray = config->GetClasses();
3454 Int_t nclasses = classesArray.GetEntriesFast();
3455 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3456 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3458 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3459 fesd->SetTriggerClass(trclass->GetName(),trindex);
3460 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3461 if (trmask & (1ull << trindex)) {
3463 trclasses += trclass->GetName();
3465 clustmask |= trclass->GetCluster()->GetClusterMask();
3469 fEventInfo.SetTriggerClasses(trclasses);
3471 // Set the information in ESD
3472 fesd->SetTriggerMask(trmask);
3473 fesd->SetTriggerCluster(clustmask);
3475 if (!aCTP->CheckTriggeredDetectors()) {
3476 if (fRawReader) delete aCTP;
3480 if (fRawReader) delete aCTP;
3482 // We have to fill also the HLT decision here!!
3488 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3490 // Match the detector list found in the rec.C or the default 'ALL'
3491 // to the list found in the GRP (stored there by the shuttle PP which
3492 // gets the information from ECS)
3493 static TString resultList;
3494 TString detList = detectorList;
3498 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3499 if ((detectorMask >> iDet) & 0x1) {
3500 TString det = AliDAQ::OfflineModuleName(iDet);
3501 if ((detList.CompareTo("ALL") == 0) ||
3502 ((detList.BeginsWith("ALL ") ||
3503 detList.EndsWith(" ALL") ||
3504 detList.Contains(" ALL ")) &&
3505 !(detList.BeginsWith("-"+det+" ") ||
3506 detList.EndsWith(" -"+det) ||
3507 detList.Contains(" -"+det+" "))) ||
3508 (detList.CompareTo(det) == 0) ||
3509 detList.BeginsWith(det+" ") ||
3510 detList.EndsWith(" "+det) ||
3511 detList.Contains( " "+det+" " )) {
3512 if (!resultList.EndsWith(det + " ")) {
3521 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3522 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3523 if ((detList.CompareTo("ALL") == 0) ||
3524 ((detList.BeginsWith("ALL ") ||
3525 detList.EndsWith(" ALL") ||
3526 detList.Contains(" ALL ")) &&
3527 !(detList.BeginsWith("-"+hltDet+" ") ||
3528 detList.EndsWith(" -"+hltDet) ||
3529 detList.Contains(" -"+hltDet+" "))) ||
3530 (detList.CompareTo(hltDet) == 0) ||
3531 detList.BeginsWith(hltDet+" ") ||
3532 detList.EndsWith(" "+hltDet) ||
3533 detList.Contains( " "+hltDet+" " )) {
3534 resultList += hltDet;
3538 return resultList.Data();
3542 //______________________________________________________________________________
3543 void AliReconstruction::Abort(const char *method, EAbort what)
3545 // Abort processing. If what = kAbortProcess, the Process() loop will be
3546 // aborted. If what = kAbortFile, the current file in a chain will be
3547 // aborted and the processing will continue with the next file, if there
3548 // is no next file then Process() will be aborted. Abort() can also be
3549 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3550 // the SlaveTerminate() and Terminate() are always called. The abort flag
3551 // can be checked in these methods using GetAbort().
3553 // The method is overwritten in AliReconstruction for better handling of
3554 // reco specific errors
3556 if (!fStopOnError) return;
3560 TString whyMess = method;
3561 whyMess += " failed! Aborting...";
3563 AliError(whyMess.Data());
3566 TString mess = "Abort";
3567 if (fAbort == kAbortProcess)
3568 mess = "AbortProcess";
3569 else if (fAbort == kAbortFile)
3572 Info(mess, whyMess.Data());
3575 //______________________________________________________________________________
3576 Bool_t AliReconstruction::ProcessEvent(void* event)
3578 // Method that is used in case the event loop
3579 // is steered from outside, for example by AMORE
3580 // 'event' is a pointer to the DATE event in the memory
3582 if (fRawReader) delete fRawReader;
3583 fRawReader = new AliRawReaderDate(event);
3584 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
3591 //______________________________________________________________________________
3592 Bool_t AliReconstruction::ParseOutput()
3594 // The method parses the output file
3595 // location string in order to steer
3596 // properly the selector
3598 TPMERegexp re1("(\\w+\\.zip#\\w+\\.root):([,*\\w+\\.root,*]+)@dataset://(\\w++)");
3599 TPMERegexp re2("(\\w+\\.root)?@?dataset://(\\w++)");
3601 if (re1.Match(fESDOutput) == 4) {
3602 // root archive with output files stored and regustered
3604 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",re1[1].Data()));
3605 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re1[3].Data()));
3606 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3607 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_ARCHIVE",re1[2].Data()));
3608 AliInfo(Form("%s files will be stored within %s in dataset %s",
3613 else if (re2.Match(fESDOutput) == 3) {
3614 // output file stored and registered
3616 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",(re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data()));
3617 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re2[2].Data()));
3618 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3619 AliInfo(Form("%s will be stored in dataset %s",
3620 (re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data(),
3624 if (fESDOutput.IsNull()) {
3625 // Output location not given.
3626 // Assuming xrootd has been already started and
3627 // the output file has to be sent back
3628 // to the client machine
3629 TString esdUrl(Form("root://%s/%s/",
3630 TUrl(gSystem->HostName()).GetHostFQDN(),
3632 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE","AliESDs.root"));
3633 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",esdUrl.Data()));
3634 AliInfo(Form("AliESDs.root will be stored in %s",
3638 // User specified an output location.
3639 // Ones has just to parse it here
3640 TUrl outputUrl(fESDOutput.Data());
3641 TString outputFile(gSystem->BaseName(outputUrl.GetFile()));
3642 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.IsNull() ? "AliESDs.root" : outputFile.Data()));
3643 TString outputLocation(outputUrl.GetUrl());
3644 outputLocation.ReplaceAll(outputFile.Data(),"");
3645 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",outputLocation.Data()));
3646 AliInfo(Form("%s will be stored in %s",
3647 outputFile.IsNull() ? "AliESDs.root" : outputFile.Data(),
3648 outputLocation.Data()));