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");
989 // LHC: "multiply by 120 to get the energy in MeV"
992 TString runType = fGRPData->GetRunType();
993 if (runType==AliGRPObject::GetInvalidString()) {
994 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
998 Int_t activeDetectors = fGRPData->GetDetectorMask();
999 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
1000 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
1001 activeDetectors = 1074790399;
1004 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1008 // Process the list of active detectors
1009 if (activeDetectors) {
1010 UInt_t detMask = activeDetectors;
1011 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1012 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1013 fFillESD = MatchDetectorList(fFillESD,detMask);
1014 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1015 fLoadCDB.Form("%s %s %s %s",
1016 fRunLocalReconstruction.Data(),
1017 fRunTracking.Data(),
1019 fQADetectors.Data());
1020 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1021 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1) &&
1022 !((detMask >> AliDAQ::DetectorID("ITSSDD")) & 0x1) &&
1023 !((detMask >> AliDAQ::DetectorID("ITSSSD")) & 0x1) ) {
1024 // switch off the vertexer
1025 AliInfo("SPD,SDD,SSD is not in the list of active detectors. Vertexer switched off.");
1026 fRunVertexFinder = kFALSE;
1028 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1029 // switch off the reading of CTP raw-data payload
1030 if (fFillTriggerESD) {
1031 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1032 fFillTriggerESD = kFALSE;
1037 AliInfo("===================================================================================");
1038 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1039 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1040 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1041 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1042 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1043 AliInfo("===================================================================================");
1045 //*** Dealing with the magnetic field map
1046 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {
1047 if (TGeoGlobalMagField::Instance()->GetField()->TestBit(AliMagF::kOverrideGRP)) {
1048 AliInfo("ExpertMode!!! GRP information will be ignored !");
1049 AliInfo("ExpertMode!!! Running with the externally locked B field !");
1052 AliInfo("Destroying existing B field instance!");
1053 delete TGeoGlobalMagField::Instance();
1056 if ( !TGeoGlobalMagField::Instance()->IsLocked() ) {
1057 // Construct the field map out of the information retrieved from GRP.
1060 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1061 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1062 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1066 Char_t l3Polarity = fGRPData->GetL3Polarity();
1067 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1068 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1073 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1074 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1075 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1079 Char_t diPolarity = fGRPData->GetDipolePolarity();
1080 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1081 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1085 // read special bits for the polarity convention and map type
1086 Int_t polConvention = fGRPData->IsPolarityConventionLHC() ? AliMagF::kConvLHC : AliMagF::kConvDCS2008;
1087 Bool_t uniformB = fGRPData->IsUniformBMap();
1090 AliMagF* fld = AliMagF::CreateFieldMap(TMath::Abs(l3Current) * (l3Polarity ? -1:1),
1091 TMath::Abs(diCurrent) * (diPolarity ? -1:1),
1092 polConvention,uniformB,beamEnergy, beamType.Data());
1094 TGeoGlobalMagField::Instance()->SetField( fld );
1095 TGeoGlobalMagField::Instance()->Lock();
1096 AliInfo("Running with the B field constructed out of GRP !");
1098 else AliFatal("Failed to create a B field map !");
1100 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1103 //*** Get the diamond profiles from OCDB
1104 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1106 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1108 AliError("No SPD diamond profile found in OCDB!");
1111 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1113 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1115 AliError("No diamond profile found in OCDB!");
1118 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1120 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1122 AliError("No TPC diamond profile found in OCDB!");
1125 entry = AliCDBManager::Instance()->Get("GRP/Calib/CosmicTriggers");
1127 fListOfCosmicTriggers = dynamic_cast<THashTable*>(entry->GetObject());
1129 AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
1132 if (!fListOfCosmicTriggers) {
1133 AliWarning("Can not get list of cosmic triggers from OCDB! Cosmic event specie will be effectively disabled!");
1139 //_____________________________________________________________________________
1140 Bool_t AliReconstruction::LoadCDB()
1142 AliCodeTimerAuto("",0);
1144 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1146 TString detStr = fLoadCDB;
1147 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1148 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1149 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1152 // Temporary fix - one has to define the correct policy in order
1153 // to load the trigger OCDB entries only for the detectors that
1154 // in the trigger or that are needed in order to put correct
1155 // information in ESD
1156 AliCDBManager::Instance()->GetAll("TRIGGER/*/*");
1160 //_____________________________________________________________________________
1161 Bool_t AliReconstruction::LoadTriggerScalersCDB()
1163 AliCodeTimerAuto("",0);
1165 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
1169 AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
1170 fRunScalers = dynamic_cast<AliTriggerRunScalers*> (entry->GetObject());
1172 if (fRunScalers->CorrectScalersOverflow() == 0) AliInfo("32bit Trigger counters corrected for overflow");
1177 //_____________________________________________________________________________
1178 Bool_t AliReconstruction::LoadCTPTimeParamsCDB()
1180 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
1184 AliInfo("Found an AliCTPTimeParams in GRP/CTP/CTPtiming, reading it");
1185 fCTPTimeParams = dynamic_cast<AliCTPTimeParams*> (entry->GetObject());
1192 //_____________________________________________________________________________
1193 Bool_t AliReconstruction::Run(const char* input)
1196 AliCodeTimerAuto("",0);
1199 if (GetAbort() != TSelector::kContinue) return kFALSE;
1201 TChain *chain = NULL;
1202 if (fRawReader && (chain = fRawReader->GetChain())) {
1203 Long64_t nEntries = (fLastEvent < 0) ? (TChain::kBigNumber) : (fLastEvent - fFirstEvent + 1);
1206 // Temporary fix for long raw-data runs (until socket timeout handling in PROOF is revised)
1207 gProof->Exec("gEnv->SetValue(\"Proof.SocketActivityTimeout\",-1)", kTRUE);
1210 gProof->Exec("TGrid::Connect(\"alien://\")",kTRUE);
1212 TMessage::EnableSchemaEvolutionForAll(kTRUE);
1213 gProof->Exec("TMessage::EnableSchemaEvolutionForAll(kTRUE)",kTRUE);
1215 gProof->AddInput(this);
1217 if (!ParseOutput()) return kFALSE;
1219 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1221 chain->Process("AliReconstruction","",nEntries,fFirstEvent);
1224 chain->Process(this,"",nEntries,fFirstEvent);
1229 if (GetAbort() != TSelector::kContinue) return kFALSE;
1231 if (GetAbort() != TSelector::kContinue) return kFALSE;
1232 //******* The loop over events
1233 AliInfo("Starting looping over events");
1235 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1236 (fRawReader && fRawReader->NextEvent())) {
1237 if (!ProcessEvent(iEvent)) {
1238 Abort("ProcessEvent",TSelector::kAbortFile);
1244 if (GetAbort() != TSelector::kContinue) return kFALSE;
1246 if (GetAbort() != TSelector::kContinue) return kFALSE;
1252 //_____________________________________________________________________________
1253 void AliReconstruction::InitRawReader(const char* input)
1255 AliCodeTimerAuto("",0);
1257 // Init raw-reader and
1258 // set the input in case of raw data
1259 if (input) fRawInput = input;
1260 fRawReader = AliRawReader::Create(fRawInput.Data());
1262 if (fRawInput.IsNull()) {
1263 AliInfo("Reconstruction will run over digits");
1266 AliFatal("Can not create raw-data reader ! Exiting...");
1270 if (!fEquipIdMap.IsNull() && fRawReader)
1271 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1273 if (!fUseHLTData.IsNull()) {
1274 // create the RawReaderHLT which performs redirection of HLT input data for
1275 // the specified detectors
1276 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1278 fParentRawReader=fRawReader;
1279 fRawReader=pRawReader;
1281 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1284 AliSysInfo::AddStamp("CreateRawReader");
1287 //_____________________________________________________________________________
1288 void AliReconstruction::InitRun(const char* input)
1290 // Initialization of raw-reader,
1291 // run number, CDB etc.
1292 AliCodeTimerAuto("",0);
1293 AliSysInfo::AddStamp("Start");
1295 // Initialize raw-reader if any
1296 InitRawReader(input);
1298 // Initialize the CDB storage
1301 // Set run number in CDBManager (if it is not already set by the user)
1302 if (!SetRunNumberFromData()) {
1303 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1307 // Set CDB lock: from now on it is forbidden to reset the run number
1308 // or the default storage or to activate any further storage!
1313 //_____________________________________________________________________________
1314 void AliReconstruction::Begin(TTree *)
1316 // Initialize AlReconstruction before
1317 // going into the event loop
1318 // Should follow the TSelector convention
1319 // i.e. initialize only the object on the client side
1320 AliCodeTimerAuto("",0);
1322 AliReconstruction *reco = NULL;
1324 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1327 AliSysInfo::AddStamp("ReadInputInBegin");
1330 // Import ideal TGeo geometry and apply misalignment
1332 TString geom(gSystem->DirName(fGAliceFileName));
1333 geom += "/geometry.root";
1334 AliGeomManager::LoadGeometry(geom.Data());
1336 Abort("LoadGeometry", TSelector::kAbortProcess);
1339 AliSysInfo::AddStamp("LoadGeom");
1340 TString detsToCheck=fRunLocalReconstruction;
1341 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1342 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1345 AliSysInfo::AddStamp("CheckGeom");
1348 if (!MisalignGeometry(fLoadAlignData)) {
1349 Abort("MisalignGeometry", TSelector::kAbortProcess);
1352 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1353 AliSysInfo::AddStamp("MisalignGeom");
1356 Abort("InitGRP", TSelector::kAbortProcess);
1359 AliSysInfo::AddStamp("InitGRP");
1362 Abort("LoadCDB", TSelector::kAbortProcess);
1365 AliSysInfo::AddStamp("LoadCDB");
1367 if (!LoadTriggerScalersCDB()) {
1368 Abort("LoadTriggerScalersCDB", TSelector::kAbortProcess);
1371 AliSysInfo::AddStamp("LoadTriggerScalersCDB");
1373 if (!LoadCTPTimeParamsCDB()) {
1374 Abort("LoadCTPTimeParamsCDB", TSelector::kAbortProcess);
1377 AliSysInfo::AddStamp("LoadCTPTimeParamsCDB");
1379 // Read the reconstruction parameters from OCDB
1380 if (!InitRecoParams()) {
1381 AliWarning("Not all detectors have correct RecoParam objects initialized");
1383 AliSysInfo::AddStamp("InitRecoParams");
1385 if (fInput && gProof) {
1386 if (reco) *reco = *this;
1388 gGeoManager->SetName("Geometry");
1389 gProof->AddInputData(gGeoManager,kTRUE);
1391 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1392 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1393 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1394 magFieldMap->SetName("MagneticFieldMap");
1395 gProof->AddInputData(magFieldMap,kTRUE);
1400 //_____________________________________________________________________________
1401 void AliReconstruction::SlaveBegin(TTree*)
1403 // Initialization related to run-loader,
1404 // vertexer, trackers, recontructors
1405 // In proof mode it is executed on the slave
1406 AliCodeTimerAuto("",0);
1408 TProofOutputFile *outProofFile = NULL;
1410 if (AliDebugLevel() > 0) fInput->Print();
1411 if (AliDebugLevel() > 10) fInput->Dump();
1412 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1415 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1417 AliGeomManager::SetGeometry(tgeo);
1419 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1420 Int_t runNumber = -1;
1421 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1422 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1423 man->SetCacheFlag(kTRUE);
1424 man->SetLock(kTRUE);
1428 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1429 AliMagF *newMap = new AliMagF(*map);
1430 if (!newMap->LoadParameterization()) {
1431 Abort("AliMagF::LoadParameterization", TSelector::kAbortProcess);
1434 TGeoGlobalMagField::Instance()->SetField(newMap);
1435 TGeoGlobalMagField::Instance()->Lock();
1437 if (TNamed *outputFileName = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE"))
1438 fProofOutputFileName = outputFileName->GetTitle();
1439 if (TNamed *outputLocation = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_LOCATION"))
1440 fProofOutputLocation = outputLocation->GetTitle();
1441 if (fInput->FindObject("PROOF_OUTPUTFILE_DATASET"))
1442 fProofOutputDataset = kTRUE;
1443 if (TNamed *archiveList = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_ARCHIVE"))
1444 fProofOutputArchive = archiveList->GetTitle();
1445 if (!fProofOutputFileName.IsNull() &&
1446 !fProofOutputLocation.IsNull() &&
1447 fProofOutputArchive.IsNull()) {
1448 if (!fProofOutputDataset) {
1449 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"M");
1450 outProofFile->SetOutputFileName(Form("%s%s",fProofOutputLocation.Data(),fProofOutputFileName.Data()));
1453 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"DROV",fProofOutputLocation.Data());
1455 if (AliDebugLevel() > 0) outProofFile->Dump();
1456 fOutput->Add(outProofFile);
1458 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1461 // get the run loader
1462 if (!InitRunLoader()) {
1463 Abort("InitRunLoader", TSelector::kAbortProcess);
1466 AliSysInfo::AddStamp("LoadLoader");
1468 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1471 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1472 Abort("CreateTrackers", TSelector::kAbortProcess);
1475 AliSysInfo::AddStamp("CreateTrackers");
1477 // create the ESD output file and tree
1478 if (!outProofFile) {
1479 ffile = TFile::Open("AliESDs.root", "RECREATE");
1480 ffile->SetCompressionLevel(2);
1481 if (!ffile->IsOpen()) {
1482 Abort("OpenESDFile", TSelector::kAbortProcess);
1487 AliInfo(Form("Opening output PROOF file: %s/%s",
1488 outProofFile->GetDir(), outProofFile->GetFileName()));
1489 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1490 Abort(Form("Problems opening output PROOF file: %s/%s",
1491 outProofFile->GetDir(), outProofFile->GetFileName()),
1492 TSelector::kAbortProcess);
1497 ftree = new TTree("esdTree", "Tree with ESD objects");
1498 fesd = new AliESDEvent();
1499 fesd->CreateStdContent();
1501 // add a so far non-std object to the ESD, this will
1502 // become part of the std content
1503 fesd->AddObject(new AliESDHLTDecision);
1505 fesd->WriteToTree(ftree);
1506 if (fWriteESDfriend) {
1508 // Since we add the branch manually we must
1509 // book and add it after WriteToTree
1510 // otherwise it is created twice,
1511 // once via writetotree and once here.
1512 // The case for AliESDfriend is now
1513 // caught also in AlIESDEvent::WriteToTree but
1514 // be careful when changing the name (AliESDfriend is not
1515 // a TNamed so we had to hardwire it)
1516 fesdf = new AliESDfriend();
1517 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1518 br->SetFile("AliESDfriends.root");
1519 fesd->AddObject(fesdf);
1521 ftree->GetUserInfo()->Add(fesd);
1523 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1524 fhltesd = new AliESDEvent();
1525 fhltesd->CreateStdContent();
1527 // read the ESD template from CDB
1528 // HLT is allowed to put non-std content to its ESD, the non-std
1529 // objects need to be created before invocation of WriteToTree in
1530 // order to create all branches. Initialization is done from an
1531 // ESD layout template in CDB
1532 AliCDBManager* man = AliCDBManager::Instance();
1533 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1534 AliCDBEntry* hltESDConfig=NULL;
1535 if (man->GetId(hltESDConfigPath)!=NULL &&
1536 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1537 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1539 // init all internal variables from the list of objects
1540 pESDLayout->GetStdContent();
1542 // copy content and create non-std objects
1543 *fhltesd=*pESDLayout;
1546 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1547 hltESDConfigPath.GetPath().Data()));
1551 fhltesd->WriteToTree(fhlttree);
1552 fhlttree->GetUserInfo()->Add(fhltesd);
1554 ProcInfo_t procInfo;
1555 gSystem->GetProcInfo(&procInfo);
1556 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1559 //Initialize the QA and start of cycle
1560 if (fRunQA || fRunGlobalQA)
1563 //Initialize the Plane Efficiency framework
1564 if (fRunPlaneEff && !InitPlaneEff()) {
1565 Abort("InitPlaneEff", TSelector::kAbortProcess);
1569 if (strcmp(gProgName,"alieve") == 0)
1570 fRunAliEVE = InitAliEVE();
1575 //_____________________________________________________________________________
1576 Bool_t AliReconstruction::Process(Long64_t entry)
1578 // run the reconstruction over a single entry
1579 // from the chain with raw data
1580 AliCodeTimerAuto("",0);
1582 TTree *currTree = fChain->GetTree();
1583 AliRawVEvent *event = NULL;
1584 currTree->SetBranchAddress("rawevent",&event);
1585 currTree->GetEntry(entry);
1586 fRawReader = new AliRawReaderRoot(event);
1587 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1595 //_____________________________________________________________________________
1596 void AliReconstruction::Init(TTree *tree)
1599 AliError("The input tree is not found!");
1605 //_____________________________________________________________________________
1606 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1608 // run the reconstruction over a single event
1609 // The event loop is steered in Run method
1611 AliCodeTimerAuto("",0);
1613 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1614 fRunLoader->SetEventNumber(iEvent);
1615 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1617 fRunLoader->TreeE()->Fill();
1618 if (fRawReader && fRawReader->UseAutoSaveESD())
1619 fRunLoader->TreeE()->AutoSave("SaveSelf");
1622 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1626 AliInfo(Form("processing event %d", iEvent));
1628 fRunLoader->GetEvent(iEvent);
1630 // Fill Event-info object
1632 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo,fListOfCosmicTriggers);
1633 AliInfo(Form("Current event specie: %s",fRecoParam.PrintEventSpecie()));
1635 // Set the reco-params
1637 TString detStr = fLoadCDB;
1638 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1639 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1640 AliReconstructor *reconstructor = GetReconstructor(iDet);
1641 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1642 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1643 reconstructor->SetRecoParam(par);
1644 reconstructor->SetEventInfo(&fEventInfo);
1646 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1647 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(par->GetEventSpecie())) ;
1651 const AliDetectorRecoParam *grppar = fRecoParam.GetDetRecoParam(kNDetectors);
1652 AliQAManager::QAManager()->SetRecoParam(AliQAv1::kGLOBAL, grppar) ;
1653 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(grppar->GetEventSpecie())) ;
1657 if (fRunQA && IsInTasks(AliQAv1::kRAWS)) {
1658 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1659 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1661 // local single event reconstruction
1662 if (!fRunLocalReconstruction.IsNull()) {
1663 TString detectors=fRunLocalReconstruction;
1664 // run HLT event reconstruction first
1665 // ;-( IsSelected changes the string
1666 if (IsSelected("HLT", detectors) &&
1667 !RunLocalEventReconstruction("HLT")) {
1668 if (fStopOnError) {CleanUp(); return kFALSE;}
1670 detectors=fRunLocalReconstruction;
1671 detectors.ReplaceAll("HLT", "");
1672 if (!RunLocalEventReconstruction(detectors)) {
1680 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1681 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1682 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1683 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1685 fesd->SetEventSpecie(fRecoParam.GetEventSpecie());
1686 fhltesd->SetEventSpecie(fRecoParam.GetEventSpecie());
1688 // Set magnetic field from the tracker
1689 fesd->SetMagneticField(AliTracker::GetBz());
1690 fhltesd->SetMagneticField(AliTracker::GetBz());
1692 AliMagF* fld = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1693 if (fld) { // set info needed for field initialization
1694 fesd->SetCurrentL3(fld->GetCurrentSol());
1695 fesd->SetCurrentDip(fld->GetCurrentDip());
1696 fesd->SetBeamEnergy(fld->GetBeamEnergy());
1697 fesd->SetBeamType(fld->GetBeamTypeText());
1698 fesd->SetUniformBMap(fld->IsUniform());
1699 fesd->SetBInfoStored();
1701 fhltesd->SetCurrentL3(fld->GetCurrentSol());
1702 fhltesd->SetCurrentDip(fld->GetCurrentDip());
1703 fhltesd->SetBeamEnergy(fld->GetBeamEnergy());
1704 fhltesd->SetBeamType(fld->GetBeamTypeText());
1705 fhltesd->SetUniformBMap(fld->IsUniform());
1706 fhltesd->SetBInfoStored();
1709 // Set most probable pt, for B=0 tracking
1710 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1711 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1712 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1714 // Fill raw-data error log into the ESD
1715 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1718 if (fRunVertexFinder) {
1719 if (!RunVertexFinder(fesd)) {
1720 if (fStopOnError) {CleanUp(); return kFALSE;}
1724 // For Plane Efficiency: run the SPD trackleter
1725 if (fRunPlaneEff && fSPDTrackleter) {
1726 if (!RunSPDTrackleting(fesd)) {
1727 if (fStopOnError) {CleanUp(); return kFALSE;}
1732 if (!fRunTracking.IsNull()) {
1733 if (fRunMuonTracking) {
1734 if (!RunMuonTracking(fesd)) {
1735 if (fStopOnError) {CleanUp(); return kFALSE;}
1741 if (!fRunTracking.IsNull()) {
1742 if (!RunTracking(fesd)) {
1743 if (fStopOnError) {CleanUp(); return kFALSE;}
1748 if (!fFillESD.IsNull()) {
1749 TString detectors=fFillESD;
1750 // run HLT first and on hltesd
1751 // ;-( IsSelected changes the string
1752 if (IsSelected("HLT", detectors) &&
1753 !FillESD(fhltesd, "HLT")) {
1754 if (fStopOnError) {CleanUp(); return kFALSE;}
1757 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1758 if (detectors.Contains("ALL")) {
1760 for (Int_t idet=0; idet<kNDetectors; ++idet){
1761 detectors += fgkDetectorName[idet];
1765 detectors.ReplaceAll("HLT", "");
1766 if (!FillESD(fesd, detectors)) {
1767 if (fStopOnError) {CleanUp(); return kFALSE;}
1771 // fill Event header information from the RawEventHeader
1772 if (fRawReader){FillRawEventHeaderESD(fesd);}
1773 if (fRawReader){FillRawEventHeaderESD(fhltesd);}
1776 AliESDpid::MakePID(fesd);
1778 if (fFillTriggerESD) {
1779 if (!FillTriggerESD(fesd)) {
1780 if (fStopOnError) {CleanUp(); return kFALSE;}
1783 // Always fill scalers
1784 if (!FillTriggerScalers(fesd)) {
1785 if (fStopOnError) {CleanUp(); return kFALSE;}
1792 // Propagate track to the beam pipe (if not already done by ITS)
1794 const Int_t ntracks = fesd->GetNumberOfTracks();
1795 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1798 UShort_t *selectedIdx=new UShort_t[ntracks];
1800 for (Int_t itrack=0; itrack<ntracks; itrack++){
1801 const Double_t kMaxStep = 1; //max step over the material
1804 AliESDtrack *track = fesd->GetTrack(itrack);
1805 if (!track) continue;
1807 AliExternalTrackParam *tpcTrack =
1808 (AliExternalTrackParam *)track->GetTPCInnerParam();
1812 PropagateTrackToBxByBz(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1815 Int_t n=trkArray.GetEntriesFast();
1816 selectedIdx[n]=track->GetID();
1817 trkArray.AddLast(tpcTrack);
1820 //Tracks refitted by ITS should already be at the SPD vertex
1821 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1824 PropagateTrackToBxByBz(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1825 Double_t x[3]; track->GetXYZ(x);
1826 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1827 track->RelateToVertexBxByBz(fesd->GetPrimaryVertexSPD(), b, kVeryBig);
1832 // Improve the reconstructed primary vertex position using the tracks
1834 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1835 if(fesd->GetPrimaryVertexSPD()) {
1836 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1837 if(vtitle.Contains("cosmics")) {
1838 runVertexFinderTracks=kFALSE;
1842 if (runVertexFinderTracks) {
1843 // TPC + ITS primary vertex
1844 ftVertexer->SetITSMode();
1845 ftVertexer->SetConstraintOff();
1846 // get cuts for vertexer from AliGRPRecoParam
1848 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1849 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1850 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1851 ftVertexer->SetCuts(cutsVertexer);
1852 delete [] cutsVertexer; cutsVertexer = NULL;
1853 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1854 ftVertexer->SetVtxStart(fDiamondProfile);
1856 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1858 if (pvtx->GetStatus()) {
1859 fesd->SetPrimaryVertexTracks(pvtx);
1860 for (Int_t i=0; i<ntracks; i++) {
1861 AliESDtrack *t = fesd->GetTrack(i);
1862 Double_t x[3]; t->GetXYZ(x);
1863 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1864 t->RelateToVertexBxByBz(pvtx, b, kVeryBig);
1867 delete pvtx; pvtx=NULL;
1870 // TPC-only primary vertex
1871 ftVertexer->SetTPCMode();
1872 ftVertexer->SetConstraintOff();
1873 // get cuts for vertexer from AliGRPRecoParam
1875 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1876 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1877 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1878 ftVertexer->SetCuts(cutsVertexer);
1879 delete [] cutsVertexer; cutsVertexer = NULL;
1880 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1881 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1883 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1885 if (pvtx->GetStatus()) {
1886 fesd->SetPrimaryVertexTPC(pvtx);
1887 for (Int_t i=0; i<ntracks; i++) {
1888 AliESDtrack *t = fesd->GetTrack(i);
1889 Double_t x[3]; t->GetXYZ(x);
1890 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1891 t->RelateToVertexTPCBxByBz(pvtx, b, kVeryBig);
1894 delete pvtx; pvtx=NULL;
1898 delete[] selectedIdx;
1900 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1905 AliV0vertexer vtxer;
1906 vtxer.Tracks2V0vertices(fesd);
1908 if (fRunCascadeFinder) {
1910 AliCascadeVertexer cvtxer;
1911 cvtxer.V0sTracks2CascadeVertices(fesd);
1916 if (fCleanESD) CleanESD(fesd);
1918 if (fRunQA && IsInTasks(AliQAv1::kESDS)) {
1919 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1920 AliQAManager::QAManager()->RunOneEvent(fesd, fhltesd) ;
1923 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1924 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1925 if (qadm && IsInTasks(AliQAv1::kESDS))
1926 qadm->Exec(AliQAv1::kESDS, fesd);
1929 // copy HLT decision from HLTesd to esd
1930 // the most relevant information is stored in a reduced container in the esd,
1931 // while the full information can be found in the HLTesd
1932 TObject* pHLTSrc=fhltesd->FindListObject(AliESDHLTDecision::Name());
1933 TObject* pHLTTgt=fesd->FindListObject(AliESDHLTDecision::Name());
1934 if (pHLTSrc && pHLTTgt) {
1935 pHLTSrc->Copy(*pHLTTgt);
1938 if (fWriteESDfriend) {
1939 // fesdf->~AliESDfriend();
1940 // new (fesdf) AliESDfriend(); // Reset...
1941 fesd->GetESDfriend(fesdf);
1945 // Auto-save the ESD tree in case of prompt reco @P2
1946 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1947 ftree->AutoSave("SaveSelf");
1948 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1949 if (friendfile) friendfile->Save();
1956 if (fRunAliEVE) RunAliEVE();
1960 if (fWriteESDfriend) {
1961 fesdf->~AliESDfriend();
1962 new (fesdf) AliESDfriend(); // Reset...
1965 ProcInfo_t procInfo;
1966 gSystem->GetProcInfo(&procInfo);
1967 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, procInfo.fMemResident, procInfo.fMemVirtual));
1970 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1971 if (fReconstructor[iDet]) {
1972 fReconstructor[iDet]->SetRecoParam(NULL);
1973 fReconstructor[iDet]->SetEventInfo(NULL);
1975 if (fTracker[iDet]) fTracker[iDet]->SetEventInfo(NULL);
1978 if (fRunQA || fRunGlobalQA)
1979 AliQAManager::QAManager()->Increment() ;
1984 //_____________________________________________________________________________
1985 void AliReconstruction::SlaveTerminate()
1987 // Finalize the run on the slave side
1988 // Called after the exit
1989 // from the event loop
1990 AliCodeTimerAuto("",0);
1992 if (fIsNewRunLoader) { // galice.root didn't exist
1993 fRunLoader->WriteHeader("OVERWRITE");
1994 fRunLoader->CdGAFile();
1995 fRunLoader->Write(0, TObject::kOverwrite);
1998 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1999 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
2001 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
2002 cdbMapCopy->SetOwner(1);
2003 cdbMapCopy->SetName("cdbMap");
2004 TIter iter(cdbMap->GetTable());
2007 while((pair = dynamic_cast<TPair*> (iter.Next()))){
2008 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
2009 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
2010 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
2013 TList *cdbListCopy = new TList();
2014 cdbListCopy->SetOwner(1);
2015 cdbListCopy->SetName("cdbList");
2017 TIter iter2(cdbList);
2020 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
2021 cdbListCopy->Add(new TObjString(id->ToString().Data()));
2024 ftree->GetUserInfo()->Add(cdbMapCopy);
2025 ftree->GetUserInfo()->Add(cdbListCopy);
2030 if (fWriteESDfriend)
2031 ftree->SetBranchStatus("ESDfriend*",0);
2032 // we want to have only one tree version number
2033 ftree->Write(ftree->GetName(),TObject::kOverwrite);
2034 fhlttree->Write(fhlttree->GetName(),TObject::kOverwrite);
2036 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
2037 if (fRunPlaneEff && !FinishPlaneEff()) {
2038 AliWarning("Finish PlaneEff evaluation failed");
2041 // End of cycle for the in-loop
2043 if (fRunQA || fRunGlobalQA) {
2044 AliQAManager::QAManager()->EndOfCycle() ;
2046 !fProofOutputLocation.IsNull() &&
2047 fProofOutputArchive.IsNull() &&
2048 !fProofOutputDataset) {
2049 TString qaOutputFile(Form("%sMerged.%s.Data.root",
2050 fProofOutputLocation.Data(),
2051 AliQAv1::GetQADataFileName()));
2052 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",
2053 AliQAv1::GetQADataFileName()));
2054 qaProofFile->SetOutputFileName(qaOutputFile.Data());
2055 if (AliDebugLevel() > 0) qaProofFile->Dump();
2056 fOutput->Add(qaProofFile);
2057 MergeQA(qaProofFile->GetFileName());
2068 if (!fProofOutputFileName.IsNull() &&
2069 !fProofOutputLocation.IsNull() &&
2070 fProofOutputDataset &&
2071 !fProofOutputArchive.IsNull()) {
2072 TProofOutputFile *zipProofFile = new TProofOutputFile(fProofOutputFileName.Data(),
2074 fProofOutputLocation.Data());
2075 if (AliDebugLevel() > 0) zipProofFile->Dump();
2076 fOutput->Add(zipProofFile);
2077 TString fileList(fProofOutputArchive.Data());
2078 fileList.ReplaceAll(","," ");
2080 #if ROOT_SVN_REVISION >= 30174
2081 command.Form("zip -n root %s/%s %s",zipProofFile->GetDir(kTRUE),zipProofFile->GetFileName(),fileList.Data());
2083 command.Form("zip -n root %s/%s %s",zipProofFile->GetDir(),zipProofFile->GetFileName(),fileList.Data());
2085 AliInfo(Form("Executing: %s",command.Data()));
2086 gSystem->Exec(command.Data());
2091 //_____________________________________________________________________________
2092 void AliReconstruction::Terminate()
2094 // Create tags for the events in the ESD tree (the ESD tree is always present)
2095 // In case of empty events the tags will contain dummy values
2096 AliCodeTimerAuto("",0);
2098 // Do not call the ESD tag creator in case of PROOF-based reconstruction
2100 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
2101 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
2102 delete esdtagCreator;
2105 // Cleanup of CDB manager: cache and active storages!
2106 AliCDBManager::Instance()->ClearCache();
2109 //_____________________________________________________________________________
2110 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
2112 // run the local reconstruction
2114 static Int_t eventNr=0;
2115 AliCodeTimerAuto("",0)
2117 TString detStr = detectors;
2118 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2119 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2120 AliReconstructor* reconstructor = GetReconstructor(iDet);
2121 if (!reconstructor) continue;
2122 AliLoader* loader = fLoader[iDet];
2123 // Matthias April 2008: temporary fix to run HLT reconstruction
2124 // although the HLT loader is missing
2125 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2127 reconstructor->Reconstruct(fRawReader, NULL);
2130 reconstructor->Reconstruct(dummy, NULL);
2135 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2138 // conversion of digits
2139 if (fRawReader && reconstructor->HasDigitConversion()) {
2140 AliInfo(Form("converting raw data digits into root objects for %s",
2141 fgkDetectorName[iDet]));
2142 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2143 // fgkDetectorName[iDet]),0);
2144 loader->LoadDigits("update");
2145 loader->CleanDigits();
2146 loader->MakeDigitsContainer();
2147 TTree* digitsTree = loader->TreeD();
2148 reconstructor->ConvertDigits(fRawReader, digitsTree);
2149 loader->WriteDigits("OVERWRITE");
2150 loader->UnloadDigits();
2152 // local reconstruction
2153 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2154 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]),0);
2155 loader->LoadRecPoints("update");
2156 loader->CleanRecPoints();
2157 loader->MakeRecPointsContainer();
2158 TTree* clustersTree = loader->TreeR();
2159 if (fRawReader && !reconstructor->HasDigitConversion()) {
2160 reconstructor->Reconstruct(fRawReader, clustersTree);
2162 loader->LoadDigits("read");
2163 TTree* digitsTree = loader->TreeD();
2165 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2169 reconstructor->Reconstruct(digitsTree, clustersTree);
2170 if (fRunQA && IsInTasks(AliQAv1::kDIGITSR)) {
2171 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2172 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2175 loader->UnloadDigits();
2177 if (fRunQA && IsInTasks(AliQAv1::kRECPOINTS)) {
2178 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2179 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2181 loader->WriteRecPoints("OVERWRITE");
2182 loader->UnloadRecPoints();
2183 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2185 IsSelected("CTP", detStr);
2186 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2187 AliError(Form("the following detectors were not found: %s",
2195 //_____________________________________________________________________________
2196 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2198 // run the SPD trackleting (for SPD efficiency purpouses)
2200 AliCodeTimerAuto("",0)
2202 Double_t vtxPos[3] = {0, 0, 0};
2203 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2205 TArrayF mcVertex(3);
2207 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2208 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2209 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2212 const AliESDVertex *vertex = esd->GetVertex();
2214 AliWarning("Vertex not found");
2217 vertex->GetXYZ(vtxPos);
2218 vertex->GetSigmaXYZ(vtxErr);
2219 if (fSPDTrackleter) {
2220 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2223 fLoader[0]->LoadRecPoints("read");
2224 TTree* tree = fLoader[0]->TreeR();
2226 AliError("Can't get the ITS cluster tree");
2229 fSPDTrackleter->LoadClusters(tree);
2230 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2232 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2233 AliError("AliITSTrackleterSPDEff Clusters2Tracks failed");
2234 // fLoader[0]->UnloadRecPoints();
2237 //fSPDTrackleter->UnloadRecPoints();
2239 AliWarning("SPDTrackleter not available");
2245 //_____________________________________________________________________________
2246 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2248 // run the barrel tracking
2250 AliCodeTimerAuto("",0)
2252 AliVertexer *vertexer = CreateVertexer();
2253 if (!vertexer) return kFALSE;
2255 AliInfo("running the ITS vertex finder");
2256 AliESDVertex* vertex = NULL;
2258 fLoader[0]->LoadRecPoints();
2259 TTree* cltree = fLoader[0]->TreeR();
2261 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2262 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2265 AliError("Can't get the ITS cluster tree");
2267 fLoader[0]->UnloadRecPoints();
2270 AliError("Can't get the ITS loader");
2273 AliWarning("Vertex not found");
2274 vertex = new AliESDVertex();
2275 vertex->SetName("default");
2278 vertex->SetName("reconstructed");
2283 vertex->GetXYZ(vtxPos);
2284 vertex->GetSigmaXYZ(vtxErr);
2286 esd->SetPrimaryVertexSPD(vertex);
2287 AliESDVertex *vpileup = NULL;
2288 Int_t novertices = 0;
2289 vpileup = vertexer->GetAllVertices(novertices);
2291 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2293 // if SPD multiplicity has been determined, it is stored in the ESD
2294 AliMultiplicity *mult = vertexer->GetMultiplicity();
2295 if(mult)esd->SetMultiplicity(mult);
2297 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2298 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2307 //_____________________________________________________________________________
2308 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2310 // run the HLT barrel tracking
2312 AliCodeTimerAuto("",0)
2315 AliError("Missing runLoader!");
2319 AliInfo("running HLT tracking");
2321 // Get a pointer to the HLT reconstructor
2322 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2323 if (!reconstructor) return kFALSE;
2326 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2327 TString detName = fgkDetectorName[iDet];
2328 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2329 reconstructor->SetOption(detName.Data());
2330 AliTracker *tracker = reconstructor->CreateTracker();
2332 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2333 if (fStopOnError) return kFALSE;
2337 Double_t vtxErr[3]={0.005,0.005,0.010};
2338 const AliESDVertex *vertex = esd->GetVertex();
2339 vertex->GetXYZ(vtxPos);
2340 tracker->SetVertex(vtxPos,vtxErr);
2342 fLoader[iDet]->LoadRecPoints("read");
2343 TTree* tree = fLoader[iDet]->TreeR();
2345 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2348 tracker->LoadClusters(tree);
2350 if (tracker->Clusters2Tracks(esd) != 0) {
2351 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2355 tracker->UnloadClusters();
2363 //_____________________________________________________________________________
2364 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2366 // run the muon spectrometer tracking
2368 AliCodeTimerAuto("",0)
2371 AliError("Missing runLoader!");
2374 Int_t iDet = 7; // for MUON
2376 AliInfo("is running...");
2378 // Get a pointer to the MUON reconstructor
2379 AliReconstructor *reconstructor = GetReconstructor(iDet);
2380 if (!reconstructor) return kFALSE;
2383 TString detName = fgkDetectorName[iDet];
2384 AliDebug(1, Form("%s tracking", detName.Data()));
2385 AliTracker *tracker = reconstructor->CreateTracker();
2387 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2392 fLoader[iDet]->LoadRecPoints("read");
2394 tracker->LoadClusters(fLoader[iDet]->TreeR());
2396 Int_t rv = tracker->Clusters2Tracks(esd);
2398 fLoader[iDet]->UnloadRecPoints();
2400 tracker->UnloadClusters();
2406 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2414 //_____________________________________________________________________________
2415 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2417 // run the barrel tracking
2418 static Int_t eventNr=0;
2419 AliCodeTimerAuto("",0)
2421 AliInfo("running tracking");
2423 // Set the event info which is used
2424 // by the trackers in order to obtain
2425 // information about read-out detectors,
2427 AliDebug(1, "Setting event info");
2428 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2429 if (!fTracker[iDet]) continue;
2430 fTracker[iDet]->SetEventInfo(&fEventInfo);
2433 //Fill the ESD with the T0 info (will be used by the TOF)
2434 if (fReconstructor[11] && fLoader[11]) {
2435 fLoader[11]->LoadRecPoints("READ");
2436 TTree *treeR = fLoader[11]->TreeR();
2438 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2442 // pass 1: TPC + ITS inwards
2443 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2444 if (!fTracker[iDet]) continue;
2445 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2448 fLoader[iDet]->LoadRecPoints("read");
2449 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2450 TTree* tree = fLoader[iDet]->TreeR();
2452 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2455 fTracker[iDet]->LoadClusters(tree);
2456 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2458 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2459 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2462 // preliminary PID in TPC needed by the ITS tracker
2464 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2465 AliESDpid::MakePID(esd);
2467 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2470 // pass 2: ALL backwards
2472 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2473 if (!fTracker[iDet]) continue;
2474 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2477 if (iDet > 1) { // all except ITS, TPC
2479 fLoader[iDet]->LoadRecPoints("read");
2480 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2481 tree = fLoader[iDet]->TreeR();
2483 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2486 fTracker[iDet]->LoadClusters(tree);
2487 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2491 if (iDet>1) // start filling residuals for the "outer" detectors
2493 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2494 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2496 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2497 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2498 if ( elem && (! elem->At(0)) ) {
2499 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2500 if (qadm) qadm->InitRecPointsForTracker() ;
2504 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2505 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2510 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2511 fTracker[iDet]->UnloadClusters();
2512 fLoader[iDet]->UnloadRecPoints();
2514 // updated PID in TPC needed by the ITS tracker -MI
2516 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2517 AliESDpid::MakePID(esd);
2519 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2521 //stop filling residuals for the "outer" detectors
2522 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2524 // pass 3: TRD + TPC + ITS refit inwards
2526 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2527 if (!fTracker[iDet]) continue;
2528 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2531 if (iDet<2) // start filling residuals for TPC and ITS
2533 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2534 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2536 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2537 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2538 if ( elem && (! elem->At(0)) ) {
2539 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2540 if (qadm) qadm->InitRecPointsForTracker() ;
2545 if (fTracker[iDet]->RefitInward(esd) != 0) {
2546 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2549 // run postprocessing
2550 if (fTracker[iDet]->PostProcess(esd) != 0) {
2551 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2554 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2557 // write space-points to the ESD in case alignment data output
2559 if (fWriteAlignmentData)
2560 WriteAlignmentData(esd);
2562 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2563 if (!fTracker[iDet]) continue;
2565 fTracker[iDet]->UnloadClusters();
2566 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2567 fLoader[iDet]->UnloadRecPoints();
2568 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2570 // stop filling residuals for TPC and ITS
2571 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2577 //_____________________________________________________________________________
2578 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2580 // Remove the data which are not needed for the physics analysis.
2583 Int_t nTracks=esd->GetNumberOfTracks();
2584 Int_t nV0s=esd->GetNumberOfV0s();
2586 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2588 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2589 Bool_t rc=esd->Clean(cleanPars);
2591 nTracks=esd->GetNumberOfTracks();
2592 nV0s=esd->GetNumberOfV0s();
2594 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2599 //_____________________________________________________________________________
2600 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2602 // fill the event summary data
2604 AliCodeTimerAuto("",0)
2605 static Int_t eventNr=0;
2606 TString detStr = detectors;
2608 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2609 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2610 AliReconstructor* reconstructor = GetReconstructor(iDet);
2611 if (!reconstructor) continue;
2612 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2613 TTree* clustersTree = NULL;
2614 if (fLoader[iDet]) {
2615 fLoader[iDet]->LoadRecPoints("read");
2616 clustersTree = fLoader[iDet]->TreeR();
2617 if (!clustersTree) {
2618 AliError(Form("Can't get the %s clusters tree",
2619 fgkDetectorName[iDet]));
2620 if (fStopOnError) return kFALSE;
2623 if (fRawReader && !reconstructor->HasDigitConversion()) {
2624 reconstructor->FillESD(fRawReader, clustersTree, esd);
2626 TTree* digitsTree = NULL;
2627 if (fLoader[iDet]) {
2628 fLoader[iDet]->LoadDigits("read");
2629 digitsTree = fLoader[iDet]->TreeD();
2631 AliError(Form("Can't get the %s digits tree",
2632 fgkDetectorName[iDet]));
2633 if (fStopOnError) return kFALSE;
2636 reconstructor->FillESD(digitsTree, clustersTree, esd);
2637 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2639 if (fLoader[iDet]) {
2640 fLoader[iDet]->UnloadRecPoints();
2644 IsSelected("CTP", detStr);
2645 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2646 AliError(Form("the following detectors were not found: %s",
2648 if (fStopOnError) return kFALSE;
2650 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2655 //_____________________________________________________________________________
2656 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2658 // Reads the trigger decision which is
2659 // stored in Trigger.root file and fills
2660 // the corresponding esd entries
2662 AliCodeTimerAuto("",0)
2664 AliInfo("Filling trigger information into the ESD");
2667 AliCTPRawStream input(fRawReader);
2668 if (!input.Next()) {
2669 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2672 if (esd->GetTriggerMask() != input.GetClassMask())
2673 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2674 input.GetClassMask(),esd->GetTriggerMask()));
2675 if (esd->GetOrbitNumber() != input.GetOrbitID())
2676 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2677 input.GetOrbitID(),esd->GetOrbitNumber()));
2678 if (esd->GetBunchCrossNumber() != input.GetBCID())
2679 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2680 input.GetBCID(),esd->GetBunchCrossNumber()));
2681 AliESDHeader* esdheader = esd->GetHeader();
2682 esdheader->SetL0TriggerInputs(input.GetL0Inputs());
2683 esdheader->SetL1TriggerInputs(input.GetL1Inputs());
2684 esdheader->SetL2TriggerInputs(input.GetL2Inputs());
2686 UInt_t orbit=input.GetOrbitID();
2687 for(Int_t i=0 ; i<input.GetNIRs() ; i++ )
2688 if(TMath::Abs(Int_t(orbit-(input.GetIR(i))->GetOrbit()))<=1){
2689 esdheader->AddTriggerIR(input.GetIR(i));
2695 //_____________________________________________________________________________
2696 Bool_t AliReconstruction::FillTriggerScalers(AliESDEvent*& esd)
2699 //fRunScalers->Print();
2700 if(fRunScalers && fRunScalers->CheckRunScalers()){
2701 AliTimeStamp* timestamp = new AliTimeStamp(esd->GetOrbitNumber(), esd->GetPeriodNumber(), esd->GetBunchCrossNumber());
2702 //AliTimeStamp* timestamp = new AliTimeStamp(10308000, 0, (ULong64_t)486238);
2703 AliESDHeader* esdheader = fesd->GetHeader();
2704 for(Int_t i=0;i<50;i++){
2705 if((1ull<<i) & esd->GetTriggerMask()){
2706 AliTriggerScalersESD* scalesd = fRunScalers->GetScalersForEventClass( timestamp, i+1);
2707 if(scalesd)esdheader->SetTriggerScalersRecord(scalesd);
2713 //_____________________________________________________________________________
2714 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2717 // Filling information from RawReader Header
2720 if (!fRawReader) return kFALSE;
2722 AliInfo("Filling information from RawReader Header");
2724 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2725 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2726 esd->SetPeriodNumber(fRawReader->GetPeriod());
2728 esd->SetTimeStamp(fRawReader->GetTimestamp());
2729 esd->SetEventType(fRawReader->GetType());
2735 //_____________________________________________________________________________
2736 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2738 // check whether detName is contained in detectors
2739 // if yes, it is removed from detectors
2741 // check if all detectors are selected
2742 if ((detectors.CompareTo("ALL") == 0) ||
2743 detectors.BeginsWith("ALL ") ||
2744 detectors.EndsWith(" ALL") ||
2745 detectors.Contains(" ALL ")) {
2750 // search for the given detector
2751 Bool_t result = kFALSE;
2752 if ((detectors.CompareTo(detName) == 0) ||
2753 detectors.BeginsWith(detName+" ") ||
2754 detectors.EndsWith(" "+detName) ||
2755 detectors.Contains(" "+detName+" ")) {
2756 detectors.ReplaceAll(detName, "");
2760 // clean up the detectors string
2761 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2762 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2763 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2768 //_____________________________________________________________________________
2769 Bool_t AliReconstruction::InitRunLoader()
2771 // get or create the run loader
2773 if (gAlice) delete gAlice;
2776 TFile *gafile = TFile::Open(fGAliceFileName.Data());
2777 // if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2778 if (gafile) { // galice.root exists
2782 // load all base libraries to get the loader classes
2783 TString libs = gSystem->GetLibraries();
2784 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2785 TString detName = fgkDetectorName[iDet];
2786 if (detName == "HLT") continue;
2787 if (libs.Contains("lib" + detName + "base.so")) continue;
2788 gSystem->Load("lib" + detName + "base.so");
2790 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2792 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2797 fRunLoader->CdGAFile();
2798 fRunLoader->LoadgAlice();
2800 //PH This is a temporary fix to give access to the kinematics
2801 //PH that is needed for the labels of ITS clusters
2802 fRunLoader->LoadHeader();
2803 fRunLoader->LoadKinematics();
2805 } else { // galice.root does not exist
2807 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2809 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2810 AliConfig::GetDefaultEventFolderName(),
2813 AliError(Form("could not create run loader in file %s",
2814 fGAliceFileName.Data()));
2818 fIsNewRunLoader = kTRUE;
2819 fRunLoader->MakeTree("E");
2821 if (fNumberOfEventsPerFile > 0)
2822 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2824 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2830 //_____________________________________________________________________________
2831 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2833 // get the reconstructor object and the loader for a detector
2835 if (fReconstructor[iDet]) {
2836 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2837 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2838 fReconstructor[iDet]->SetRecoParam(par);
2839 fReconstructor[iDet]->SetRunInfo(fRunInfo);
2841 return fReconstructor[iDet];
2844 // load the reconstructor object
2845 TPluginManager* pluginManager = gROOT->GetPluginManager();
2846 TString detName = fgkDetectorName[iDet];
2847 TString recName = "Ali" + detName + "Reconstructor";
2849 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2851 AliReconstructor* reconstructor = NULL;
2852 // first check if a plugin is defined for the reconstructor
2853 TPluginHandler* pluginHandler =
2854 pluginManager->FindHandler("AliReconstructor", detName);
2855 // if not, add a plugin for it
2856 if (!pluginHandler) {
2857 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2858 TString libs = gSystem->GetLibraries();
2859 if (libs.Contains("lib" + detName + "base.so") ||
2860 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2861 pluginManager->AddHandler("AliReconstructor", detName,
2862 recName, detName + "rec", recName + "()");
2864 pluginManager->AddHandler("AliReconstructor", detName,
2865 recName, detName, recName + "()");
2867 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2869 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2870 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2872 if (reconstructor) {
2873 TObject* obj = fOptions.FindObject(detName.Data());
2874 if (obj) reconstructor->SetOption(obj->GetTitle());
2875 reconstructor->SetRunInfo(fRunInfo);
2876 reconstructor->Init();
2877 fReconstructor[iDet] = reconstructor;
2880 // get or create the loader
2881 if (detName != "HLT") {
2882 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2883 if (!fLoader[iDet]) {
2884 AliConfig::Instance()
2885 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2887 // first check if a plugin is defined for the loader
2889 pluginManager->FindHandler("AliLoader", detName);
2890 // if not, add a plugin for it
2891 if (!pluginHandler) {
2892 TString loaderName = "Ali" + detName + "Loader";
2893 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2894 pluginManager->AddHandler("AliLoader", detName,
2895 loaderName, detName + "base",
2896 loaderName + "(const char*, TFolder*)");
2897 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2899 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2901 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2902 fRunLoader->GetEventFolder());
2904 if (!fLoader[iDet]) { // use default loader
2905 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2907 if (!fLoader[iDet]) {
2908 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2909 if (fStopOnError) return NULL;
2911 fRunLoader->AddLoader(fLoader[iDet]);
2912 fRunLoader->CdGAFile();
2913 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2914 fRunLoader->Write(0, TObject::kOverwrite);
2919 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2920 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2921 reconstructor->SetRecoParam(par);
2922 reconstructor->SetRunInfo(fRunInfo);
2924 return reconstructor;
2927 //_____________________________________________________________________________
2928 AliVertexer* AliReconstruction::CreateVertexer()
2930 // create the vertexer
2931 // Please note that the caller is the owner of the
2934 AliVertexer* vertexer = NULL;
2935 AliReconstructor* itsReconstructor = GetReconstructor(0);
2936 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) || fRunTracking.Contains("ITS"))) {
2937 vertexer = itsReconstructor->CreateVertexer();
2940 AliWarning("couldn't create a vertexer for ITS");
2946 //_____________________________________________________________________________
2947 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2949 // create the trackers
2950 AliInfo("Creating trackers");
2952 TString detStr = detectors;
2953 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2954 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2955 AliReconstructor* reconstructor = GetReconstructor(iDet);
2956 if (!reconstructor) continue;
2957 TString detName = fgkDetectorName[iDet];
2958 if (detName == "HLT") {
2959 fRunHLTTracking = kTRUE;
2962 if (detName == "MUON") {
2963 fRunMuonTracking = kTRUE;
2968 fTracker[iDet] = reconstructor->CreateTracker();
2969 if (!fTracker[iDet] && (iDet < 7)) {
2970 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2971 if (fStopOnError) return kFALSE;
2973 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2979 //_____________________________________________________________________________
2980 void AliReconstruction::CleanUp()
2982 // delete trackers and the run loader and close and delete the file
2984 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2985 delete fReconstructor[iDet];
2986 fReconstructor[iDet] = NULL;
2987 fLoader[iDet] = NULL;
2988 delete fTracker[iDet];
2989 fTracker[iDet] = NULL;
2994 delete fSPDTrackleter;
2995 fSPDTrackleter = NULL;
3004 delete fParentRawReader;
3005 fParentRawReader=NULL;
3013 if (AliQAManager::QAManager())
3014 AliQAManager::QAManager()->ShowQA() ;
3015 AliQAManager::Destroy() ;
3019 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
3021 // Write space-points which are then used in the alignment procedures
3022 // For the moment only ITS, TPC, TRD and TOF
3024 Int_t ntracks = esd->GetNumberOfTracks();
3025 for (Int_t itrack = 0; itrack < ntracks; itrack++)
3027 AliESDtrack *track = esd->GetTrack(itrack);
3030 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
3031 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
3032 nsp += (iDet==GetDetIndex("TRD")) ? track->GetTRDntracklets():track->GetNcls(iDet);
3034 if (iDet==GetDetIndex("ITS")) { // ITS "extra" clusters
3035 track->GetClusters(iDet,idx);
3036 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
3041 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
3042 track->SetTrackPointArray(sp);
3044 for (Int_t iDet = 5; iDet >= 0; iDet--) {
3045 AliTracker *tracker = fTracker[iDet];
3046 if (!tracker) continue;
3047 Int_t nspdet = (iDet==GetDetIndex("TRD")) ? track->GetTRDtracklets(idx):track->GetClusters(iDet,idx);
3049 if (iDet==GetDetIndex("ITS")) // ITS "extra" clusters
3050 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
3052 if (nspdet <= 0) continue;
3056 while (isp2 < nspdet) {
3057 Bool_t isvalid=kTRUE;
3059 Int_t index=idx[isp++];
3060 if (index < 0) continue;
3062 TString dets = fgkDetectorName[iDet];
3063 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
3064 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
3065 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
3066 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
3067 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
3069 isvalid = tracker->GetTrackPoint(index,p);
3072 if (!isvalid) continue;
3073 if (iDet==GetDetIndex("ITS") && (isp-1)>=6) p.SetExtra();
3074 sp->AddPoint(isptrack,&p); isptrack++;
3081 //_____________________________________________________________________________
3082 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
3084 // The method reads the raw-data error log
3085 // accumulated within the rawReader.
3086 // It extracts the raw-data errors related to
3087 // the current event and stores them into
3088 // a TClonesArray inside the esd object.
3090 if (!fRawReader) return;
3092 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
3094 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
3096 if (iEvent != log->GetEventNumber()) continue;
3098 esd->AddRawDataErrorLog(log);
3103 //_____________________________________________________________________________
3104 void AliReconstruction::CheckQA()
3106 // check the QA of SIM for this run and remove the detectors
3107 // with status Fatal
3109 // TString newRunLocalReconstruction ;
3110 // TString newRunTracking ;
3111 // TString newFillESD ;
3113 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
3114 // TString detName(AliQAv1::GetDetName(iDet)) ;
3115 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
3116 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
3117 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
3118 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
3120 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
3121 // fRunLocalReconstruction.Contains("ALL") ) {
3122 // newRunLocalReconstruction += detName ;
3123 // newRunLocalReconstruction += " " ;
3125 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
3126 // fRunTracking.Contains("ALL") ) {
3127 // newRunTracking += detName ;
3128 // newRunTracking += " " ;
3130 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
3131 // fFillESD.Contains("ALL") ) {
3132 // newFillESD += detName ;
3133 // newFillESD += " " ;
3137 // fRunLocalReconstruction = newRunLocalReconstruction ;
3138 // fRunTracking = newRunTracking ;
3139 // fFillESD = newFillESD ;
3142 //_____________________________________________________________________________
3143 Int_t AliReconstruction::GetDetIndex(const char* detector)
3145 // return the detector index corresponding to detector
3147 for (index = 0; index < kNDetectors ; index++) {
3148 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3153 //_____________________________________________________________________________
3154 Bool_t AliReconstruction::FinishPlaneEff() {
3156 // Here execute all the necessary operationis, at the end of the tracking phase,
3157 // in case that evaluation of PlaneEfficiencies was required for some detector.
3158 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3160 // This Preliminary version works only FOR ITS !!!!!
3161 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3164 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3167 TString detStr = fLoadCDB;
3168 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3169 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3170 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3171 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
3172 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3173 TString name=planeeff->GetName();
3175 TFile* pefile = TFile::Open(name, "RECREATE");
3176 ret=(Bool_t)planeeff->Write();
3178 if(planeeff->GetCreateHistos()) {
3179 TString hname=planeeff->GetName();
3180 hname+="Histo.root";
3181 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3184 if(fSPDTrackleter) {
3185 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3186 TString name="AliITSPlaneEffSPDtracklet.root";
3187 TFile* pefile = TFile::Open(name, "RECREATE");
3188 ret=(Bool_t)planeeff->Write();
3190 AliESDEvent *dummy=NULL;
3191 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3196 //_____________________________________________________________________________
3197 Bool_t AliReconstruction::InitPlaneEff() {
3199 // Here execute all the necessary operations, before of the tracking phase,
3200 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3201 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3202 // which should be updated/recalculated.
3204 // This Preliminary version will work only FOR ITS !!!!!
3205 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3208 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3210 AliWarning(Form("Implementation of this method not yet completed !! Method return kTRUE"));
3212 fSPDTrackleter = NULL;
3213 TString detStr = fLoadCDB;
3214 if (IsSelected(fgkDetectorName[0], detStr)) {
3215 AliReconstructor* itsReconstructor = GetReconstructor(0);
3216 if (itsReconstructor) {
3217 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3219 if (fSPDTrackleter) {
3220 AliInfo("Trackleter for SPD has been created");
3226 //_____________________________________________________________________________
3227 Bool_t AliReconstruction::InitAliEVE()
3229 // This method should be called only in case
3230 // AliReconstruction is run
3231 // within the alieve environment.
3232 // It will initialize AliEVE in a way
3233 // so that it can visualize event processed
3234 // by AliReconstruction.
3235 // The return flag shows whenever the
3236 // AliEVE initialization was successful or not.
3239 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3240 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3241 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3243 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3244 gROOT->ProcessLine("alieve_online_init()");
3249 //_____________________________________________________________________________
3250 void AliReconstruction::RunAliEVE()
3252 // Runs AliEVE visualisation of
3253 // the current event.
3254 // Should be executed only after
3255 // successful initialization of AliEVE.
3257 AliInfo("Running AliEVE...");
3258 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3262 //_____________________________________________________________________________
3263 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3265 // Allows to run QA for a selected set of detectors
3266 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3267 // all selected detectors run the same selected tasks
3269 if (!detAndAction.Contains(":")) {
3270 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3274 Int_t colon = detAndAction.Index(":") ;
3275 fQADetectors = detAndAction(0, colon) ;
3276 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3277 if (fQATasks.Contains("ALL") ) {
3278 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3280 fQATasks.ToUpper() ;
3282 if ( fQATasks.Contains("RAW") )
3283 tempo = Form("%d ", AliQAv1::kRAWS) ;
3284 if ( fQATasks.Contains("DIGIT") )
3285 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3286 if ( fQATasks.Contains("RECPOINT") )
3287 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3288 if ( fQATasks.Contains("ESD") )
3289 tempo += Form("%d ", AliQAv1::kESDS) ;
3291 if (fQATasks.IsNull()) {
3292 AliInfo("No QA requested\n") ;
3297 TString tempo(fQATasks) ;
3298 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3299 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3300 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3301 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3302 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3307 //_____________________________________________________________________________
3308 Bool_t AliReconstruction::InitRecoParams()
3310 // The method accesses OCDB and retrieves all
3311 // the available reco-param objects from there.
3313 Bool_t isOK = kTRUE;
3315 if (fRecoParam.GetDetRecoParamArray(kNDetectors)) {
3316 AliInfo("Using custom GRP reconstruction parameters");
3319 AliInfo("Loading GRP reconstruction parameter objects");
3321 AliCDBPath path("GRP","Calib","RecoParam");
3322 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3324 AliWarning("Couldn't find GRP RecoParam entry in OCDB");
3328 TObject *recoParamObj = entry->GetObject();
3329 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3330 // GRP has a normal TobjArray of AliDetectorRecoParam objects
3331 // Registering them in AliRecoParam
3332 fRecoParam.AddDetRecoParamArray(kNDetectors,dynamic_cast<TObjArray*>(recoParamObj));
3334 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3335 // GRP has only onse set of reco parameters
3336 // Registering it in AliRecoParam
3337 AliInfo("Single set of GRP reconstruction parameters found");
3338 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3339 fRecoParam.AddDetRecoParam(kNDetectors,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3342 AliError("No valid GRP RecoParam object found in the OCDB");
3349 TString detStr = fLoadCDB;
3350 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3352 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3354 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3355 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3359 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3361 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3362 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3364 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3368 TObject *recoParamObj = entry->GetObject();
3369 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3370 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3371 // Registering them in AliRecoParam
3372 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3374 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3375 // The detector has only onse set of reco parameters
3376 // Registering it in AliRecoParam
3377 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3378 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3379 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3382 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3386 // FIX ME: We have to disable the unloading of reco-param CDB
3387 // entries because QA framework is using them. Has to be fix in
3388 // a way that the QA takes the objects already constructed in
3390 // AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3394 if (AliDebugLevel() > 0) fRecoParam.Print();
3399 //_____________________________________________________________________________
3400 Bool_t AliReconstruction::GetEventInfo()
3402 // Fill the event info object
3404 AliCodeTimerAuto("",0)
3406 AliCentralTrigger *aCTP = NULL;
3408 fEventInfo.SetEventType(fRawReader->GetType());
3410 ULong64_t mask = fRawReader->GetClassMask();
3411 fEventInfo.SetTriggerMask(mask);
3412 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3413 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3415 aCTP = new AliCentralTrigger();
3416 TString configstr("");
3417 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3418 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3422 aCTP->SetClassMask(mask);
3423 aCTP->SetClusterMask(clmask);
3426 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3428 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3429 aCTP = fRunLoader->GetTrigger();
3430 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3431 // get inputs from actp - just get
3432 AliESDHeader* esdheader = fesd->GetHeader();
3433 esdheader->SetL0TriggerInputs(aCTP->GetL0TriggerInputs());
3434 esdheader->SetL1TriggerInputs(aCTP->GetL1TriggerInputs());
3435 esdheader->SetL2TriggerInputs(aCTP->GetL2TriggerInputs());
3436 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3439 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3444 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3446 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3447 if (fRawReader) delete aCTP;
3451 UChar_t clustmask = 0;
3453 ULong64_t trmask = fEventInfo.GetTriggerMask();
3454 const TObjArray& classesArray = config->GetClasses();
3455 Int_t nclasses = classesArray.GetEntriesFast();
3456 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3457 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3459 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3460 fesd->SetTriggerClass(trclass->GetName(),trindex);
3461 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3462 if (trmask & (1ull << trindex)) {
3464 trclasses += trclass->GetName();
3466 clustmask |= trclass->GetCluster()->GetClusterMask();
3470 fEventInfo.SetTriggerClasses(trclasses);
3472 // Set the information in ESD
3473 fesd->SetTriggerMask(trmask);
3474 fesd->SetTriggerCluster(clustmask);
3476 if (!aCTP->CheckTriggeredDetectors()) {
3477 if (fRawReader) delete aCTP;
3481 if (fRawReader) delete aCTP;
3483 // We have to fill also the HLT decision here!!
3489 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3491 // Match the detector list found in the rec.C or the default 'ALL'
3492 // to the list found in the GRP (stored there by the shuttle PP which
3493 // gets the information from ECS)
3494 static TString resultList;
3495 TString detList = detectorList;
3499 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3500 if ((detectorMask >> iDet) & 0x1) {
3501 TString det = AliDAQ::OfflineModuleName(iDet);
3502 if ((detList.CompareTo("ALL") == 0) ||
3503 ((detList.BeginsWith("ALL ") ||
3504 detList.EndsWith(" ALL") ||
3505 detList.Contains(" ALL ")) &&
3506 !(detList.BeginsWith("-"+det+" ") ||
3507 detList.EndsWith(" -"+det) ||
3508 detList.Contains(" -"+det+" "))) ||
3509 (detList.CompareTo(det) == 0) ||
3510 detList.BeginsWith(det+" ") ||
3511 detList.EndsWith(" "+det) ||
3512 detList.Contains( " "+det+" " )) {
3513 if (!resultList.EndsWith(det + " ")) {
3522 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3523 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3524 if ((detList.CompareTo("ALL") == 0) ||
3525 ((detList.BeginsWith("ALL ") ||
3526 detList.EndsWith(" ALL") ||
3527 detList.Contains(" ALL ")) &&
3528 !(detList.BeginsWith("-"+hltDet+" ") ||
3529 detList.EndsWith(" -"+hltDet) ||
3530 detList.Contains(" -"+hltDet+" "))) ||
3531 (detList.CompareTo(hltDet) == 0) ||
3532 detList.BeginsWith(hltDet+" ") ||
3533 detList.EndsWith(" "+hltDet) ||
3534 detList.Contains( " "+hltDet+" " )) {
3535 resultList += hltDet;
3539 return resultList.Data();
3543 //______________________________________________________________________________
3544 void AliReconstruction::Abort(const char *method, EAbort what)
3546 // Abort processing. If what = kAbortProcess, the Process() loop will be
3547 // aborted. If what = kAbortFile, the current file in a chain will be
3548 // aborted and the processing will continue with the next file, if there
3549 // is no next file then Process() will be aborted. Abort() can also be
3550 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3551 // the SlaveTerminate() and Terminate() are always called. The abort flag
3552 // can be checked in these methods using GetAbort().
3554 // The method is overwritten in AliReconstruction for better handling of
3555 // reco specific errors
3557 if (!fStopOnError) return;
3561 TString whyMess = method;
3562 whyMess += " failed! Aborting...";
3564 AliError(whyMess.Data());
3567 TString mess = "Abort";
3568 if (fAbort == kAbortProcess)
3569 mess = "AbortProcess";
3570 else if (fAbort == kAbortFile)
3573 Info(mess, whyMess.Data());
3576 //______________________________________________________________________________
3577 Bool_t AliReconstruction::ProcessEvent(void* event)
3579 // Method that is used in case the event loop
3580 // is steered from outside, for example by AMORE
3581 // 'event' is a pointer to the DATE event in the memory
3583 if (fRawReader) delete fRawReader;
3584 fRawReader = new AliRawReaderDate(event);
3585 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
3592 //______________________________________________________________________________
3593 Bool_t AliReconstruction::ParseOutput()
3595 // The method parses the output file
3596 // location string in order to steer
3597 // properly the selector
3599 TPMERegexp re1("(\\w+\\.zip#\\w+\\.root):([,*\\w+\\.root,*]+)@dataset://(\\w++)");
3600 TPMERegexp re2("(\\w+\\.root)?@?dataset://(\\w++)");
3602 if (re1.Match(fESDOutput) == 4) {
3603 // root archive with output files stored and regustered
3605 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",re1[1].Data()));
3606 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re1[3].Data()));
3607 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3608 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_ARCHIVE",re1[2].Data()));
3609 AliInfo(Form("%s files will be stored within %s in dataset %s",
3614 else if (re2.Match(fESDOutput) == 3) {
3615 // output file stored and registered
3617 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",(re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data()));
3618 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re2[2].Data()));
3619 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3620 AliInfo(Form("%s will be stored in dataset %s",
3621 (re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data(),
3625 if (fESDOutput.IsNull()) {
3626 // Output location not given.
3627 // Assuming xrootd has been already started and
3628 // the output file has to be sent back
3629 // to the client machine
3630 TString esdUrl(Form("root://%s/%s/",
3631 TUrl(gSystem->HostName()).GetHostFQDN(),
3633 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE","AliESDs.root"));
3634 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",esdUrl.Data()));
3635 AliInfo(Form("AliESDs.root will be stored in %s",
3639 // User specified an output location.
3640 // Ones has just to parse it here
3641 TUrl outputUrl(fESDOutput.Data());
3642 TString outputFile(gSystem->BaseName(outputUrl.GetFile()));
3643 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.IsNull() ? "AliESDs.root" : outputFile.Data()));
3644 TString outputLocation(outputUrl.GetUrl());
3645 outputLocation.ReplaceAll(outputFile.Data(),"");
3646 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",outputLocation.Data()));
3647 AliInfo(Form("%s will be stored in %s",
3648 outputFile.IsNull() ? "AliESDs.root" : outputFile.Data(),
3649 outputLocation.Data()));