1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // class for running the reconstruction //
22 // Clusters and tracks are created for all detectors and all events by //
25 // AliReconstruction rec; //
28 // The Run method returns kTRUE in case of successful execution. //
30 // If the input to the reconstruction are not simulated digits but raw data, //
31 // this can be specified by an argument of the Run method or by the method //
33 // rec.SetInput("..."); //
35 // The input formats and the corresponding argument are: //
36 // - DDL raw data files: directory name, ends with "/" //
37 // - raw data root file: root file name, extension ".root" //
38 // - raw data DATE file: DATE file name, any other non-empty string //
39 // - MC root files : empty string, default //
41 // By default all events are reconstructed. The reconstruction can be //
42 // limited to a range of events by giving the index of the first and the //
43 // last event as an argument to the Run method or by calling //
45 // rec.SetEventRange(..., ...); //
47 // The index -1 (default) can be used for the last event to indicate no //
48 // upper limit of the event range. //
50 // In case of raw-data reconstruction the user can modify the default //
51 // number of events per digits/clusters/tracks file. In case the option //
52 // is not used the number is set 1. In case the user provides 0, than //
53 // the number of events is equal to the number of events inside the //
54 // raw-data file (i.e. one digits/clusters/tracks file): //
56 // rec.SetNumberOfEventsPerFile(...); //
59 // The name of the galice file can be changed from the default //
60 // "galice.root" by passing it as argument to the AliReconstruction //
61 // constructor or by //
63 // rec.SetGAliceFile("..."); //
65 // The local reconstruction can be switched on or off for individual //
68 // rec.SetRunLocalReconstruction("..."); //
70 // The argument is a (case sensitive) string with the names of the //
71 // detectors separated by a space. The special string "ALL" selects all //
72 // available detectors. This is the default. //
74 // The reconstruction of the primary vertex position can be switched off by //
76 // rec.SetRunVertexFinder(kFALSE); //
78 // The tracking and the creation of ESD tracks can be switched on for //
79 // selected detectors by //
81 // rec.SetRunTracking("..."); //
83 // Uniform/nonuniform field tracking switches (default: uniform field) //
85 // rec.SetUniformFieldTracking(); ( rec.SetUniformFieldTracking(kFALSE); ) //
87 // The filling of additional ESD information can be steered by //
89 // rec.SetFillESD("..."); //
91 // Again, for both methods the string specifies the list of detectors. //
92 // The default is "ALL". //
94 // The call of the shortcut method //
96 // rec.SetRunReconstruction("..."); //
98 // is equivalent to calling SetRunLocalReconstruction, SetRunTracking and //
99 // SetFillESD with the same detector selecting string as argument. //
101 // The reconstruction requires digits or raw data as input. For the creation //
102 // of digits and raw data have a look at the class AliSimulation. //
104 // The input data of a detector can be replaced by the corresponding HLT //
105 // data by calling (usual detector string) //
106 // SetUseHLTData("..."); //
109 ///////////////////////////////////////////////////////////////////////////////
116 #include <TGeoGlobalMagField.h>
117 #include <TGeoManager.h>
119 #include <TLorentzVector.h>
121 #include <TObjArray.h>
122 #include <TPRegexp.h>
123 #include <TParameter.h>
124 #include <TPluginManager.h>
126 #include <TProofOutputFile.h>
129 #include <THashTable.h>
131 #include <TMessage.h>
133 #include "AliAlignObj.h"
134 #include "AliCDBEntry.h"
135 #include "AliCDBManager.h"
136 #include "AliCDBStorage.h"
137 #include "AliCTPRawStream.h"
138 #include "AliCascadeVertexer.h"
139 #include "AliCentralTrigger.h"
140 #include "AliCodeTimer.h"
142 #include "AliDetectorRecoParam.h"
143 #include "AliESDCaloCells.h"
144 #include "AliESDCaloCluster.h"
145 #include "AliESDEvent.h"
146 #include "AliESDMuonTrack.h"
147 #include "AliESDPmdTrack.h"
148 #include "AliESDTagCreator.h"
149 #include "AliESDVertex.h"
150 #include "AliESDcascade.h"
151 #include "AliESDfriend.h"
152 #include "AliESDkink.h"
153 #include "AliESDpid.h"
154 #include "AliESDtrack.h"
155 #include "AliESDtrack.h"
156 #include "AliEventInfo.h"
157 #include "AliGRPObject.h"
158 #include "AliGRPRecoParam.h"
159 #include "AliGenEventHeader.h"
160 #include "AliGeomManager.h"
161 #include "AliGlobalQADataMaker.h"
162 #include "AliHeader.h"
165 #include "AliMultiplicity.h"
167 #include "AliPlaneEff.h"
169 #include "AliQADataMakerRec.h"
170 #include "AliQAManager.h"
171 #include "AliRawVEvent.h"
172 #include "AliRawEventHeaderBase.h"
173 #include "AliRawHLTManager.h"
174 #include "AliRawReaderDate.h"
175 #include "AliRawReaderFile.h"
176 #include "AliRawReaderRoot.h"
177 #include "AliReconstruction.h"
178 #include "AliReconstructor.h"
180 #include "AliRunInfo.h"
181 #include "AliRunLoader.h"
182 #include "AliSysInfo.h" // memory snapshots
183 #include "AliTrackPointArray.h"
184 #include "AliTracker.h"
185 #include "AliTriggerClass.h"
186 #include "AliTriggerCluster.h"
187 #include "AliTriggerIR.h"
188 #include "AliTriggerConfiguration.h"
189 #include "AliV0vertexer.h"
190 #include "AliVertexer.h"
191 #include "AliVertexerTracks.h"
192 #include "AliTriggerRunScalers.h"
193 #include "AliCTPTimeParams.h"
194 #include "AliESDHLTDecision.h"
196 ClassImp(AliReconstruction)
198 //_____________________________________________________________________________
199 const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
201 //_____________________________________________________________________________
202 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
204 fRunVertexFinder(kTRUE),
205 fRunVertexFinderTracks(kTRUE),
206 fRunHLTTracking(kFALSE),
207 fRunMuonTracking(kFALSE),
209 fRunCascadeFinder(kTRUE),
211 fWriteAlignmentData(kFALSE),
212 fWriteESDfriend(kFALSE),
213 fFillTriggerESD(kTRUE),
221 fRunLocalReconstruction("ALL"),
225 fUseTrackingErrorsForAlignment(""),
226 fGAliceFileName(gAliceFilename),
229 fProofOutputFileName(""),
230 fProofOutputLocation(""),
231 fProofOutputDataset(kFALSE),
232 fProofOutputArchive(""),
236 fNumberOfEventsPerFile((UInt_t)-1),
238 fLoadAlignFromCDB(kTRUE),
239 fLoadAlignData("ALL"),
244 fCTPTimeParams(NULL),
248 fParentRawReader(NULL),
252 fSPDTrackleter(NULL),
254 fDiamondProfileSPD(NULL),
255 fDiamondProfile(NULL),
256 fDiamondProfileTPC(NULL),
257 fListOfCosmicTriggers(NULL),
261 fAlignObjArray(NULL),
265 fInitCDBCalled(kFALSE),
266 fSetRunNumberFromDataCalled(kFALSE),
271 fSameQACycle(kFALSE),
272 fInitQACalled(kFALSE),
273 fWriteQAExpertData(kTRUE),
274 fRunPlaneEff(kFALSE),
283 fIsNewRunLoader(kFALSE),
287 // create reconstruction object with default parameters
290 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
291 fReconstructor[iDet] = NULL;
292 fLoader[iDet] = NULL;
293 fTracker[iDet] = NULL;
295 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
296 fQACycles[iDet] = 999999 ;
297 fQAWriteExpert[iDet] = kFALSE ;
303 //_____________________________________________________________________________
304 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
306 fRunVertexFinder(rec.fRunVertexFinder),
307 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
308 fRunHLTTracking(rec.fRunHLTTracking),
309 fRunMuonTracking(rec.fRunMuonTracking),
310 fRunV0Finder(rec.fRunV0Finder),
311 fRunCascadeFinder(rec.fRunCascadeFinder),
312 fStopOnError(rec.fStopOnError),
313 fWriteAlignmentData(rec.fWriteAlignmentData),
314 fWriteESDfriend(rec.fWriteESDfriend),
315 fFillTriggerESD(rec.fFillTriggerESD),
317 fCleanESD(rec.fCleanESD),
318 fV0DCAmax(rec.fV0DCAmax),
319 fV0CsPmin(rec.fV0CsPmin),
323 fRunLocalReconstruction(rec.fRunLocalReconstruction),
324 fRunTracking(rec.fRunTracking),
325 fFillESD(rec.fFillESD),
326 fLoadCDB(rec.fLoadCDB),
327 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
328 fGAliceFileName(rec.fGAliceFileName),
329 fRawInput(rec.fRawInput),
330 fESDOutput(rec.fESDOutput),
331 fProofOutputFileName(rec.fProofOutputFileName),
332 fProofOutputLocation(rec.fProofOutputLocation),
333 fProofOutputDataset(rec.fProofOutputDataset),
334 fProofOutputArchive(rec.fProofOutputArchive),
335 fEquipIdMap(rec.fEquipIdMap),
336 fFirstEvent(rec.fFirstEvent),
337 fLastEvent(rec.fLastEvent),
338 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
340 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
341 fLoadAlignData(rec.fLoadAlignData),
342 fUseHLTData(rec.fUseHLTData),
346 fCTPTimeParams(NULL),
350 fParentRawReader(NULL),
352 fRecoParam(rec.fRecoParam),
354 fSPDTrackleter(NULL),
356 fDiamondProfileSPD(rec.fDiamondProfileSPD),
357 fDiamondProfile(rec.fDiamondProfile),
358 fDiamondProfileTPC(rec.fDiamondProfileTPC),
359 fListOfCosmicTriggers(NULL),
363 fAlignObjArray(rec.fAlignObjArray),
364 fCDBUri(rec.fCDBUri),
365 fQARefUri(rec.fQARefUri),
367 fInitCDBCalled(rec.fInitCDBCalled),
368 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
369 fQADetectors(rec.fQADetectors),
370 fQATasks(rec.fQATasks),
372 fRunGlobalQA(rec.fRunGlobalQA),
373 fSameQACycle(rec.fSameQACycle),
374 fInitQACalled(rec.fInitQACalled),
375 fWriteQAExpertData(rec.fWriteQAExpertData),
376 fRunPlaneEff(rec.fRunPlaneEff),
385 fIsNewRunLoader(rec.fIsNewRunLoader),
391 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
392 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
394 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
395 fReconstructor[iDet] = NULL;
396 fLoader[iDet] = NULL;
397 fTracker[iDet] = NULL;
400 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
401 fQACycles[iDet] = rec.fQACycles[iDet];
402 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
405 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
406 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
410 //_____________________________________________________________________________
411 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
413 // assignment operator
414 // Used in PROOF mode
415 // Be very careful while modifing it!
416 // Simple rules to follow:
417 // for persistent data members - use their assignment operators
418 // for non-persistent ones - do nothing or take the default values from constructor
419 // TSelector members should not be touched
420 if(&rec == this) return *this;
422 fRunVertexFinder = rec.fRunVertexFinder;
423 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
424 fRunHLTTracking = rec.fRunHLTTracking;
425 fRunMuonTracking = rec.fRunMuonTracking;
426 fRunV0Finder = rec.fRunV0Finder;
427 fRunCascadeFinder = rec.fRunCascadeFinder;
428 fStopOnError = rec.fStopOnError;
429 fWriteAlignmentData = rec.fWriteAlignmentData;
430 fWriteESDfriend = rec.fWriteESDfriend;
431 fFillTriggerESD = rec.fFillTriggerESD;
433 fCleanESD = rec.fCleanESD;
434 fV0DCAmax = rec.fV0DCAmax;
435 fV0CsPmin = rec.fV0CsPmin;
439 fRunLocalReconstruction = rec.fRunLocalReconstruction;
440 fRunTracking = rec.fRunTracking;
441 fFillESD = rec.fFillESD;
442 fLoadCDB = rec.fLoadCDB;
443 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
444 fGAliceFileName = rec.fGAliceFileName;
445 fRawInput = rec.fRawInput;
446 fESDOutput = rec.fESDOutput;
447 fProofOutputFileName = rec.fProofOutputFileName;
448 fProofOutputLocation = rec.fProofOutputLocation;
449 fProofOutputDataset = rec.fProofOutputDataset;
450 fProofOutputArchive = rec.fProofOutputArchive;
451 fEquipIdMap = rec.fEquipIdMap;
452 fFirstEvent = rec.fFirstEvent;
453 fLastEvent = rec.fLastEvent;
454 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
456 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
457 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
460 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
461 fLoadAlignData = rec.fLoadAlignData;
462 fUseHLTData = rec.fUseHLTData;
464 delete fRunInfo; fRunInfo = NULL;
465 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
467 fEventInfo = rec.fEventInfo;
469 delete fRunScalers; fRunScalers = NULL;
470 if (rec.fRunScalers) fRunScalers = new AliTriggerRunScalers(*rec.fRunScalers);
472 delete fCTPTimeParams; fCTPTimeParams = NULL;
473 if (rec.fCTPTimeParams) fCTPTimeParams = new AliCTPTimeParams(*rec.fCTPTimeParams);
477 fParentRawReader = NULL;
479 fRecoParam = rec.fRecoParam;
481 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
482 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
483 delete fLoader[iDet]; fLoader[iDet] = NULL;
484 delete fTracker[iDet]; fTracker[iDet] = NULL;
487 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
488 fQACycles[iDet] = rec.fQACycles[iDet];
489 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
492 delete fSPDTrackleter; fSPDTrackleter = NULL;
494 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
495 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
496 delete fDiamondProfile; fDiamondProfile = NULL;
497 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
498 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
499 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
501 delete fListOfCosmicTriggers; fListOfCosmicTriggers = NULL;
502 if (rec.fListOfCosmicTriggers) fListOfCosmicTriggers = (THashTable*)((rec.fListOfCosmicTriggers)->Clone());
504 delete fGRPData; fGRPData = NULL;
505 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
506 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
508 delete fAlignObjArray; fAlignObjArray = NULL;
511 fQARefUri = rec.fQARefUri;
512 fSpecCDBUri.Delete();
513 fInitCDBCalled = rec.fInitCDBCalled;
514 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
515 fQADetectors = rec.fQADetectors;
516 fQATasks = rec.fQATasks;
518 fRunGlobalQA = rec.fRunGlobalQA;
519 fSameQACycle = rec.fSameQACycle;
520 fInitQACalled = rec.fInitQACalled;
521 fWriteQAExpertData = rec.fWriteQAExpertData;
522 fRunPlaneEff = rec.fRunPlaneEff;
531 fIsNewRunLoader = rec.fIsNewRunLoader;
538 //_____________________________________________________________________________
539 AliReconstruction::~AliReconstruction()
544 if (fListOfCosmicTriggers) {
545 fListOfCosmicTriggers->Delete();
546 delete fListOfCosmicTriggers;
550 delete fCTPTimeParams;
552 if (fAlignObjArray) {
553 fAlignObjArray->Delete();
554 delete fAlignObjArray;
556 fSpecCDBUri.Delete();
558 AliCodeTimer::Instance()->Print();
561 //_____________________________________________________________________________
562 void AliReconstruction::InitQA()
564 //Initialize the QA and start of cycle
565 AliCodeTimerAuto("",0);
567 if (fInitQACalled) return;
568 fInitQACalled = kTRUE;
570 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kRECMODE) ;
571 if (fWriteQAExpertData)
572 qam->SetWriteExpert() ;
574 if (qam->IsDefaultStorageSet()) {
575 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
576 AliWarning("Default QA reference storage has been already set !");
577 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fQARefUri.Data()));
578 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
579 fQARefUri = qam->GetDefaultStorage()->GetURI();
581 if (fQARefUri.Length() > 0) {
582 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
583 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
584 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
586 fQARefUri="local://$ALICE_ROOT/QAref";
587 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
588 AliWarning("Default QA refeference storage not yet set !!!!");
589 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
590 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
593 qam->SetDefaultStorage(fQARefUri);
597 qam->SetActiveDetectors(fQADetectors) ;
598 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) {
599 qam->SetCycleLength(AliQAv1::DETECTORINDEX_t(det), fQACycles[det]) ;
600 qam->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
602 if (!fRawReader && !fInput && IsInTasks(AliQAv1::kRAWS))
603 fQATasks.ReplaceAll(Form("%d",AliQAv1::kRAWS), "") ;
604 qam->SetTasks(fQATasks) ;
605 qam->InitQADataMaker(AliCDBManager::Instance()->GetRun()) ;
608 Bool_t sameCycle = kFALSE ;
609 AliQADataMaker *qadm = qam->GetQADataMaker(AliQAv1::kGLOBAL);
610 AliInfo(Form("Initializing the global QA data maker"));
611 if (IsInTasks(AliQAv1::kRECPOINTS)) {
612 qadm->StartOfCycle(AliQAv1::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
613 TObjArray **arr=qadm->Init(AliQAv1::kRECPOINTS);
614 AliTracker::SetResidualsArray(arr);
617 if (IsInTasks(AliQAv1::kESDS)) {
618 qadm->StartOfCycle(AliQAv1::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
619 qadm->Init(AliQAv1::kESDS);
622 AliSysInfo::AddStamp("InitQA") ;
625 //_____________________________________________________________________________
626 void AliReconstruction::MergeQA(const char *fileName)
628 //Initialize the QA and start of cycle
629 AliCodeTimerAuto("",0) ;
630 AliQAManager::QAManager()->Merge(AliCDBManager::Instance()->GetRun(),fileName) ;
631 AliSysInfo::AddStamp("MergeQA") ;
634 //_____________________________________________________________________________
635 void AliReconstruction::InitCDB()
637 // activate a default CDB storage
638 // First check if we have any CDB storage set, because it is used
639 // to retrieve the calibration and alignment constants
640 AliCodeTimerAuto("",0);
642 if (fInitCDBCalled) return;
643 fInitCDBCalled = kTRUE;
645 AliCDBManager* man = AliCDBManager::Instance();
646 if (man->IsDefaultStorageSet())
648 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
649 AliWarning("Default CDB storage has been already set !");
650 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
651 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
652 fCDBUri = man->GetDefaultStorage()->GetURI();
655 if (fCDBUri.Length() > 0)
657 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
658 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
659 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
660 man->SetDefaultStorage(fCDBUri);
662 else if (!man->GetRaw()){
663 fCDBUri="local://$ALICE_ROOT/OCDB";
664 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
665 AliWarning("Default CDB storage not yet set !!!!");
666 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
667 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
668 man->SetDefaultStorage(fCDBUri);
671 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
672 AliWarning("Default storage will be set after setting the Run Number!!!");
673 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
677 // Now activate the detector specific CDB storage locations
678 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
679 TObject* obj = fSpecCDBUri[i];
681 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
682 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
683 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
684 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
686 AliSysInfo::AddStamp("InitCDB");
689 //_____________________________________________________________________________
690 void AliReconstruction::SetDefaultStorage(const char* uri) {
691 // Store the desired default CDB storage location
692 // Activate it later within the Run() method
698 //_____________________________________________________________________________
699 void AliReconstruction::SetQARefDefaultStorage(const char* uri) {
700 // Store the desired default CDB storage location
701 // Activate it later within the Run() method
704 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
707 //_____________________________________________________________________________
708 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
709 // Store a detector-specific CDB storage location
710 // Activate it later within the Run() method
712 AliCDBPath aPath(calibType);
713 if(!aPath.IsValid()){
714 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
715 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
716 if(!strcmp(calibType, fgkDetectorName[iDet])) {
717 aPath.SetPath(Form("%s/*", calibType));
718 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
722 if(!aPath.IsValid()){
723 AliError(Form("Not a valid path or detector: %s", calibType));
728 // // check that calibType refers to a "valid" detector name
729 // Bool_t isDetector = kFALSE;
730 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
731 // TString detName = fgkDetectorName[iDet];
732 // if(aPath.GetLevel0() == detName) {
733 // isDetector = kTRUE;
739 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
743 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
744 if (obj) fSpecCDBUri.Remove(obj);
745 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
749 //_____________________________________________________________________________
750 Bool_t AliReconstruction::SetRunNumberFromData()
752 // The method is called in Run() in order
753 // to set a correct run number.
754 // In case of raw data reconstruction the
755 // run number is taken from the raw data header
757 if (fSetRunNumberFromDataCalled) return kTRUE;
758 fSetRunNumberFromDataCalled = kTRUE;
760 AliCDBManager* man = AliCDBManager::Instance();
763 if(fRawReader->NextEvent()) {
764 if(man->GetRun() > 0) {
765 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
767 man->SetRun(fRawReader->GetRunNumber());
768 fRawReader->RewindEvents();
771 if(man->GetRun() > 0) {
772 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
775 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
781 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
783 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
788 // read run number from gAlice
789 if(rl->GetHeader()) {
790 man->SetRun(rl->GetHeader()->GetRun());
795 AliError("Neither run-loader header nor RawReader objects are found !");
807 //_____________________________________________________________________________
808 void AliReconstruction::SetCDBLock() {
809 // Set CDB lock: from now on it is forbidden to reset the run number
810 // or the default storage or to activate any further storage!
812 AliCDBManager::Instance()->SetLock(1);
815 //_____________________________________________________________________________
816 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
818 // Read the alignment objects from CDB.
819 // Each detector is supposed to have the
820 // alignment objects in DET/Align/Data CDB path.
821 // All the detector objects are then collected,
822 // sorted by geometry level (starting from ALIC) and
823 // then applied to the TGeo geometry.
824 // Finally an overlaps check is performed.
826 // Load alignment data from CDB and fill fAlignObjArray
827 if(fLoadAlignFromCDB){
829 TString detStr = detectors;
830 TString loadAlObjsListOfDets = "";
832 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
833 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
834 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
836 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
838 loadAlObjsListOfDets += fgkDetectorName[iDet];
839 loadAlObjsListOfDets += " ";
841 } // end loop over detectors
843 if(AliGeomManager::GetNalignable("GRP") != 0)
844 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
845 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
846 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
848 // Check if the array with alignment objects was
849 // provided by the user. If yes, apply the objects
850 // to the present TGeo geometry
851 if (fAlignObjArray) {
852 if (gGeoManager && gGeoManager->IsClosed()) {
853 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
854 AliError("The misalignment of one or more volumes failed!"
855 "Compare the list of simulated detectors and the list of detector alignment data!");
860 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
866 if (fAlignObjArray) {
867 fAlignObjArray->Delete();
868 delete fAlignObjArray; fAlignObjArray=NULL;
874 //_____________________________________________________________________________
875 void AliReconstruction::SetGAliceFile(const char* fileName)
877 // set the name of the galice file
879 fGAliceFileName = fileName;
882 //_____________________________________________________________________________
883 void AliReconstruction::SetInput(const char* input)
885 // In case the input string starts with 'mem://', we run in an online mode
886 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
887 // file is assumed. One can give as an input:
888 // mem://: - events taken from DAQ monitoring libs online
890 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
891 if (input) fRawInput = input;
894 //_____________________________________________________________________________
895 void AliReconstruction::SetOutput(const char* output)
897 // Set the output ESD filename
898 // 'output' is a normalt ROOT url
899 // The method is used in case of raw-data reco with PROOF
900 if (output) fESDOutput = output;
903 //_____________________________________________________________________________
904 void AliReconstruction::SetOption(const char* detector, const char* option)
906 // set options for the reconstruction of a detector
908 TObject* obj = fOptions.FindObject(detector);
909 if (obj) fOptions.Remove(obj);
910 fOptions.Add(new TNamed(detector, option));
913 //_____________________________________________________________________________
914 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
916 // Set custom reconstruction parameters for a given detector
917 // Single set of parameters for all the events
919 // First check if the reco-params are global
920 if(!strcmp(detector, "GRP")) {
922 fRecoParam.AddDetRecoParam(kNDetectors,par);
926 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
927 if(!strcmp(detector, fgkDetectorName[iDet])) {
929 fRecoParam.AddDetRecoParam(iDet,par);
936 //_____________________________________________________________________________
937 Bool_t AliReconstruction::InitGRP() {
938 //------------------------------------
939 // Initialization of the GRP entry
940 //------------------------------------
941 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
945 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
948 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
950 fGRPData = new AliGRPObject();
951 fGRPData->ReadValuesFromMap(m);
955 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
956 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
960 // FIX ME: The unloading of GRP entry is temporarily disabled
961 // because ZDC and VZERO are using it in order to initialize
962 // their reconstructor objects. In the future one has to think
963 // of propagating AliRunInfo to the reconstructors.
964 // AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
968 AliError("No GRP entry found in OCDB!");
972 TString lhcState = fGRPData->GetLHCState();
973 if (lhcState==AliGRPObject::GetInvalidString()) {
974 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
975 lhcState = "UNKNOWN";
978 TString beamType = fGRPData->GetBeamType();
979 if (beamType==AliGRPObject::GetInvalidString()) {
980 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
981 beamType = "UNKNOWN";
984 Float_t beamEnergy = fGRPData->GetBeamEnergy();
985 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
986 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
990 TString runType = fGRPData->GetRunType();
991 if (runType==AliGRPObject::GetInvalidString()) {
992 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
996 Int_t activeDetectors = fGRPData->GetDetectorMask();
997 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
998 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
999 activeDetectors = 1074790399;
1002 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1006 // Process the list of active detectors
1007 if (activeDetectors) {
1008 UInt_t detMask = activeDetectors;
1009 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1010 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1011 fFillESD = MatchDetectorList(fFillESD,detMask);
1012 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1013 fLoadCDB.Form("%s %s %s %s",
1014 fRunLocalReconstruction.Data(),
1015 fRunTracking.Data(),
1017 fQADetectors.Data());
1018 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1019 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1) &&
1020 !((detMask >> AliDAQ::DetectorID("ITSSDD")) & 0x1) &&
1021 !((detMask >> AliDAQ::DetectorID("ITSSSD")) & 0x1) ) {
1022 // switch off the vertexer
1023 AliInfo("SPD,SDD,SSD is not in the list of active detectors. Vertexer switched off.");
1024 fRunVertexFinder = kFALSE;
1026 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1027 // switch off the reading of CTP raw-data payload
1028 if (fFillTriggerESD) {
1029 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1030 fFillTriggerESD = kFALSE;
1035 AliInfo("===================================================================================");
1036 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1037 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1038 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1039 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1040 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1041 AliInfo("===================================================================================");
1043 //*** Dealing with the magnetic field map
1044 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {
1045 if (TGeoGlobalMagField::Instance()->GetField()->TestBit(AliMagF::kOverrideGRP)) {
1046 AliInfo("ExpertMode!!! GRP information will be ignored !");
1047 AliInfo("ExpertMode!!! Running with the externally locked B field !");
1050 AliInfo("Destroying existing B field instance!");
1051 delete TGeoGlobalMagField::Instance();
1054 if ( !TGeoGlobalMagField::Instance()->IsLocked() ) {
1055 // Construct the field map out of the information retrieved from GRP.
1058 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1059 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1060 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1064 Char_t l3Polarity = fGRPData->GetL3Polarity();
1065 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1066 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1071 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1072 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1073 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1077 Char_t diPolarity = fGRPData->GetDipolePolarity();
1078 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1079 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1083 // read special bits for the polarity convention and map type
1084 Int_t polConvention = fGRPData->IsPolarityConventionLHC() ? AliMagF::kConvLHC : AliMagF::kConvDCS2008;
1085 Bool_t uniformB = fGRPData->IsUniformBMap();
1088 AliMagF* fld = AliMagF::CreateFieldMap(TMath::Abs(l3Current) * (l3Polarity ? -1:1),
1089 TMath::Abs(diCurrent) * (diPolarity ? -1:1),
1090 polConvention,uniformB,beamEnergy, beamType.Data());
1092 TGeoGlobalMagField::Instance()->SetField( fld );
1093 TGeoGlobalMagField::Instance()->Lock();
1094 AliInfo("Running with the B field constructed out of GRP !");
1096 else AliFatal("Failed to create a B field map !");
1098 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1101 //*** Get the diamond profiles from OCDB
1102 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1104 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1106 AliError("No SPD diamond profile found in OCDB!");
1109 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1111 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1113 AliError("No diamond profile found in OCDB!");
1116 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1118 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1120 AliError("No TPC diamond profile found in OCDB!");
1123 entry = AliCDBManager::Instance()->Get("GRP/Calib/CosmicTriggers");
1125 fListOfCosmicTriggers = dynamic_cast<THashTable*>(entry->GetObject());
1127 AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
1130 if (!fListOfCosmicTriggers) {
1131 AliWarning("Can not get list of cosmic triggers from OCDB! Cosmic event specie will be effectively disabled!");
1137 //_____________________________________________________________________________
1138 Bool_t AliReconstruction::LoadCDB()
1140 AliCodeTimerAuto("",0);
1142 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1144 TString detStr = fLoadCDB;
1145 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1146 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1147 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1150 // Temporary fix - one has to define the correct policy in order
1151 // to load the trigger OCDB entries only for the detectors that
1152 // in the trigger or that are needed in order to put correct
1153 // information in ESD
1154 AliCDBManager::Instance()->GetAll("TRIGGER/*/*");
1158 //_____________________________________________________________________________
1159 Bool_t AliReconstruction::LoadTriggerScalersCDB()
1161 AliCodeTimerAuto("",0);
1163 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
1167 AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
1168 fRunScalers = dynamic_cast<AliTriggerRunScalers*> (entry->GetObject());
1170 if (fRunScalers->CorrectScalersOverflow() == 0) AliInfo("32bit Trigger counters corrected for overflow");
1175 //_____________________________________________________________________________
1176 Bool_t AliReconstruction::LoadCTPTimeParamsCDB()
1178 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
1182 AliInfo("Found an AliCTPTimeParams in GRP/CTP/CTPtiming, reading it");
1183 fCTPTimeParams = dynamic_cast<AliCTPTimeParams*> (entry->GetObject());
1190 //_____________________________________________________________________________
1191 Bool_t AliReconstruction::Run(const char* input)
1194 AliCodeTimerAuto("",0);
1197 if (GetAbort() != TSelector::kContinue) return kFALSE;
1199 TChain *chain = NULL;
1200 if (fRawReader && (chain = fRawReader->GetChain())) {
1201 Long64_t nEntries = (fLastEvent < 0) ? (TChain::kBigNumber) : (fLastEvent - fFirstEvent + 1);
1204 // Temporary fix for long raw-data runs (until socket timeout handling in PROOF is revised)
1205 gProof->Exec("gEnv->SetValue(\"Proof.SocketActivityTimeout\",-1)", kTRUE);
1208 gProof->Exec("TGrid::Connect(\"alien://\")",kTRUE);
1210 TMessage::EnableSchemaEvolutionForAll(kTRUE);
1211 gProof->Exec("TMessage::EnableSchemaEvolutionForAll(kTRUE)",kTRUE);
1213 gProof->AddInput(this);
1215 if (!ParseOutput()) return kFALSE;
1217 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1219 chain->Process("AliReconstruction","",nEntries,fFirstEvent);
1222 chain->Process(this,"",nEntries,fFirstEvent);
1227 if (GetAbort() != TSelector::kContinue) return kFALSE;
1229 if (GetAbort() != TSelector::kContinue) return kFALSE;
1230 //******* The loop over events
1231 AliInfo("Starting looping over events");
1233 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1234 (fRawReader && fRawReader->NextEvent())) {
1235 if (!ProcessEvent(iEvent)) {
1236 Abort("ProcessEvent",TSelector::kAbortFile);
1242 if (GetAbort() != TSelector::kContinue) return kFALSE;
1244 if (GetAbort() != TSelector::kContinue) return kFALSE;
1250 //_____________________________________________________________________________
1251 void AliReconstruction::InitRawReader(const char* input)
1253 AliCodeTimerAuto("",0);
1255 // Init raw-reader and
1256 // set the input in case of raw data
1257 if (input) fRawInput = input;
1258 fRawReader = AliRawReader::Create(fRawInput.Data());
1260 if (fRawInput.IsNull()) {
1261 AliInfo("Reconstruction will run over digits");
1264 AliFatal("Can not create raw-data reader ! Exiting...");
1268 if (!fEquipIdMap.IsNull() && fRawReader)
1269 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1271 if (!fUseHLTData.IsNull()) {
1272 // create the RawReaderHLT which performs redirection of HLT input data for
1273 // the specified detectors
1274 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1276 fParentRawReader=fRawReader;
1277 fRawReader=pRawReader;
1279 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1282 AliSysInfo::AddStamp("CreateRawReader");
1285 //_____________________________________________________________________________
1286 void AliReconstruction::InitRun(const char* input)
1288 // Initialization of raw-reader,
1289 // run number, CDB etc.
1290 AliCodeTimerAuto("",0);
1291 AliSysInfo::AddStamp("Start");
1293 // Initialize raw-reader if any
1294 InitRawReader(input);
1296 // Initialize the CDB storage
1299 // Set run number in CDBManager (if it is not already set by the user)
1300 if (!SetRunNumberFromData()) {
1301 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1305 // Set CDB lock: from now on it is forbidden to reset the run number
1306 // or the default storage or to activate any further storage!
1311 //_____________________________________________________________________________
1312 void AliReconstruction::Begin(TTree *)
1314 // Initialize AlReconstruction before
1315 // going into the event loop
1316 // Should follow the TSelector convention
1317 // i.e. initialize only the object on the client side
1318 AliCodeTimerAuto("",0);
1320 AliReconstruction *reco = NULL;
1322 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1325 AliSysInfo::AddStamp("ReadInputInBegin");
1328 // Import ideal TGeo geometry and apply misalignment
1330 TString geom(gSystem->DirName(fGAliceFileName));
1331 geom += "/geometry.root";
1332 AliGeomManager::LoadGeometry(geom.Data());
1334 Abort("LoadGeometry", TSelector::kAbortProcess);
1337 AliSysInfo::AddStamp("LoadGeom");
1338 TString detsToCheck=fRunLocalReconstruction;
1339 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1340 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1343 AliSysInfo::AddStamp("CheckGeom");
1346 if (!MisalignGeometry(fLoadAlignData)) {
1347 Abort("MisalignGeometry", TSelector::kAbortProcess);
1350 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1351 AliSysInfo::AddStamp("MisalignGeom");
1354 Abort("InitGRP", TSelector::kAbortProcess);
1357 AliSysInfo::AddStamp("InitGRP");
1360 Abort("LoadCDB", TSelector::kAbortProcess);
1363 AliSysInfo::AddStamp("LoadCDB");
1365 if (!LoadTriggerScalersCDB()) {
1366 Abort("LoadTriggerScalersCDB", TSelector::kAbortProcess);
1369 AliSysInfo::AddStamp("LoadTriggerScalersCDB");
1371 if (!LoadCTPTimeParamsCDB()) {
1372 Abort("LoadCTPTimeParamsCDB", TSelector::kAbortProcess);
1375 AliSysInfo::AddStamp("LoadCTPTimeParamsCDB");
1377 // Read the reconstruction parameters from OCDB
1378 if (!InitRecoParams()) {
1379 AliWarning("Not all detectors have correct RecoParam objects initialized");
1381 AliSysInfo::AddStamp("InitRecoParams");
1383 if (fInput && gProof) {
1384 if (reco) *reco = *this;
1386 gGeoManager->SetName("Geometry");
1387 gProof->AddInputData(gGeoManager,kTRUE);
1389 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1390 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1391 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1392 magFieldMap->SetName("MagneticFieldMap");
1393 gProof->AddInputData(magFieldMap,kTRUE);
1398 //_____________________________________________________________________________
1399 void AliReconstruction::SlaveBegin(TTree*)
1401 // Initialization related to run-loader,
1402 // vertexer, trackers, recontructors
1403 // In proof mode it is executed on the slave
1404 AliCodeTimerAuto("",0);
1406 TProofOutputFile *outProofFile = NULL;
1408 if (AliDebugLevel() > 0) fInput->Print();
1409 if (AliDebugLevel() > 10) fInput->Dump();
1410 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1413 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1415 AliGeomManager::SetGeometry(tgeo);
1417 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1418 Int_t runNumber = -1;
1419 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1420 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1421 man->SetCacheFlag(kTRUE);
1422 man->SetLock(kTRUE);
1426 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1427 AliMagF *newMap = new AliMagF(*map);
1428 if (!newMap->LoadParameterization()) {
1429 Abort("AliMagF::LoadParameterization", TSelector::kAbortProcess);
1432 TGeoGlobalMagField::Instance()->SetField(newMap);
1433 TGeoGlobalMagField::Instance()->Lock();
1435 if (TNamed *outputFileName = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE"))
1436 fProofOutputFileName = outputFileName->GetTitle();
1437 if (TNamed *outputLocation = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_LOCATION"))
1438 fProofOutputLocation = outputLocation->GetTitle();
1439 if (fInput->FindObject("PROOF_OUTPUTFILE_DATASET"))
1440 fProofOutputDataset = kTRUE;
1441 if (TNamed *archiveList = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_ARCHIVE"))
1442 fProofOutputArchive = archiveList->GetTitle();
1443 if (!fProofOutputFileName.IsNull() &&
1444 !fProofOutputLocation.IsNull() &&
1445 fProofOutputArchive.IsNull()) {
1446 if (!fProofOutputDataset) {
1447 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"M");
1448 outProofFile->SetOutputFileName(Form("%s%s",fProofOutputLocation.Data(),fProofOutputFileName.Data()));
1451 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"DROV",fProofOutputLocation.Data());
1453 if (AliDebugLevel() > 0) outProofFile->Dump();
1454 fOutput->Add(outProofFile);
1456 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1459 // get the run loader
1460 if (!InitRunLoader()) {
1461 Abort("InitRunLoader", TSelector::kAbortProcess);
1464 AliSysInfo::AddStamp("LoadLoader");
1466 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1469 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1470 Abort("CreateTrackers", TSelector::kAbortProcess);
1473 AliSysInfo::AddStamp("CreateTrackers");
1475 // create the ESD output file and tree
1476 if (!outProofFile) {
1477 ffile = TFile::Open("AliESDs.root", "RECREATE");
1478 ffile->SetCompressionLevel(2);
1479 if (!ffile->IsOpen()) {
1480 Abort("OpenESDFile", TSelector::kAbortProcess);
1485 AliInfo(Form("Opening output PROOF file: %s/%s",
1486 outProofFile->GetDir(), outProofFile->GetFileName()));
1487 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1488 Abort(Form("Problems opening output PROOF file: %s/%s",
1489 outProofFile->GetDir(), outProofFile->GetFileName()),
1490 TSelector::kAbortProcess);
1495 ftree = new TTree("esdTree", "Tree with ESD objects");
1496 fesd = new AliESDEvent();
1497 fesd->CreateStdContent();
1498 // add a so far non-std object to the ESD, this will
1499 // become part of the std content
1500 fesd->AddObject(new AliESDHLTDecision);
1502 fesd->WriteToTree(ftree);
1503 if (fWriteESDfriend) {
1505 // Since we add the branch manually we must
1506 // book and add it after WriteToTree
1507 // otherwise it is created twice,
1508 // once via writetotree and once here.
1509 // The case for AliESDfriend is now
1510 // caught also in AlIESDEvent::WriteToTree but
1511 // be careful when changing the name (AliESDfriend is not
1512 // a TNamed so we had to hardwire it)
1513 fesdf = new AliESDfriend();
1514 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1515 br->SetFile("AliESDfriends.root");
1516 fesd->AddObject(fesdf);
1518 ftree->GetUserInfo()->Add(fesd);
1520 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1521 fhltesd = new AliESDEvent();
1522 fhltesd->CreateStdContent();
1523 // read the ESD template from CDB
1524 // HLT is allowed to put non-std content to its ESD, the non-std
1525 // objects need to be created before invocation of WriteToTree in
1526 // order to create all branches. Initialization is done from an
1527 // ESD layout template in CDB
1528 AliCDBManager* man = AliCDBManager::Instance();
1529 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1530 AliCDBEntry* hltESDConfig=NULL;
1531 if (man->GetId(hltESDConfigPath)!=NULL &&
1532 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1533 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1535 // init all internal variables from the list of objects
1536 pESDLayout->GetStdContent();
1538 // copy content and create non-std objects
1539 *fhltesd=*pESDLayout;
1542 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1543 hltESDConfigPath.GetPath().Data()));
1547 fhltesd->WriteToTree(fhlttree);
1548 fhlttree->GetUserInfo()->Add(fhltesd);
1550 ProcInfo_t procInfo;
1551 gSystem->GetProcInfo(&procInfo);
1552 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1555 //Initialize the QA and start of cycle
1556 if (fRunQA || fRunGlobalQA)
1559 //Initialize the Plane Efficiency framework
1560 if (fRunPlaneEff && !InitPlaneEff()) {
1561 Abort("InitPlaneEff", TSelector::kAbortProcess);
1565 if (strcmp(gProgName,"alieve") == 0)
1566 fRunAliEVE = InitAliEVE();
1571 //_____________________________________________________________________________
1572 Bool_t AliReconstruction::Process(Long64_t entry)
1574 // run the reconstruction over a single entry
1575 // from the chain with raw data
1576 AliCodeTimerAuto("",0);
1578 TTree *currTree = fChain->GetTree();
1579 AliRawVEvent *event = NULL;
1580 currTree->SetBranchAddress("rawevent",&event);
1581 currTree->GetEntry(entry);
1582 fRawReader = new AliRawReaderRoot(event);
1583 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1591 //_____________________________________________________________________________
1592 void AliReconstruction::Init(TTree *tree)
1595 AliError("The input tree is not found!");
1601 //_____________________________________________________________________________
1602 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1604 // run the reconstruction over a single event
1605 // The event loop is steered in Run method
1608 static Long_t oldMres=0;
1609 static Long_t oldMvir=0;
1610 static Float_t oldCPU=0;
1612 AliCodeTimerAuto("",0);
1616 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1617 fRunLoader->SetEventNumber(iEvent);
1618 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1620 fRunLoader->TreeE()->Fill();
1621 if (fRawReader && fRawReader->UseAutoSaveESD())
1622 fRunLoader->TreeE()->AutoSave("SaveSelf");
1625 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1630 fRunLoader->GetEvent(iEvent);
1632 // Fill Event-info object
1634 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo,fListOfCosmicTriggers);
1635 AliInfo(Form("================================= Processing event %d of type %-10s ==================================", iEvent,fRecoParam.PrintEventSpecie()));
1637 // Set the reco-params
1639 TString detStr = fLoadCDB;
1640 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1641 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1642 AliReconstructor *reconstructor = GetReconstructor(iDet);
1643 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1644 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1645 reconstructor->SetRecoParam(par);
1646 reconstructor->GetPidSettings(&PID);
1647 reconstructor->SetEventInfo(&fEventInfo);
1649 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1650 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(par->GetEventSpecie())) ;
1654 const AliDetectorRecoParam *grppar = fRecoParam.GetDetRecoParam(kNDetectors);
1655 AliQAManager::QAManager()->SetRecoParam(AliQAv1::kGLOBAL, grppar) ;
1656 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(grppar->GetEventSpecie())) ;
1660 if (fRunQA && IsInTasks(AliQAv1::kRAWS)) {
1661 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1662 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1664 // local single event reconstruction
1665 if (!fRunLocalReconstruction.IsNull()) {
1666 TString detectors=fRunLocalReconstruction;
1667 // run HLT event reconstruction first
1668 // ;-( IsSelected changes the string
1669 if (IsSelected("HLT", detectors) &&
1670 !RunLocalEventReconstruction("HLT")) {
1671 if (fStopOnError) {CleanUp(); return kFALSE;}
1673 detectors=fRunLocalReconstruction;
1674 detectors.ReplaceAll("HLT", "");
1675 if (!RunLocalEventReconstruction(detectors)) {
1684 // fill Event header information from the RawEventHeader
1685 if (fRawReader){FillRawEventHeaderESD(fesd);}
1686 if (fRawReader){FillRawEventHeaderESD(fhltesd);}
1688 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1689 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1690 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1691 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1693 fesd->SetEventSpecie(fRecoParam.GetEventSpecie());
1694 fhltesd->SetEventSpecie(fRecoParam.GetEventSpecie());
1696 // Set magnetic field from the tracker
1697 fesd->SetMagneticField(AliTracker::GetBz());
1698 fhltesd->SetMagneticField(AliTracker::GetBz());
1700 ((AliESDRun*)fesd->GetESDRun())->SetBeamEnergyIsSqrtSHalfGeV();
1701 ((AliESDRun*)fhltesd->GetESDRun())->SetBeamEnergyIsSqrtSHalfGeV();
1703 AliMagF* fld = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1704 if (fld) { // set info needed for field initialization
1705 fesd->SetCurrentL3(fld->GetCurrentSol());
1706 fesd->SetCurrentDip(fld->GetCurrentDip());
1707 fesd->SetBeamEnergy(fld->GetBeamEnergy());
1708 fesd->SetBeamType(fld->GetBeamTypeText());
1709 fesd->SetUniformBMap(fld->IsUniform());
1710 fesd->SetBInfoStored();
1712 fhltesd->SetCurrentL3(fld->GetCurrentSol());
1713 fhltesd->SetCurrentDip(fld->GetCurrentDip());
1714 fhltesd->SetBeamEnergy(fld->GetBeamEnergy());
1715 fhltesd->SetBeamType(fld->GetBeamTypeText());
1716 fhltesd->SetUniformBMap(fld->IsUniform());
1717 fhltesd->SetBInfoStored();
1720 // Set most probable pt, for B=0 tracking
1721 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1722 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1723 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1725 // Fill raw-data error log into the ESD
1726 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1729 if (fRunVertexFinder) {
1730 if (!RunVertexFinder(fesd)) {
1731 if (fStopOnError) {CleanUp(); return kFALSE;}
1735 // For Plane Efficiency: run the SPD trackleter
1736 if (fRunPlaneEff && fSPDTrackleter) {
1737 if (!RunSPDTrackleting(fesd)) {
1738 if (fStopOnError) {CleanUp(); return kFALSE;}
1743 if (!fRunTracking.IsNull()) {
1744 if (fRunMuonTracking) {
1745 if (!RunMuonTracking(fesd)) {
1746 if (fStopOnError) {CleanUp(); return kFALSE;}
1752 if (!fRunTracking.IsNull()) {
1753 if (!RunTracking(fesd,PID)) {
1754 if (fStopOnError) {CleanUp(); return kFALSE;}
1759 if (!fFillESD.IsNull()) {
1760 TString detectors=fFillESD;
1761 // run HLT first and on hltesd
1762 // ;-( IsSelected changes the string
1763 if (IsSelected("HLT", detectors) &&
1764 !FillESD(fhltesd, "HLT")) {
1765 if (fStopOnError) {CleanUp(); return kFALSE;}
1768 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1769 if (detectors.Contains("ALL")) {
1771 for (Int_t idet=0; idet<kNDetectors; ++idet){
1772 detectors += fgkDetectorName[idet];
1776 detectors.ReplaceAll("HLT", "");
1777 if (!FillESD(fesd, detectors)) {
1778 if (fStopOnError) {CleanUp(); return kFALSE;}
1785 if (fFillTriggerESD) {
1786 if (!FillTriggerESD(fesd)) {
1787 if (fStopOnError) {CleanUp(); return kFALSE;}
1790 // Always fill scalers
1791 if (!FillTriggerScalers(fesd)) {
1792 if (fStopOnError) {CleanUp(); return kFALSE;}
1799 // Propagate track to the beam pipe (if not already done by ITS)
1801 const Int_t ntracks = fesd->GetNumberOfTracks();
1802 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1805 UShort_t *selectedIdx=new UShort_t[ntracks];
1807 for (Int_t itrack=0; itrack<ntracks; itrack++){
1808 const Double_t kMaxStep = 1; //max step over the material
1811 AliESDtrack *track = fesd->GetTrack(itrack);
1812 if (!track) continue;
1814 AliExternalTrackParam *tpcTrack =
1815 (AliExternalTrackParam *)track->GetTPCInnerParam();
1819 PropagateTrackToBxByBz(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1822 Int_t n=trkArray.GetEntriesFast();
1823 selectedIdx[n]=track->GetID();
1824 trkArray.AddLast(tpcTrack);
1827 //Tracks refitted by ITS should already be at the SPD vertex
1828 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1831 PropagateTrackToBxByBz(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1832 Double_t x[3]; track->GetXYZ(x);
1833 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1834 track->RelateToVertexBxByBz(fesd->GetPrimaryVertexSPD(), b, kVeryBig);
1839 // Improve the reconstructed primary vertex position using the tracks
1841 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1842 if(fesd->GetPrimaryVertexSPD()) {
1843 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1844 if(vtitle.Contains("cosmics")) {
1845 runVertexFinderTracks=kFALSE;
1849 if (runVertexFinderTracks) {
1850 // TPC + ITS primary vertex
1851 ftVertexer->SetITSMode();
1852 ftVertexer->SetConstraintOff();
1853 // get cuts for vertexer from AliGRPRecoParam
1855 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1856 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1857 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1858 ftVertexer->SetCuts(cutsVertexer);
1859 delete [] cutsVertexer; cutsVertexer = NULL;
1860 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1861 ftVertexer->SetVtxStart(fDiamondProfile);
1863 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1865 if (pvtx->GetStatus()) {
1866 fesd->SetPrimaryVertexTracks(pvtx);
1867 for (Int_t i=0; i<ntracks; i++) {
1868 AliESDtrack *t = fesd->GetTrack(i);
1869 Double_t x[3]; t->GetXYZ(x);
1870 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1871 t->RelateToVertexBxByBz(pvtx, b, kVeryBig);
1874 delete pvtx; pvtx=NULL;
1877 // TPC-only primary vertex
1878 ftVertexer->SetTPCMode();
1879 ftVertexer->SetConstraintOff();
1880 // get cuts for vertexer from AliGRPRecoParam
1882 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1883 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1884 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1885 ftVertexer->SetCuts(cutsVertexer);
1886 delete [] cutsVertexer; cutsVertexer = NULL;
1887 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1888 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1890 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1892 if (pvtx->GetStatus()) {
1893 fesd->SetPrimaryVertexTPC(pvtx);
1894 for (Int_t i=0; i<ntracks; i++) {
1895 AliESDtrack *t = fesd->GetTrack(i);
1896 Double_t x[3]; t->GetXYZ(x);
1897 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1898 t->RelateToVertexTPCBxByBz(pvtx, b, kVeryBig);
1901 delete pvtx; pvtx=NULL;
1905 delete[] selectedIdx;
1907 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1912 AliV0vertexer vtxer;
1913 vtxer.Tracks2V0vertices(fesd);
1915 if (fRunCascadeFinder) {
1917 AliCascadeVertexer cvtxer;
1918 cvtxer.V0sTracks2CascadeVertices(fesd);
1923 if (fCleanESD) CleanESD(fesd);
1925 if (fRunQA && IsInTasks(AliQAv1::kESDS)) {
1926 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1927 AliQAManager::QAManager()->RunOneEvent(fesd, fhltesd) ;
1930 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
1931 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1932 if (qadm && IsInTasks(AliQAv1::kESDS))
1933 qadm->Exec(AliQAv1::kESDS, fesd);
1936 // copy HLT decision from HLTesd to esd
1937 // the most relevant information is stored in a reduced container in the esd,
1938 // while the full information can be found in the HLTesd
1939 TObject* pHLTSrc=fhltesd->FindListObject(AliESDHLTDecision::Name());
1940 TObject* pHLTTgt=fesd->FindListObject(AliESDHLTDecision::Name());
1941 if (pHLTSrc && pHLTTgt) {
1942 pHLTSrc->Copy(*pHLTTgt);
1945 if (fWriteESDfriend) {
1946 // fesdf->~AliESDfriend();
1947 // new (fesdf) AliESDfriend(); // Reset...
1948 fesd->GetESDfriend(fesdf);
1952 // Auto-save the ESD tree in case of prompt reco @P2
1953 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1954 ftree->AutoSave("SaveSelf");
1955 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1956 if (friendfile) friendfile->Save();
1963 if (fRunAliEVE) RunAliEVE();
1967 if (fWriteESDfriend) {
1968 fesdf->~AliESDfriend();
1969 new (fesdf) AliESDfriend(); // Reset...
1972 ProcInfo_t procInfo;
1973 gSystem->GetProcInfo(&procInfo);
1974 AliInfo(Form("========================= End Event %d: memory res %d(%3d) vir %d(%3d) CPU %5.2f =====================",
1975 iEvent, procInfo.fMemResident/1024, (procInfo.fMemResident-oldMres)/1024,
1976 procInfo.fMemVirtual/1024,(procInfo.fMemVirtual-oldMvir)/1024,procInfo.fCpuUser+procInfo.fCpuSys-oldCPU));
1977 oldMres=procInfo.fMemResident;
1978 oldMvir=procInfo.fMemVirtual;
1979 oldCPU=procInfo.fCpuUser+procInfo.fCpuSys;
1982 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1983 if (fReconstructor[iDet]) {
1984 fReconstructor[iDet]->SetRecoParam(NULL);
1985 fReconstructor[iDet]->SetEventInfo(NULL);
1987 if (fTracker[iDet]) fTracker[iDet]->SetEventInfo(NULL);
1990 if (fRunQA || fRunGlobalQA)
1991 AliQAManager::QAManager()->Increment() ;
1996 //_____________________________________________________________________________
1997 void AliReconstruction::SlaveTerminate()
1999 // Finalize the run on the slave side
2000 // Called after the exit
2001 // from the event loop
2002 AliCodeTimerAuto("",0);
2004 if (fIsNewRunLoader) { // galice.root didn't exist
2005 fRunLoader->WriteHeader("OVERWRITE");
2006 fRunLoader->CdGAFile();
2007 fRunLoader->Write(0, TObject::kOverwrite);
2010 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
2011 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
2013 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
2014 cdbMapCopy->SetOwner(1);
2015 cdbMapCopy->SetName("cdbMap");
2016 TIter iter(cdbMap->GetTable());
2019 while((pair = dynamic_cast<TPair*> (iter.Next()))){
2020 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
2021 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
2022 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
2025 TList *cdbListCopy = new TList();
2026 cdbListCopy->SetOwner(1);
2027 cdbListCopy->SetName("cdbList");
2029 TIter iter2(cdbList);
2032 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
2033 cdbListCopy->Add(new TObjString(id->ToString().Data()));
2036 ftree->GetUserInfo()->Add(cdbMapCopy);
2037 ftree->GetUserInfo()->Add(cdbListCopy);
2042 if (fWriteESDfriend)
2043 ftree->SetBranchStatus("ESDfriend*",0);
2044 // we want to have only one tree version number
2045 ftree->Write(ftree->GetName(),TObject::kOverwrite);
2046 fhlttree->Write(fhlttree->GetName(),TObject::kOverwrite);
2048 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
2049 if (fRunPlaneEff && !FinishPlaneEff()) {
2050 AliWarning("Finish PlaneEff evaluation failed");
2053 // End of cycle for the in-loop
2055 if (fRunQA || fRunGlobalQA) {
2056 AliQAManager::QAManager()->EndOfCycle() ;
2058 !fProofOutputLocation.IsNull() &&
2059 fProofOutputArchive.IsNull() &&
2060 !fProofOutputDataset) {
2061 TString qaOutputFile(Form("%sMerged.%s.Data.root",
2062 fProofOutputLocation.Data(),
2063 AliQAv1::GetQADataFileName()));
2064 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",
2065 AliQAv1::GetQADataFileName()));
2066 qaProofFile->SetOutputFileName(qaOutputFile.Data());
2067 if (AliDebugLevel() > 0) qaProofFile->Dump();
2068 fOutput->Add(qaProofFile);
2069 MergeQA(qaProofFile->GetFileName());
2080 if (!fProofOutputFileName.IsNull() &&
2081 !fProofOutputLocation.IsNull() &&
2082 fProofOutputDataset &&
2083 !fProofOutputArchive.IsNull()) {
2084 TProofOutputFile *zipProofFile = new TProofOutputFile(fProofOutputFileName.Data(),
2086 fProofOutputLocation.Data());
2087 if (AliDebugLevel() > 0) zipProofFile->Dump();
2088 fOutput->Add(zipProofFile);
2089 TString fileList(fProofOutputArchive.Data());
2090 fileList.ReplaceAll(","," ");
2092 #if ROOT_SVN_REVISION >= 30174
2093 command.Form("zip -n root %s/%s %s",zipProofFile->GetDir(kTRUE),zipProofFile->GetFileName(),fileList.Data());
2095 command.Form("zip -n root %s/%s %s",zipProofFile->GetDir(),zipProofFile->GetFileName(),fileList.Data());
2097 AliInfo(Form("Executing: %s",command.Data()));
2098 gSystem->Exec(command.Data());
2103 //_____________________________________________________________________________
2104 void AliReconstruction::Terminate()
2106 // Create tags for the events in the ESD tree (the ESD tree is always present)
2107 // In case of empty events the tags will contain dummy values
2108 AliCodeTimerAuto("",0);
2110 // Do not call the ESD tag creator in case of PROOF-based reconstruction
2112 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
2113 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
2114 delete esdtagCreator;
2117 // Cleanup of CDB manager: cache and active storages!
2118 AliCDBManager::Instance()->ClearCache();
2121 //_____________________________________________________________________________
2122 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
2124 // run the local reconstruction
2126 static Int_t eventNr=0;
2127 AliCodeTimerAuto("",0)
2129 TString detStr = detectors;
2130 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2131 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2132 AliReconstructor* reconstructor = GetReconstructor(iDet);
2133 if (!reconstructor) continue;
2134 AliLoader* loader = fLoader[iDet];
2135 // Matthias April 2008: temporary fix to run HLT reconstruction
2136 // although the HLT loader is missing
2137 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2139 reconstructor->Reconstruct(fRawReader, NULL);
2142 reconstructor->Reconstruct(dummy, NULL);
2147 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2150 // conversion of digits
2151 if (fRawReader && reconstructor->HasDigitConversion()) {
2152 AliInfo(Form("converting raw data digits into root objects for %s",
2153 fgkDetectorName[iDet]));
2154 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2155 // fgkDetectorName[iDet]),0);
2156 loader->LoadDigits("update");
2157 loader->CleanDigits();
2158 loader->MakeDigitsContainer();
2159 TTree* digitsTree = loader->TreeD();
2160 reconstructor->ConvertDigits(fRawReader, digitsTree);
2161 loader->WriteDigits("OVERWRITE");
2162 loader->UnloadDigits();
2164 // local reconstruction
2165 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2166 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]),0);
2167 loader->LoadRecPoints("update");
2168 loader->CleanRecPoints();
2169 loader->MakeRecPointsContainer();
2170 TTree* clustersTree = loader->TreeR();
2171 if (fRawReader && !reconstructor->HasDigitConversion()) {
2172 reconstructor->Reconstruct(fRawReader, clustersTree);
2174 loader->LoadDigits("read");
2175 TTree* digitsTree = loader->TreeD();
2177 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2181 reconstructor->Reconstruct(digitsTree, clustersTree);
2182 if (fRunQA && IsInTasks(AliQAv1::kDIGITSR)) {
2183 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2184 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2187 loader->UnloadDigits();
2189 if (fRunQA && IsInTasks(AliQAv1::kRECPOINTS)) {
2190 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2191 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2193 loader->WriteRecPoints("OVERWRITE");
2194 loader->UnloadRecPoints();
2195 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2197 IsSelected("CTP", detStr);
2198 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2199 AliError(Form("the following detectors were not found: %s",
2207 //_____________________________________________________________________________
2208 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2210 // run the SPD trackleting (for SPD efficiency purpouses)
2212 AliCodeTimerAuto("",0)
2214 Double_t vtxPos[3] = {0, 0, 0};
2215 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2217 TArrayF mcVertex(3);
2219 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2220 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2221 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2224 const AliESDVertex *vertex = esd->GetVertex();
2226 AliWarning("Vertex not found");
2229 vertex->GetXYZ(vtxPos);
2230 vertex->GetSigmaXYZ(vtxErr);
2231 if (fSPDTrackleter) {
2232 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2235 fLoader[0]->LoadRecPoints("read");
2236 TTree* tree = fLoader[0]->TreeR();
2238 AliError("Can't get the ITS cluster tree");
2241 fSPDTrackleter->LoadClusters(tree);
2242 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2244 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2245 AliWarning("AliITSTrackleterSPDEff Clusters2Tracks failed");
2246 // fLoader[0]->UnloadRecPoints();
2249 //fSPDTrackleter->UnloadRecPoints();
2251 AliWarning("SPDTrackleter not available");
2257 //_____________________________________________________________________________
2258 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2260 // run the barrel tracking
2262 AliCodeTimerAuto("",0)
2264 AliVertexer *vertexer = CreateVertexer();
2265 if (!vertexer) return kFALSE;
2267 AliInfo("running the ITS vertex finder");
2268 AliESDVertex* vertex = NULL;
2270 fLoader[0]->LoadRecPoints();
2271 TTree* cltree = fLoader[0]->TreeR();
2273 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2274 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2277 AliError("Can't get the ITS cluster tree");
2279 fLoader[0]->UnloadRecPoints();
2282 AliError("Can't get the ITS loader");
2285 AliWarning("Vertex not found");
2286 vertex = new AliESDVertex();
2287 vertex->SetName("default");
2290 vertex->SetName("reconstructed");
2295 vertex->GetXYZ(vtxPos);
2296 vertex->GetSigmaXYZ(vtxErr);
2298 esd->SetPrimaryVertexSPD(vertex);
2299 AliESDVertex *vpileup = NULL;
2300 Int_t novertices = 0;
2301 vpileup = vertexer->GetAllVertices(novertices);
2303 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2305 // if SPD multiplicity has been determined, it is stored in the ESD
2306 AliMultiplicity *mult = vertexer->GetMultiplicity();
2307 if(mult)esd->SetMultiplicity(mult);
2309 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2310 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2319 //_____________________________________________________________________________
2320 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2322 // run the HLT barrel tracking
2324 AliCodeTimerAuto("",0)
2327 AliError("Missing runLoader!");
2331 AliInfo("running HLT tracking");
2333 // Get a pointer to the HLT reconstructor
2334 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2335 if (!reconstructor) return kFALSE;
2338 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2339 TString detName = fgkDetectorName[iDet];
2340 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2341 reconstructor->SetOption(detName.Data());
2342 AliTracker *tracker = reconstructor->CreateTracker();
2344 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2345 if (fStopOnError) return kFALSE;
2349 Double_t vtxErr[3]={0.005,0.005,0.010};
2350 const AliESDVertex *vertex = esd->GetVertex();
2351 vertex->GetXYZ(vtxPos);
2352 tracker->SetVertex(vtxPos,vtxErr);
2354 fLoader[iDet]->LoadRecPoints("read");
2355 TTree* tree = fLoader[iDet]->TreeR();
2357 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2360 tracker->LoadClusters(tree);
2362 if (tracker->Clusters2Tracks(esd) != 0) {
2363 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2367 tracker->UnloadClusters();
2375 //_____________________________________________________________________________
2376 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2378 // run the muon spectrometer tracking
2380 AliCodeTimerAuto("",0)
2383 AliError("Missing runLoader!");
2386 Int_t iDet = 7; // for MUON
2388 AliInfo("is running...");
2390 // Get a pointer to the MUON reconstructor
2391 AliReconstructor *reconstructor = GetReconstructor(iDet);
2392 if (!reconstructor) return kFALSE;
2395 TString detName = fgkDetectorName[iDet];
2396 AliDebug(1, Form("%s tracking", detName.Data()));
2397 AliTracker *tracker = reconstructor->CreateTracker();
2399 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2404 fLoader[iDet]->LoadRecPoints("read");
2406 tracker->LoadClusters(fLoader[iDet]->TreeR());
2408 Int_t rv = tracker->Clusters2Tracks(esd);
2410 fLoader[iDet]->UnloadRecPoints();
2412 tracker->UnloadClusters();
2418 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2426 //_____________________________________________________________________________
2427 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd,AliESDpid &PID)
2429 // run the barrel tracking
2430 static Int_t eventNr=0;
2431 AliCodeTimerAuto("",0)
2433 AliInfo("running tracking");
2435 // Set the event info which is used
2436 // by the trackers in order to obtain
2437 // information about read-out detectors,
2439 AliDebug(1, "Setting event info");
2440 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2441 if (!fTracker[iDet]) continue;
2442 fTracker[iDet]->SetEventInfo(&fEventInfo);
2445 //Fill the ESD with the T0 info (will be used by the TOF)
2446 if (fReconstructor[11] && fLoader[11]) {
2447 fLoader[11]->LoadRecPoints("READ");
2448 TTree *treeR = fLoader[11]->TreeR();
2450 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2454 // pass 1: TPC + ITS inwards
2455 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2456 if (!fTracker[iDet]) continue;
2457 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2460 fLoader[iDet]->LoadRecPoints("read");
2461 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2462 TTree* tree = fLoader[iDet]->TreeR();
2464 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2467 fTracker[iDet]->LoadClusters(tree);
2468 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2470 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2471 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2474 // preliminary PID in TPC needed by the ITS tracker
2476 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2477 PID.MakePID(esd,kTRUE);
2479 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2482 // pass 2: ALL backwards
2484 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2485 if (!fTracker[iDet]) continue;
2486 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2489 if (iDet > 1) { // all except ITS, TPC
2491 fLoader[iDet]->LoadRecPoints("read");
2492 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2493 tree = fLoader[iDet]->TreeR();
2495 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2498 fTracker[iDet]->LoadClusters(tree);
2499 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2503 if (iDet>1) // start filling residuals for the "outer" detectors
2505 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2506 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2508 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2509 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2510 if ( elem && (! elem->At(0)) ) {
2511 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2512 if (qadm) qadm->InitRecPointsForTracker() ;
2516 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2517 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2522 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2523 fTracker[iDet]->UnloadClusters();
2524 fLoader[iDet]->UnloadRecPoints();
2526 // updated PID in TPC needed by the ITS tracker -MI
2528 //GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2529 //AliESDpid::MakePID(esd);
2530 PID.MakePID(esd,kTRUE);
2532 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2534 //stop filling residuals for the "outer" detectors
2535 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2537 // pass 3: TRD + TPC + ITS refit inwards
2539 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2540 if (!fTracker[iDet]) continue;
2541 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2544 if (iDet<2) // start filling residuals for TPC and ITS
2546 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2547 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2549 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2550 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2551 if ( elem && (! elem->At(0)) ) {
2552 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2553 if (qadm) qadm->InitRecPointsForTracker() ;
2558 if (fTracker[iDet]->RefitInward(esd) != 0) {
2559 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2562 // run postprocessing
2563 if (fTracker[iDet]->PostProcess(esd) != 0) {
2564 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2567 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2570 // write space-points to the ESD in case alignment data output
2572 if (fWriteAlignmentData)
2573 WriteAlignmentData(esd);
2575 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2576 if (!fTracker[iDet]) continue;
2578 fTracker[iDet]->UnloadClusters();
2579 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2580 fLoader[iDet]->UnloadRecPoints();
2581 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2583 // stop filling residuals for TPC and ITS
2584 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2590 //_____________________________________________________________________________
2591 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2593 // Remove the data which are not needed for the physics analysis.
2596 Int_t nTracks=esd->GetNumberOfTracks();
2597 Int_t nV0s=esd->GetNumberOfV0s();
2599 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2601 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2602 Bool_t rc=esd->Clean(cleanPars);
2604 nTracks=esd->GetNumberOfTracks();
2605 nV0s=esd->GetNumberOfV0s();
2607 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2612 //_____________________________________________________________________________
2613 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2615 // fill the event summary data
2617 AliCodeTimerAuto("",0)
2618 static Int_t eventNr=0;
2619 TString detStr = detectors;
2621 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2622 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2623 AliReconstructor* reconstructor = GetReconstructor(iDet);
2624 if (!reconstructor) continue;
2625 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2626 TTree* clustersTree = NULL;
2627 if (fLoader[iDet]) {
2628 fLoader[iDet]->LoadRecPoints("read");
2629 clustersTree = fLoader[iDet]->TreeR();
2630 if (!clustersTree) {
2631 AliError(Form("Can't get the %s clusters tree",
2632 fgkDetectorName[iDet]));
2633 if (fStopOnError) return kFALSE;
2636 if (fRawReader && !reconstructor->HasDigitConversion()) {
2637 reconstructor->FillESD(fRawReader, clustersTree, esd);
2639 TTree* digitsTree = NULL;
2640 if (fLoader[iDet]) {
2641 fLoader[iDet]->LoadDigits("read");
2642 digitsTree = fLoader[iDet]->TreeD();
2644 AliError(Form("Can't get the %s digits tree",
2645 fgkDetectorName[iDet]));
2646 if (fStopOnError) return kFALSE;
2649 reconstructor->FillESD(digitsTree, clustersTree, esd);
2650 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2652 if (fLoader[iDet]) {
2653 fLoader[iDet]->UnloadRecPoints();
2657 IsSelected("CTP", detStr);
2658 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2659 AliError(Form("the following detectors were not found: %s",
2661 if (fStopOnError) return kFALSE;
2663 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2668 //_____________________________________________________________________________
2669 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2671 // Reads the trigger decision which is
2672 // stored in Trigger.root file and fills
2673 // the corresponding esd entries
2675 AliCodeTimerAuto("",0)
2677 AliInfo("Filling trigger information into the ESD");
2680 AliCTPRawStream input(fRawReader);
2681 if (!input.Next()) {
2682 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2685 if (esd->GetTriggerMask() != input.GetClassMask())
2686 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2687 input.GetClassMask(),esd->GetTriggerMask()));
2688 if (esd->GetOrbitNumber() != input.GetOrbitID())
2689 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2690 input.GetOrbitID(),esd->GetOrbitNumber()));
2691 if (esd->GetBunchCrossNumber() != input.GetBCID())
2692 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2693 input.GetBCID(),esd->GetBunchCrossNumber()));
2694 AliESDHeader* esdheader = esd->GetHeader();
2695 esdheader->SetL0TriggerInputs(input.GetL0Inputs());
2696 esdheader->SetL1TriggerInputs(input.GetL1Inputs());
2697 esdheader->SetL2TriggerInputs(input.GetL2Inputs());
2699 UInt_t orbit=input.GetOrbitID();
2700 for(Int_t i=0 ; i<input.GetNIRs() ; i++ )
2701 if(TMath::Abs(Int_t(orbit-(input.GetIR(i))->GetOrbit()))<=1){
2702 esdheader->AddTriggerIR(input.GetIR(i));
2708 //_____________________________________________________________________________
2709 Bool_t AliReconstruction::FillTriggerScalers(AliESDEvent*& esd)
2712 //fRunScalers->Print();
2713 if(fRunScalers && fRunScalers->CheckRunScalers()){
2714 AliTimeStamp* timestamp = new AliTimeStamp(esd->GetOrbitNumber(), esd->GetPeriodNumber(), esd->GetBunchCrossNumber());
2715 //AliTimeStamp* timestamp = new AliTimeStamp(10308000, 0, (ULong64_t)486238);
2716 AliESDHeader* esdheader = fesd->GetHeader();
2717 for(Int_t i=0;i<50;i++){
2718 if((1ull<<i) & esd->GetTriggerMask()){
2719 AliTriggerScalersESD* scalesd = fRunScalers->GetScalersForEventClass( timestamp, i+1);
2720 if(scalesd)esdheader->SetTriggerScalersRecord(scalesd);
2726 //_____________________________________________________________________________
2727 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2730 // Filling information from RawReader Header
2733 if (!fRawReader) return kFALSE;
2735 AliInfo("Filling information from RawReader Header");
2737 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2738 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2739 esd->SetPeriodNumber(fRawReader->GetPeriod());
2741 esd->SetTimeStamp(fRawReader->GetTimestamp());
2742 esd->SetEventType(fRawReader->GetType());
2748 //_____________________________________________________________________________
2749 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2751 // check whether detName is contained in detectors
2752 // if yes, it is removed from detectors
2754 // check if all detectors are selected
2755 if ((detectors.CompareTo("ALL") == 0) ||
2756 detectors.BeginsWith("ALL ") ||
2757 detectors.EndsWith(" ALL") ||
2758 detectors.Contains(" ALL ")) {
2763 // search for the given detector
2764 Bool_t result = kFALSE;
2765 if ((detectors.CompareTo(detName) == 0) ||
2766 detectors.BeginsWith(detName+" ") ||
2767 detectors.EndsWith(" "+detName) ||
2768 detectors.Contains(" "+detName+" ")) {
2769 detectors.ReplaceAll(detName, "");
2773 // clean up the detectors string
2774 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2775 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2776 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2781 //_____________________________________________________________________________
2782 Bool_t AliReconstruction::InitRunLoader()
2784 // get or create the run loader
2786 if (gAlice) delete gAlice;
2789 TFile *gafile = TFile::Open(fGAliceFileName.Data());
2790 // if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2791 if (gafile) { // galice.root exists
2795 // load all base libraries to get the loader classes
2796 TString libs = gSystem->GetLibraries();
2797 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2798 TString detName = fgkDetectorName[iDet];
2799 if (detName == "HLT") continue;
2800 if (libs.Contains("lib" + detName + "base.so")) continue;
2801 gSystem->Load("lib" + detName + "base.so");
2803 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2805 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2810 fRunLoader->CdGAFile();
2811 fRunLoader->LoadgAlice();
2813 //PH This is a temporary fix to give access to the kinematics
2814 //PH that is needed for the labels of ITS clusters
2815 fRunLoader->LoadHeader();
2816 fRunLoader->LoadKinematics();
2818 } else { // galice.root does not exist
2820 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2822 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2823 AliConfig::GetDefaultEventFolderName(),
2826 AliError(Form("could not create run loader in file %s",
2827 fGAliceFileName.Data()));
2831 fIsNewRunLoader = kTRUE;
2832 fRunLoader->MakeTree("E");
2834 if (fNumberOfEventsPerFile > 0)
2835 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2837 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2843 //_____________________________________________________________________________
2844 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2846 // get the reconstructor object and the loader for a detector
2848 if (fReconstructor[iDet]) {
2849 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2850 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2851 fReconstructor[iDet]->SetRecoParam(par);
2852 fReconstructor[iDet]->SetRunInfo(fRunInfo);
2854 return fReconstructor[iDet];
2857 // load the reconstructor object
2858 TPluginManager* pluginManager = gROOT->GetPluginManager();
2859 TString detName = fgkDetectorName[iDet];
2860 TString recName = "Ali" + detName + "Reconstructor";
2862 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2864 AliReconstructor* reconstructor = NULL;
2865 // first check if a plugin is defined for the reconstructor
2866 TPluginHandler* pluginHandler =
2867 pluginManager->FindHandler("AliReconstructor", detName);
2868 // if not, add a plugin for it
2869 if (!pluginHandler) {
2870 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2871 TString libs = gSystem->GetLibraries();
2872 if (libs.Contains("lib" + detName + "base.so") ||
2873 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2874 pluginManager->AddHandler("AliReconstructor", detName,
2875 recName, detName + "rec", recName + "()");
2877 pluginManager->AddHandler("AliReconstructor", detName,
2878 recName, detName, recName + "()");
2880 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2882 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2883 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2885 if (reconstructor) {
2886 TObject* obj = fOptions.FindObject(detName.Data());
2887 if (obj) reconstructor->SetOption(obj->GetTitle());
2888 reconstructor->SetRunInfo(fRunInfo);
2889 reconstructor->Init();
2890 fReconstructor[iDet] = reconstructor;
2893 // get or create the loader
2894 if (detName != "HLT") {
2895 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2896 if (!fLoader[iDet]) {
2897 AliConfig::Instance()
2898 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2900 // first check if a plugin is defined for the loader
2902 pluginManager->FindHandler("AliLoader", detName);
2903 // if not, add a plugin for it
2904 if (!pluginHandler) {
2905 TString loaderName = "Ali" + detName + "Loader";
2906 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2907 pluginManager->AddHandler("AliLoader", detName,
2908 loaderName, detName + "base",
2909 loaderName + "(const char*, TFolder*)");
2910 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2912 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2914 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2915 fRunLoader->GetEventFolder());
2917 if (!fLoader[iDet]) { // use default loader
2918 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2920 if (!fLoader[iDet]) {
2921 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2922 if (fStopOnError) return NULL;
2924 fRunLoader->AddLoader(fLoader[iDet]);
2925 fRunLoader->CdGAFile();
2926 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2927 fRunLoader->Write(0, TObject::kOverwrite);
2932 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2933 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2934 reconstructor->SetRecoParam(par);
2935 reconstructor->SetRunInfo(fRunInfo);
2937 return reconstructor;
2940 //_____________________________________________________________________________
2941 AliVertexer* AliReconstruction::CreateVertexer()
2943 // create the vertexer
2944 // Please note that the caller is the owner of the
2947 AliVertexer* vertexer = NULL;
2948 AliReconstructor* itsReconstructor = GetReconstructor(0);
2949 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) || fRunTracking.Contains("ITS"))) {
2950 vertexer = itsReconstructor->CreateVertexer();
2953 AliWarning("couldn't create a vertexer for ITS");
2959 //_____________________________________________________________________________
2960 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2962 // create the trackers
2963 AliInfo("Creating trackers");
2965 TString detStr = detectors;
2966 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2967 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2968 AliReconstructor* reconstructor = GetReconstructor(iDet);
2969 if (!reconstructor) continue;
2970 TString detName = fgkDetectorName[iDet];
2971 if (detName == "HLT") {
2972 fRunHLTTracking = kTRUE;
2975 if (detName == "MUON") {
2976 fRunMuonTracking = kTRUE;
2981 fTracker[iDet] = reconstructor->CreateTracker();
2982 if (!fTracker[iDet] && (iDet < 7)) {
2983 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2984 if (fStopOnError) return kFALSE;
2986 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2992 //_____________________________________________________________________________
2993 void AliReconstruction::CleanUp()
2995 // delete trackers and the run loader and close and delete the file
2997 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2998 delete fReconstructor[iDet];
2999 fReconstructor[iDet] = NULL;
3000 fLoader[iDet] = NULL;
3001 delete fTracker[iDet];
3002 fTracker[iDet] = NULL;
3007 delete fSPDTrackleter;
3008 fSPDTrackleter = NULL;
3017 delete fParentRawReader;
3018 fParentRawReader=NULL;
3026 if (AliQAManager::QAManager())
3027 AliQAManager::QAManager()->ShowQA() ;
3028 AliQAManager::Destroy() ;
3032 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
3034 // Write space-points which are then used in the alignment procedures
3035 // For the moment only ITS, TPC, TRD and TOF
3037 Int_t ntracks = esd->GetNumberOfTracks();
3038 for (Int_t itrack = 0; itrack < ntracks; itrack++)
3040 AliESDtrack *track = esd->GetTrack(itrack);
3043 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
3044 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
3045 nsp += (iDet==GetDetIndex("TRD")) ? track->GetTRDntracklets():track->GetNcls(iDet);
3047 if (iDet==GetDetIndex("ITS")) { // ITS "extra" clusters
3048 track->GetClusters(iDet,idx);
3049 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
3054 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
3055 track->SetTrackPointArray(sp);
3057 for (Int_t iDet = 5; iDet >= 0; iDet--) {
3058 AliTracker *tracker = fTracker[iDet];
3059 if (!tracker) continue;
3060 Int_t nspdet = (iDet==GetDetIndex("TRD")) ? track->GetTRDtracklets(idx):track->GetClusters(iDet,idx);
3062 if (iDet==GetDetIndex("ITS")) // ITS "extra" clusters
3063 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
3065 if (nspdet <= 0) continue;
3069 while (isp2 < nspdet) {
3070 Bool_t isvalid=kTRUE;
3072 Int_t index=idx[isp++];
3073 if (index < 0) continue;
3075 TString dets = fgkDetectorName[iDet];
3076 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
3077 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
3078 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
3079 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
3080 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
3082 isvalid = tracker->GetTrackPoint(index,p);
3085 if (!isvalid) continue;
3086 if (iDet==GetDetIndex("ITS") && (isp-1)>=6) p.SetExtra();
3087 sp->AddPoint(isptrack,&p); isptrack++;
3094 //_____________________________________________________________________________
3095 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
3097 // The method reads the raw-data error log
3098 // accumulated within the rawReader.
3099 // It extracts the raw-data errors related to
3100 // the current event and stores them into
3101 // a TClonesArray inside the esd object.
3103 if (!fRawReader) return;
3105 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
3107 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
3109 if (iEvent != log->GetEventNumber()) continue;
3111 esd->AddRawDataErrorLog(log);
3116 //_____________________________________________________________________________
3117 void AliReconstruction::CheckQA()
3119 // check the QA of SIM for this run and remove the detectors
3120 // with status Fatal
3122 // TString newRunLocalReconstruction ;
3123 // TString newRunTracking ;
3124 // TString newFillESD ;
3126 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
3127 // TString detName(AliQAv1::GetDetName(iDet)) ;
3128 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
3129 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
3130 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
3131 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
3133 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
3134 // fRunLocalReconstruction.Contains("ALL") ) {
3135 // newRunLocalReconstruction += detName ;
3136 // newRunLocalReconstruction += " " ;
3138 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
3139 // fRunTracking.Contains("ALL") ) {
3140 // newRunTracking += detName ;
3141 // newRunTracking += " " ;
3143 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
3144 // fFillESD.Contains("ALL") ) {
3145 // newFillESD += detName ;
3146 // newFillESD += " " ;
3150 // fRunLocalReconstruction = newRunLocalReconstruction ;
3151 // fRunTracking = newRunTracking ;
3152 // fFillESD = newFillESD ;
3155 //_____________________________________________________________________________
3156 Int_t AliReconstruction::GetDetIndex(const char* detector)
3158 // return the detector index corresponding to detector
3160 for (index = 0; index < kNDetectors ; index++) {
3161 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3166 //_____________________________________________________________________________
3167 Bool_t AliReconstruction::FinishPlaneEff() {
3169 // Here execute all the necessary operationis, at the end of the tracking phase,
3170 // in case that evaluation of PlaneEfficiencies was required for some detector.
3171 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3173 // This Preliminary version works only FOR ITS !!!!!
3174 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3177 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3180 TString detStr = fLoadCDB;
3181 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3182 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3183 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3184 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
3185 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3186 TString name=planeeff->GetName();
3188 TFile* pefile = TFile::Open(name, "RECREATE");
3189 ret=(Bool_t)planeeff->Write();
3191 if(planeeff->GetCreateHistos()) {
3192 TString hname=planeeff->GetName();
3193 hname+="Histo.root";
3194 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3197 if(fSPDTrackleter) {
3198 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3199 TString name="AliITSPlaneEffSPDtracklet.root";
3200 TFile* pefile = TFile::Open(name, "RECREATE");
3201 ret=(Bool_t)planeeff->Write();
3203 AliESDEvent *dummy=NULL;
3204 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3209 //_____________________________________________________________________________
3210 Bool_t AliReconstruction::InitPlaneEff() {
3212 // Here execute all the necessary operations, before of the tracking phase,
3213 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3214 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3215 // which should be updated/recalculated.
3217 // This Preliminary version will work only FOR ITS !!!!!
3218 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3221 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3223 AliWarning(Form("Implementation of this method not yet completed !! Method return kTRUE"));
3225 fSPDTrackleter = NULL;
3226 TString detStr = fLoadCDB;
3227 if (IsSelected(fgkDetectorName[0], detStr)) {
3228 AliReconstructor* itsReconstructor = GetReconstructor(0);
3229 if (itsReconstructor) {
3230 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3232 if (fSPDTrackleter) {
3233 AliInfo("Trackleter for SPD has been created");
3239 //_____________________________________________________________________________
3240 Bool_t AliReconstruction::InitAliEVE()
3242 // This method should be called only in case
3243 // AliReconstruction is run
3244 // within the alieve environment.
3245 // It will initialize AliEVE in a way
3246 // so that it can visualize event processed
3247 // by AliReconstruction.
3248 // The return flag shows whenever the
3249 // AliEVE initialization was successful or not.
3252 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3253 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3254 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3256 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3257 gROOT->ProcessLine("alieve_online_init()");
3262 //_____________________________________________________________________________
3263 void AliReconstruction::RunAliEVE()
3265 // Runs AliEVE visualisation of
3266 // the current event.
3267 // Should be executed only after
3268 // successful initialization of AliEVE.
3270 AliInfo("Running AliEVE...");
3271 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3275 //_____________________________________________________________________________
3276 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3278 // Allows to run QA for a selected set of detectors
3279 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3280 // all selected detectors run the same selected tasks
3282 if (!detAndAction.Contains(":")) {
3283 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3287 Int_t colon = detAndAction.Index(":") ;
3288 fQADetectors = detAndAction(0, colon) ;
3289 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3290 if (fQATasks.Contains("ALL") ) {
3291 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3293 fQATasks.ToUpper() ;
3295 if ( fQATasks.Contains("RAW") )
3296 tempo = Form("%d ", AliQAv1::kRAWS) ;
3297 if ( fQATasks.Contains("DIGIT") )
3298 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3299 if ( fQATasks.Contains("RECPOINT") )
3300 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3301 if ( fQATasks.Contains("ESD") )
3302 tempo += Form("%d ", AliQAv1::kESDS) ;
3304 if (fQATasks.IsNull()) {
3305 AliInfo("No QA requested\n") ;
3310 TString tempo(fQATasks) ;
3311 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3312 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3313 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3314 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3315 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3320 //_____________________________________________________________________________
3321 Bool_t AliReconstruction::InitRecoParams()
3323 // The method accesses OCDB and retrieves all
3324 // the available reco-param objects from there.
3326 Bool_t isOK = kTRUE;
3328 if (fRecoParam.GetDetRecoParamArray(kNDetectors)) {
3329 AliInfo("Using custom GRP reconstruction parameters");
3332 AliInfo("Loading GRP reconstruction parameter objects");
3334 AliCDBPath path("GRP","Calib","RecoParam");
3335 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3337 AliWarning("Couldn't find GRP RecoParam entry in OCDB");
3341 TObject *recoParamObj = entry->GetObject();
3342 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3343 // GRP has a normal TobjArray of AliDetectorRecoParam objects
3344 // Registering them in AliRecoParam
3345 fRecoParam.AddDetRecoParamArray(kNDetectors,dynamic_cast<TObjArray*>(recoParamObj));
3347 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3348 // GRP has only onse set of reco parameters
3349 // Registering it in AliRecoParam
3350 AliInfo("Single set of GRP reconstruction parameters found");
3351 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3352 fRecoParam.AddDetRecoParam(kNDetectors,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3355 AliError("No valid GRP RecoParam object found in the OCDB");
3362 TString detStr = fLoadCDB;
3363 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3365 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3367 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3368 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3372 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3374 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3375 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3377 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3381 TObject *recoParamObj = entry->GetObject();
3382 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3383 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3384 // Registering them in AliRecoParam
3385 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3387 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3388 // The detector has only onse set of reco parameters
3389 // Registering it in AliRecoParam
3390 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3391 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3392 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3395 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3399 // FIX ME: We have to disable the unloading of reco-param CDB
3400 // entries because QA framework is using them. Has to be fix in
3401 // a way that the QA takes the objects already constructed in
3403 // AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3407 if (AliDebugLevel() > 0) fRecoParam.Print();
3412 //_____________________________________________________________________________
3413 Bool_t AliReconstruction::GetEventInfo()
3415 // Fill the event info object
3417 AliCodeTimerAuto("",0)
3419 AliCentralTrigger *aCTP = NULL;
3421 fEventInfo.SetEventType(fRawReader->GetType());
3423 ULong64_t mask = fRawReader->GetClassMask();
3424 fEventInfo.SetTriggerMask(mask);
3425 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3426 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3428 aCTP = new AliCentralTrigger();
3429 TString configstr("");
3430 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3431 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3435 aCTP->SetClassMask(mask);
3436 aCTP->SetClusterMask(clmask);
3439 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3441 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3442 aCTP = fRunLoader->GetTrigger();
3443 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3444 // get inputs from actp - just get
3445 AliESDHeader* esdheader = fesd->GetHeader();
3446 esdheader->SetL0TriggerInputs(aCTP->GetL0TriggerInputs());
3447 esdheader->SetL1TriggerInputs(aCTP->GetL1TriggerInputs());
3448 esdheader->SetL2TriggerInputs(aCTP->GetL2TriggerInputs());
3449 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3452 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3457 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3459 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3460 if (fRawReader) delete aCTP;
3464 UChar_t clustmask = 0;
3466 ULong64_t trmask = fEventInfo.GetTriggerMask();
3467 const TObjArray& classesArray = config->GetClasses();
3468 Int_t nclasses = classesArray.GetEntriesFast();
3469 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3470 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3472 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3473 fesd->SetTriggerClass(trclass->GetName(),trindex);
3474 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3475 if (trmask & (1ull << trindex)) {
3477 trclasses += trclass->GetName();
3479 clustmask |= trclass->GetCluster()->GetClusterMask();
3483 fEventInfo.SetTriggerClasses(trclasses);
3485 // Set the information in ESD
3486 fesd->SetTriggerMask(trmask);
3487 fesd->SetTriggerCluster(clustmask);
3489 if (!aCTP->CheckTriggeredDetectors()) {
3490 if (fRawReader) delete aCTP;
3494 if (fRawReader) delete aCTP;
3496 // We have to fill also the HLT decision here!!
3502 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3504 // Match the detector list found in the rec.C or the default 'ALL'
3505 // to the list found in the GRP (stored there by the shuttle PP which
3506 // gets the information from ECS)
3507 static TString resultList;
3508 TString detList = detectorList;
3512 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3513 if ((detectorMask >> iDet) & 0x1) {
3514 TString det = AliDAQ::OfflineModuleName(iDet);
3515 if ((detList.CompareTo("ALL") == 0) ||
3516 ((detList.BeginsWith("ALL ") ||
3517 detList.EndsWith(" ALL") ||
3518 detList.Contains(" ALL ")) &&
3519 !(detList.BeginsWith("-"+det+" ") ||
3520 detList.EndsWith(" -"+det) ||
3521 detList.Contains(" -"+det+" "))) ||
3522 (detList.CompareTo(det) == 0) ||
3523 detList.BeginsWith(det+" ") ||
3524 detList.EndsWith(" "+det) ||
3525 detList.Contains( " "+det+" " )) {
3526 if (!resultList.EndsWith(det + " ")) {
3535 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3536 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3537 if ((detList.CompareTo("ALL") == 0) ||
3538 ((detList.BeginsWith("ALL ") ||
3539 detList.EndsWith(" ALL") ||
3540 detList.Contains(" ALL ")) &&
3541 !(detList.BeginsWith("-"+hltDet+" ") ||
3542 detList.EndsWith(" -"+hltDet) ||
3543 detList.Contains(" -"+hltDet+" "))) ||
3544 (detList.CompareTo(hltDet) == 0) ||
3545 detList.BeginsWith(hltDet+" ") ||
3546 detList.EndsWith(" "+hltDet) ||
3547 detList.Contains( " "+hltDet+" " )) {
3548 resultList += hltDet;
3552 return resultList.Data();
3556 //______________________________________________________________________________
3557 void AliReconstruction::Abort(const char *method, EAbort what)
3559 // Abort processing. If what = kAbortProcess, the Process() loop will be
3560 // aborted. If what = kAbortFile, the current file in a chain will be
3561 // aborted and the processing will continue with the next file, if there
3562 // is no next file then Process() will be aborted. Abort() can also be
3563 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3564 // the SlaveTerminate() and Terminate() are always called. The abort flag
3565 // can be checked in these methods using GetAbort().
3567 // The method is overwritten in AliReconstruction for better handling of
3568 // reco specific errors
3570 if (!fStopOnError) return;
3574 TString whyMess = method;
3575 whyMess += " failed! Aborting...";
3577 AliError(whyMess.Data());
3580 TString mess = "Abort";
3581 if (fAbort == kAbortProcess)
3582 mess = "AbortProcess";
3583 else if (fAbort == kAbortFile)
3586 Info(mess, whyMess.Data());
3589 //______________________________________________________________________________
3590 Bool_t AliReconstruction::ProcessEvent(void* event)
3592 // Method that is used in case the event loop
3593 // is steered from outside, for example by AMORE
3594 // 'event' is a pointer to the DATE event in the memory
3596 if (fRawReader) delete fRawReader;
3597 fRawReader = new AliRawReaderDate(event);
3598 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
3605 //______________________________________________________________________________
3606 Bool_t AliReconstruction::ParseOutput()
3608 // The method parses the output file
3609 // location string in order to steer
3610 // properly the selector
3612 TPMERegexp re1("(\\w+\\.zip#\\w+\\.root):([,*\\w+\\.root,*]+)@dataset://(\\w++)");
3613 TPMERegexp re2("(\\w+\\.root)?@?dataset://(\\w++)");
3615 if (re1.Match(fESDOutput) == 4) {
3616 // root archive with output files stored and regustered
3618 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",re1[1].Data()));
3619 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re1[3].Data()));
3620 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3621 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_ARCHIVE",re1[2].Data()));
3622 AliInfo(Form("%s files will be stored within %s in dataset %s",
3627 else if (re2.Match(fESDOutput) == 3) {
3628 // output file stored and registered
3630 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",(re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data()));
3631 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re2[2].Data()));
3632 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3633 AliInfo(Form("%s will be stored in dataset %s",
3634 (re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data(),
3638 if (fESDOutput.IsNull()) {
3639 // Output location not given.
3640 // Assuming xrootd has been already started and
3641 // the output file has to be sent back
3642 // to the client machine
3643 TString esdUrl(Form("root://%s/%s/",
3644 TUrl(gSystem->HostName()).GetHostFQDN(),
3646 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE","AliESDs.root"));
3647 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",esdUrl.Data()));
3648 AliInfo(Form("AliESDs.root will be stored in %s",
3652 // User specified an output location.
3653 // Ones has just to parse it here
3654 TUrl outputUrl(fESDOutput.Data());
3655 TString outputFile(gSystem->BaseName(outputUrl.GetFile()));
3656 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.IsNull() ? "AliESDs.root" : outputFile.Data()));
3657 TString outputLocation(outputUrl.GetUrl());
3658 outputLocation.ReplaceAll(outputFile.Data(),"");
3659 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",outputLocation.Data()));
3660 AliInfo(Form("%s will be stored in %s",
3661 outputFile.IsNull() ? "AliESDs.root" : outputFile.Data(),
3662 outputLocation.Data()));