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>
135 #include "AliAlignObj.h"
136 #include "AliCDBEntry.h"
137 #include "AliCDBManager.h"
138 #include "AliCDBStorage.h"
139 #include "AliCTPRawStream.h"
140 #include "AliCascadeVertexer.h"
141 #include "AliCentralTrigger.h"
142 #include "AliCodeTimer.h"
144 #include "AliDetectorRecoParam.h"
145 #include "AliESDCaloCells.h"
146 #include "AliESDCaloCluster.h"
147 #include "AliESDEvent.h"
148 #include "AliESDMuonTrack.h"
149 #include "AliESDPmdTrack.h"
150 #include "AliESDTagCreator.h"
151 #include "AliESDVertex.h"
152 #include "AliESDcascade.h"
153 #include "AliESDfriend.h"
154 #include "AliESDkink.h"
155 #include "AliESDpid.h"
156 #include "AliESDtrack.h"
157 #include "AliESDtrack.h"
158 #include "AliEventInfo.h"
159 #include "AliGRPObject.h"
160 #include "AliGRPRecoParam.h"
161 #include "AliGenEventHeader.h"
162 #include "AliGeomManager.h"
163 #include "AliGlobalQADataMaker.h"
164 #include "AliHeader.h"
167 #include "AliMultiplicity.h"
169 #include "AliPlaneEff.h"
171 #include "AliQADataMakerRec.h"
172 #include "AliQAManager.h"
173 #include "AliRawVEvent.h"
174 #include "AliRawEventHeaderBase.h"
175 #include "AliRawHLTManager.h"
176 #include "AliRawReaderDate.h"
177 #include "AliRawReaderFile.h"
178 #include "AliRawReaderRoot.h"
179 #include "AliReconstruction.h"
180 #include "AliReconstructor.h"
182 #include "AliRunInfo.h"
183 #include "AliRunLoader.h"
184 #include "AliSysInfo.h" // memory snapshots
185 #include "AliTrackPointArray.h"
186 #include "AliTracker.h"
187 #include "AliTriggerClass.h"
188 #include "AliTriggerCluster.h"
189 #include "AliTriggerIR.h"
190 #include "AliTriggerConfiguration.h"
191 #include "AliV0vertexer.h"
192 #include "AliVertexer.h"
193 #include "AliTrackleter.h"
194 #include "AliVertexerTracks.h"
195 #include "AliTriggerRunScalers.h"
196 #include "AliCTPTimeParams.h"
197 #include "AliESDHLTDecision.h"
198 #include "AliTriggerInput.h"
199 #include "AliLHCData.h"
200 ClassImp(AliReconstruction)
202 //_____________________________________________________________________________
203 const char* AliReconstruction::fgkDetectorName[AliReconstruction::kNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
205 //_____________________________________________________________________________
206 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
208 fRunVertexFinder(kTRUE),
209 fRunVertexFinderTracks(kTRUE),
210 fRunHLTTracking(kFALSE),
211 fRunMuonTracking(kFALSE),
213 fRunCascadeFinder(kTRUE),
214 fRunMultFinder(kTRUE),
216 fWriteAlignmentData(kFALSE),
217 fWriteESDfriend(kFALSE),
218 fFillTriggerESD(kTRUE),
226 fRunLocalReconstruction("ALL"),
230 fUseTrackingErrorsForAlignment(""),
231 fGAliceFileName(gAliceFilename),
234 fProofOutputFileName(""),
235 fProofOutputLocation(""),
236 fProofOutputDataset(kFALSE),
237 fProofOutputArchive(""),
241 fNumberOfEventsPerFile((UInt_t)-1),
242 fFractionFriends(0.04),
244 fLoadAlignFromCDB(kTRUE),
245 fLoadAlignData("ALL"),
250 fCTPTimeParams(NULL),
254 fParentRawReader(NULL),
258 fSPDTrackleter(NULL),
260 fDiamondProfileSPD(NULL),
261 fDiamondProfile(NULL),
262 fDiamondProfileTPC(NULL),
263 fListOfCosmicTriggers(NULL),
267 fAlignObjArray(NULL),
271 fInitCDBCalled(kFALSE),
272 fSetRunNumberFromDataCalled(kFALSE),
277 fSameQACycle(kFALSE),
278 fInitQACalled(kFALSE),
279 fWriteQAExpertData(kTRUE),
280 fRunPlaneEff(kFALSE),
291 fIsNewRunLoader(kFALSE),
295 // create reconstruction object with default parameters
298 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
299 fReconstructor[iDet] = NULL;
300 fLoader[iDet] = NULL;
301 fTracker[iDet] = NULL;
303 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
304 fQACycles[iDet] = 999999 ;
305 fQAWriteExpert[iDet] = kFALSE ;
307 fBeamInt[0][0]=fBeamInt[0][1]=fBeamInt[1][0]=fBeamInt[1][1] = -1;
312 //_____________________________________________________________________________
313 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
315 fRunVertexFinder(rec.fRunVertexFinder),
316 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
317 fRunHLTTracking(rec.fRunHLTTracking),
318 fRunMuonTracking(rec.fRunMuonTracking),
319 fRunV0Finder(rec.fRunV0Finder),
320 fRunCascadeFinder(rec.fRunCascadeFinder),
321 fRunMultFinder(rec.fRunMultFinder),
322 fStopOnError(rec.fStopOnError),
323 fWriteAlignmentData(rec.fWriteAlignmentData),
324 fWriteESDfriend(rec.fWriteESDfriend),
325 fFillTriggerESD(rec.fFillTriggerESD),
327 fCleanESD(rec.fCleanESD),
328 fV0DCAmax(rec.fV0DCAmax),
329 fV0CsPmin(rec.fV0CsPmin),
333 fRunLocalReconstruction(rec.fRunLocalReconstruction),
334 fRunTracking(rec.fRunTracking),
335 fFillESD(rec.fFillESD),
336 fLoadCDB(rec.fLoadCDB),
337 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
338 fGAliceFileName(rec.fGAliceFileName),
339 fRawInput(rec.fRawInput),
340 fESDOutput(rec.fESDOutput),
341 fProofOutputFileName(rec.fProofOutputFileName),
342 fProofOutputLocation(rec.fProofOutputLocation),
343 fProofOutputDataset(rec.fProofOutputDataset),
344 fProofOutputArchive(rec.fProofOutputArchive),
345 fEquipIdMap(rec.fEquipIdMap),
346 fFirstEvent(rec.fFirstEvent),
347 fLastEvent(rec.fLastEvent),
348 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
349 fFractionFriends(rec.fFractionFriends),
351 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
352 fLoadAlignData(rec.fLoadAlignData),
353 fUseHLTData(rec.fUseHLTData),
357 fCTPTimeParams(NULL),
361 fParentRawReader(NULL),
363 fRecoParam(rec.fRecoParam),
365 fSPDTrackleter(NULL),
367 fDiamondProfileSPD(rec.fDiamondProfileSPD),
368 fDiamondProfile(rec.fDiamondProfile),
369 fDiamondProfileTPC(rec.fDiamondProfileTPC),
370 fListOfCosmicTriggers(NULL),
374 fAlignObjArray(rec.fAlignObjArray),
375 fCDBUri(rec.fCDBUri),
376 fQARefUri(rec.fQARefUri),
378 fInitCDBCalled(rec.fInitCDBCalled),
379 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
380 fQADetectors(rec.fQADetectors),
381 fQATasks(rec.fQATasks),
383 fRunGlobalQA(rec.fRunGlobalQA),
384 fSameQACycle(rec.fSameQACycle),
385 fInitQACalled(rec.fInitQACalled),
386 fWriteQAExpertData(rec.fWriteQAExpertData),
387 fRunPlaneEff(rec.fRunPlaneEff),
398 fIsNewRunLoader(rec.fIsNewRunLoader),
404 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
405 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
407 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
408 fReconstructor[iDet] = NULL;
409 fLoader[iDet] = NULL;
410 fTracker[iDet] = NULL;
413 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
414 fQACycles[iDet] = rec.fQACycles[iDet];
415 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
418 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
419 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
422 for (int i=2;i--;) for (int j=2;j--;) fBeamInt[i][j] = rec.fBeamInt[i][j];
426 //_____________________________________________________________________________
427 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
429 // assignment operator
430 // Used in PROOF mode
431 // Be very careful while modifing it!
432 // Simple rules to follow:
433 // for persistent data members - use their assignment operators
434 // for non-persistent ones - do nothing or take the default values from constructor
435 // TSelector members should not be touched
436 if(&rec == this) return *this;
438 fRunVertexFinder = rec.fRunVertexFinder;
439 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
440 fRunHLTTracking = rec.fRunHLTTracking;
441 fRunMuonTracking = rec.fRunMuonTracking;
442 fRunV0Finder = rec.fRunV0Finder;
443 fRunCascadeFinder = rec.fRunCascadeFinder;
444 fRunMultFinder = rec.fRunMultFinder;
445 fStopOnError = rec.fStopOnError;
446 fWriteAlignmentData = rec.fWriteAlignmentData;
447 fWriteESDfriend = rec.fWriteESDfriend;
448 fFillTriggerESD = rec.fFillTriggerESD;
450 fCleanESD = rec.fCleanESD;
451 fV0DCAmax = rec.fV0DCAmax;
452 fV0CsPmin = rec.fV0CsPmin;
456 fRunLocalReconstruction = rec.fRunLocalReconstruction;
457 fRunTracking = rec.fRunTracking;
458 fFillESD = rec.fFillESD;
459 fLoadCDB = rec.fLoadCDB;
460 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
461 fGAliceFileName = rec.fGAliceFileName;
462 fRawInput = rec.fRawInput;
463 fESDOutput = rec.fESDOutput;
464 fProofOutputFileName = rec.fProofOutputFileName;
465 fProofOutputLocation = rec.fProofOutputLocation;
466 fProofOutputDataset = rec.fProofOutputDataset;
467 fProofOutputArchive = rec.fProofOutputArchive;
468 fEquipIdMap = rec.fEquipIdMap;
469 fFirstEvent = rec.fFirstEvent;
470 fLastEvent = rec.fLastEvent;
471 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
472 fFractionFriends = rec.fFractionFriends;
474 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
475 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
478 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
479 fLoadAlignData = rec.fLoadAlignData;
480 fUseHLTData = rec.fUseHLTData;
482 delete fRunInfo; fRunInfo = NULL;
483 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
485 fEventInfo = rec.fEventInfo;
487 delete fRunScalers; fRunScalers = NULL;
488 if (rec.fRunScalers) fRunScalers = new AliTriggerRunScalers(*rec.fRunScalers);
490 delete fCTPTimeParams; fCTPTimeParams = NULL;
491 if (rec.fCTPTimeParams) fCTPTimeParams = new AliCTPTimeParams(*rec.fCTPTimeParams);
495 fParentRawReader = NULL;
497 fRecoParam = rec.fRecoParam;
499 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
500 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
501 delete fLoader[iDet]; fLoader[iDet] = NULL;
502 delete fTracker[iDet]; fTracker[iDet] = NULL;
505 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
506 fQACycles[iDet] = rec.fQACycles[iDet];
507 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
510 delete fSPDTrackleter; fSPDTrackleter = NULL;
512 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
513 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
514 delete fDiamondProfile; fDiamondProfile = NULL;
515 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
516 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
517 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
519 delete fListOfCosmicTriggers; fListOfCosmicTriggers = NULL;
520 if (rec.fListOfCosmicTriggers) fListOfCosmicTriggers = (THashTable*)((rec.fListOfCosmicTriggers)->Clone());
522 delete fGRPData; fGRPData = NULL;
523 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
524 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
526 delete fAlignObjArray; fAlignObjArray = NULL;
529 fQARefUri = rec.fQARefUri;
530 fSpecCDBUri.Delete();
531 fInitCDBCalled = rec.fInitCDBCalled;
532 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
533 fQADetectors = rec.fQADetectors;
534 fQATasks = rec.fQATasks;
536 fRunGlobalQA = rec.fRunGlobalQA;
537 fSameQACycle = rec.fSameQACycle;
538 fInitQACalled = rec.fInitQACalled;
539 fWriteQAExpertData = rec.fWriteQAExpertData;
540 fRunPlaneEff = rec.fRunPlaneEff;
541 for (int i=2;i--;) for (int j=2;j--;) fBeamInt[i][j] = rec.fBeamInt[i][j];
551 fIsNewRunLoader = rec.fIsNewRunLoader;
558 //_____________________________________________________________________________
559 AliReconstruction::~AliReconstruction()
564 if (fListOfCosmicTriggers) {
565 fListOfCosmicTriggers->Delete();
566 delete fListOfCosmicTriggers;
570 delete fCTPTimeParams;
572 if (fAlignObjArray) {
573 fAlignObjArray->Delete();
574 delete fAlignObjArray;
576 fSpecCDBUri.Delete();
578 AliCodeTimer::Instance()->Print();
581 //_____________________________________________________________________________
582 void AliReconstruction::InitQA()
584 //Initialize the QA and start of cycle
585 AliCodeTimerAuto("",0);
587 if (fInitQACalled) return;
588 fInitQACalled = kTRUE;
590 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kRECMODE) ;
591 if (fWriteQAExpertData)
592 qam->SetWriteExpert() ;
594 if (qam->IsDefaultStorageSet()) {
595 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
596 AliWarning("Default QA reference storage has been already set !");
597 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fQARefUri.Data()));
598 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
599 fQARefUri = qam->GetDefaultStorage()->GetURI();
601 if (fQARefUri.Length() > 0) {
602 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
603 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
604 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
606 fQARefUri="local://$ALICE_ROOT/QAref";
607 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
608 AliWarning("Default QA refeference storage not yet set !!!!");
609 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
610 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
613 qam->SetDefaultStorage(fQARefUri);
617 qam->SetActiveDetectors(fQADetectors) ;
618 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) {
619 qam->SetCycleLength(AliQAv1::DETECTORINDEX_t(det), fQACycles[det]) ;
620 qam->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
622 if (!fRawReader && !fInput && IsInTasks(AliQAv1::kRAWS))
623 fQATasks.ReplaceAll(Form("%d",AliQAv1::kRAWS), "") ;
624 qam->SetTasks(fQATasks) ;
625 qam->InitQADataMaker(AliCDBManager::Instance()->GetRun()) ;
628 Bool_t sameCycle = kFALSE ;
629 AliQADataMaker *qadm = qam->GetQADataMaker(AliQAv1::kGLOBAL);
630 AliInfo(Form("Initializing the global QA data maker"));
631 if (IsInTasks(AliQAv1::kRECPOINTS)) {
632 qadm->StartOfCycle(AliQAv1::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
633 TObjArray **arr=qadm->Init(AliQAv1::kRECPOINTS);
634 AliTracker::SetResidualsArray(arr);
637 if (IsInTasks(AliQAv1::kESDS)) {
638 qadm->StartOfCycle(AliQAv1::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
639 qadm->Init(AliQAv1::kESDS);
642 AliSysInfo::AddStamp("InitQA") ;
645 //_____________________________________________________________________________
646 void AliReconstruction::MergeQA(const char *fileName)
648 //Initialize the QA and start of cycle
649 AliCodeTimerAuto("",0) ;
650 AliQAManager::QAManager()->Merge(AliCDBManager::Instance()->GetRun(),fileName) ;
651 AliSysInfo::AddStamp("MergeQA") ;
654 //_____________________________________________________________________________
655 void AliReconstruction::InitCDB()
657 // activate a default CDB storage
658 // First check if we have any CDB storage set, because it is used
659 // to retrieve the calibration and alignment constants
660 AliCodeTimerAuto("",0);
662 if (fInitCDBCalled) return;
663 fInitCDBCalled = kTRUE;
665 AliCDBManager* man = AliCDBManager::Instance();
666 if (man->IsDefaultStorageSet())
668 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
669 AliWarning("Default CDB storage has been already set !");
670 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
671 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
672 fCDBUri = man->GetDefaultStorage()->GetURI();
675 if (fCDBUri.Length() > 0)
677 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
678 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
679 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
680 man->SetDefaultStorage(fCDBUri);
682 else if (!man->GetRaw()){
683 fCDBUri="local://$ALICE_ROOT/OCDB";
684 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
685 AliWarning("Default CDB storage not yet set !!!!");
686 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
687 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
688 man->SetDefaultStorage(fCDBUri);
691 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
692 AliWarning("Default storage will be set after setting the Run Number!!!");
693 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
697 // Now activate the detector specific CDB storage locations
698 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
699 TObject* obj = fSpecCDBUri[i];
701 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
702 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
703 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
704 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
706 AliSysInfo::AddStamp("InitCDB");
709 //_____________________________________________________________________________
710 void AliReconstruction::SetDefaultStorage(const char* uri) {
711 // Store the desired default CDB storage location
712 // Activate it later within the Run() method
718 //_____________________________________________________________________________
719 void AliReconstruction::SetQARefDefaultStorage(const char* uri) {
720 // Store the desired default CDB storage location
721 // Activate it later within the Run() method
724 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
727 //_____________________________________________________________________________
728 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
729 // Store a detector-specific CDB storage location
730 // Activate it later within the Run() method
732 AliCDBPath aPath(calibType);
733 if(!aPath.IsValid()){
734 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
735 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
736 if(!strcmp(calibType, fgkDetectorName[iDet])) {
737 aPath.SetPath(Form("%s/*", calibType));
738 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
742 if(!aPath.IsValid()){
743 AliError(Form("Not a valid path or detector: %s", calibType));
748 // // check that calibType refers to a "valid" detector name
749 // Bool_t isDetector = kFALSE;
750 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
751 // TString detName = fgkDetectorName[iDet];
752 // if(aPath.GetLevel0() == detName) {
753 // isDetector = kTRUE;
759 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
763 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
764 if (obj) fSpecCDBUri.Remove(obj);
765 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
769 //_____________________________________________________________________________
770 Bool_t AliReconstruction::SetRunNumberFromData()
772 // The method is called in Run() in order
773 // to set a correct run number.
774 // In case of raw data reconstruction the
775 // run number is taken from the raw data header
777 if (fSetRunNumberFromDataCalled) return kTRUE;
778 fSetRunNumberFromDataCalled = kTRUE;
780 AliCDBManager* man = AliCDBManager::Instance();
783 if(fRawReader->NextEvent()) {
784 if(man->GetRun() > 0) {
785 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
787 man->SetRun(fRawReader->GetRunNumber());
788 fRawReader->RewindEvents();
791 if(man->GetRun() > 0) {
792 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
795 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
801 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
803 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
808 // read run number from gAlice
809 if(rl->GetHeader()) {
810 man->SetRun(rl->GetHeader()->GetRun());
815 AliError("Neither run-loader header nor RawReader objects are found !");
827 //_____________________________________________________________________________
828 void AliReconstruction::SetCDBLock() {
829 // Set CDB lock: from now on it is forbidden to reset the run number
830 // or the default storage or to activate any further storage!
832 AliCDBManager::Instance()->SetLock(1);
835 //_____________________________________________________________________________
836 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
838 // Read the alignment objects from CDB.
839 // Each detector is supposed to have the
840 // alignment objects in DET/Align/Data CDB path.
841 // All the detector objects are then collected,
842 // sorted by geometry level (starting from ALIC) and
843 // then applied to the TGeo geometry.
844 // Finally an overlaps check is performed.
846 // Load alignment data from CDB and fill fAlignObjArray
847 if(fLoadAlignFromCDB){
849 TString detStr = detectors;
850 TString loadAlObjsListOfDets = "";
852 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
853 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
854 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
856 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
858 loadAlObjsListOfDets += fgkDetectorName[iDet];
859 loadAlObjsListOfDets += " ";
861 } // end loop over detectors
863 if(AliGeomManager::GetNalignable("GRP") != 0)
864 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
865 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
866 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
868 // Check if the array with alignment objects was
869 // provided by the user. If yes, apply the objects
870 // to the present TGeo geometry
871 if (fAlignObjArray) {
872 if (gGeoManager && gGeoManager->IsClosed()) {
873 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
874 AliError("The misalignment of one or more volumes failed!"
875 "Compare the list of simulated detectors and the list of detector alignment data!");
880 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
886 if (fAlignObjArray) {
887 fAlignObjArray->Delete();
888 delete fAlignObjArray; fAlignObjArray=NULL;
894 //_____________________________________________________________________________
895 void AliReconstruction::SetGAliceFile(const char* fileName)
897 // set the name of the galice file
899 fGAliceFileName = fileName;
902 //_____________________________________________________________________________
903 void AliReconstruction::SetInput(const char* input)
905 // In case the input string starts with 'mem://', we run in an online mode
906 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
907 // file is assumed. One can give as an input:
908 // mem://: - events taken from DAQ monitoring libs online
910 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
911 if (input) fRawInput = input;
914 //_____________________________________________________________________________
915 void AliReconstruction::SetOutput(const char* output)
917 // Set the output ESD filename
918 // 'output' is a normalt ROOT url
919 // The method is used in case of raw-data reco with PROOF
920 if (output) fESDOutput = output;
923 //_____________________________________________________________________________
924 void AliReconstruction::SetOption(const char* detector, const char* option)
926 // set options for the reconstruction of a detector
928 TObject* obj = fOptions.FindObject(detector);
929 if (obj) fOptions.Remove(obj);
930 fOptions.Add(new TNamed(detector, option));
933 //_____________________________________________________________________________
934 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
936 // Set custom reconstruction parameters for a given detector
937 // Single set of parameters for all the events
939 // First check if the reco-params are global
940 if(!strcmp(detector, "GRP")) {
942 fRecoParam.AddDetRecoParam(kNDetectors,par);
946 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
947 if(!strcmp(detector, fgkDetectorName[iDet])) {
949 fRecoParam.AddDetRecoParam(iDet,par);
956 //_____________________________________________________________________________
957 Bool_t AliReconstruction::InitGRP() {
958 //------------------------------------
959 // Initialization of the GRP entry
960 //------------------------------------
961 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
965 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
968 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
970 fGRPData = new AliGRPObject();
971 fGRPData->ReadValuesFromMap(m);
975 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
976 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
980 // FIX ME: The unloading of GRP entry is temporarily disabled
981 // because ZDC and VZERO are using it in order to initialize
982 // their reconstructor objects. In the future one has to think
983 // of propagating AliRunInfo to the reconstructors.
984 // AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
988 AliError("No GRP entry found in OCDB!");
992 TString lhcState = fGRPData->GetLHCState();
993 if (lhcState==AliGRPObject::GetInvalidString()) {
994 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
995 lhcState = "UNKNOWN";
998 TString beamType = fGRPData->GetBeamType();
999 if (beamType==AliGRPObject::GetInvalidString()) {
1000 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
1001 beamType = "UNKNOWN";
1004 Float_t beamEnergy = fGRPData->GetBeamEnergy();
1005 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
1006 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
1010 TString runType = fGRPData->GetRunType();
1011 if (runType==AliGRPObject::GetInvalidString()) {
1012 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
1013 runType = "UNKNOWN";
1016 Int_t activeDetectors = fGRPData->GetDetectorMask();
1017 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
1018 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
1019 activeDetectors = 1074790399;
1022 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1026 // Process the list of active detectors
1027 if (activeDetectors) {
1028 UInt_t detMask = activeDetectors;
1029 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1030 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1031 fFillESD = MatchDetectorList(fFillESD,detMask);
1032 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1033 fLoadCDB.Form("%s %s %s %s",
1034 fRunLocalReconstruction.Data(),
1035 fRunTracking.Data(),
1037 fQADetectors.Data());
1038 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1039 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1) &&
1040 !((detMask >> AliDAQ::DetectorID("ITSSDD")) & 0x1) &&
1041 !((detMask >> AliDAQ::DetectorID("ITSSSD")) & 0x1) ) {
1042 // switch off the vertexer
1043 AliInfo("SPD,SDD,SSD is not in the list of active detectors. Vertexer and Trackleter are switched off.");
1044 fRunVertexFinder = kFALSE;
1045 fRunMultFinder = kFALSE;
1047 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1048 // switch off the reading of CTP raw-data payload
1049 if (fFillTriggerESD) {
1050 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1051 fFillTriggerESD = kFALSE;
1056 AliInfo("===================================================================================");
1057 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1058 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1059 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1060 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1061 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1062 AliInfo("===================================================================================");
1064 //*** Dealing with the magnetic field map
1065 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {
1066 if (TGeoGlobalMagField::Instance()->GetField()->TestBit(AliMagF::kOverrideGRP)) {
1067 AliInfo("ExpertMode!!! GRP information will be ignored !");
1068 AliInfo("ExpertMode!!! Running with the externally locked B field !");
1071 AliInfo("Destroying existing B field instance!");
1072 delete TGeoGlobalMagField::Instance();
1075 if ( !TGeoGlobalMagField::Instance()->IsLocked() ) {
1076 // Construct the field map out of the information retrieved from GRP.
1079 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1080 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1081 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1085 Char_t l3Polarity = fGRPData->GetL3Polarity();
1086 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1087 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1092 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1093 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1094 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1098 Char_t diPolarity = fGRPData->GetDipolePolarity();
1099 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1100 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1104 // read special bits for the polarity convention and map type
1105 Int_t polConvention = fGRPData->IsPolarityConventionLHC() ? AliMagF::kConvLHC : AliMagF::kConvDCS2008;
1106 Bool_t uniformB = fGRPData->IsUniformBMap();
1109 AliMagF* fld = AliMagF::CreateFieldMap(TMath::Abs(l3Current) * (l3Polarity ? -1:1),
1110 TMath::Abs(diCurrent) * (diPolarity ? -1:1),
1111 polConvention,uniformB,beamEnergy, beamType.Data());
1113 TGeoGlobalMagField::Instance()->SetField( fld );
1114 TGeoGlobalMagField::Instance()->Lock();
1115 AliInfo("Running with the B field constructed out of GRP !");
1117 else AliFatal("Failed to create a B field map !");
1119 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1122 //*** Get the diamond profiles from OCDB
1123 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1125 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1127 AliError("No SPD diamond profile found in OCDB!");
1130 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1132 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1134 AliError("No diamond profile found in OCDB!");
1137 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1139 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1141 AliError("No TPC diamond profile found in OCDB!");
1144 entry = AliCDBManager::Instance()->Get("GRP/Calib/CosmicTriggers");
1146 fListOfCosmicTriggers = dynamic_cast<THashTable*>(entry->GetObject());
1148 AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
1151 if (!fListOfCosmicTriggers) {
1152 AliWarning("Can not get list of cosmic triggers from OCDB! Cosmic event specie will be effectively disabled!");
1158 //_____________________________________________________________________________
1159 Bool_t AliReconstruction::LoadCDB()
1161 // Load CDB entries for all active detectors.
1162 // By default we load all the entries in <det>/Calib
1165 AliCodeTimerAuto("",0);
1167 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1169 TString detStr = fLoadCDB;
1170 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1171 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1172 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1175 // Temporary fix - one has to define the correct policy in order
1176 // to load the trigger OCDB entries only for the detectors that
1177 // in the trigger or that are needed in order to put correct
1178 // information in ESD
1179 AliCDBManager::Instance()->GetAll("TRIGGER/*/*");
1183 //_____________________________________________________________________________
1184 Bool_t AliReconstruction::LoadTriggerScalersCDB()
1186 // Load CTP scalers from OCDB.
1187 // The scalers are checked for consistency.
1189 AliCodeTimerAuto("",0);
1191 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
1195 AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
1196 fRunScalers = dynamic_cast<AliTriggerRunScalers*> (entry->GetObject());
1198 if (fRunScalers->CorrectScalersOverflow() == 0) AliInfo("32bit Trigger counters corrected for overflow");
1203 //_____________________________________________________________________________
1204 Bool_t AliReconstruction::LoadCTPTimeParamsCDB()
1206 // Load CTP timing information (alignment)
1209 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
1213 AliInfo("Found an AliCTPTimeParams in GRP/CTP/CTPtiming, reading it");
1214 fCTPTimeParams = dynamic_cast<AliCTPTimeParams*> (entry->GetObject());
1222 //_____________________________________________________________________________
1223 Bool_t AliReconstruction::ReadIntensityInfoCDB()
1225 // Load LHC DIP data
1226 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/LHCData");
1229 AliInfo("Found an AliLHCData in GRP/GRP/LHCData, reading it");
1230 AliLHCData* dipData = dynamic_cast<AliLHCData*> (entry->GetObject());
1231 for (int ib=2;ib--;) {
1233 if (dipData->GetMeanIntensity(ib,intI,intNI)>=0) {
1234 fBeamInt[ib][0] = intI;
1235 fBeamInt[ib][1] = intNI;
1244 //_____________________________________________________________________________
1245 Bool_t AliReconstruction::Run(const char* input)
1248 AliCodeTimerAuto("",0);
1251 if (GetAbort() != TSelector::kContinue) return kFALSE;
1253 TChain *chain = NULL;
1254 if (fRawReader && (chain = fRawReader->GetChain())) {
1255 Long64_t nEntries = (fLastEvent < 0) ? (TChain::kBigNumber) : (fLastEvent - fFirstEvent + 1);
1258 // Temporary fix for long raw-data runs (until socket timeout handling in PROOF is revised)
1259 gProof->Exec("gEnv->SetValue(\"Proof.SocketActivityTimeout\",-1)", kTRUE);
1262 gProof->Exec("TGrid::Connect(\"alien://\")",kTRUE);
1264 TMessage::EnableSchemaEvolutionForAll(kTRUE);
1265 gProof->Exec("TMessage::EnableSchemaEvolutionForAll(kTRUE)",kTRUE);
1267 gProof->AddInput(this);
1269 if (!ParseOutput()) return kFALSE;
1271 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1273 chain->Process("AliReconstruction","",nEntries,fFirstEvent);
1276 chain->Process(this,"",nEntries,fFirstEvent);
1281 if (GetAbort() != TSelector::kContinue) return kFALSE;
1283 if (GetAbort() != TSelector::kContinue) return kFALSE;
1284 //******* The loop over events
1285 AliInfo("Starting looping over events");
1287 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1288 (fRawReader && fRawReader->NextEvent())) {
1289 if (!ProcessEvent(iEvent)) {
1290 Abort("ProcessEvent",TSelector::kAbortFile);
1296 if (GetAbort() != TSelector::kContinue) return kFALSE;
1298 if (GetAbort() != TSelector::kContinue) return kFALSE;
1304 //_____________________________________________________________________________
1305 void AliReconstruction::InitRawReader(const char* input)
1307 // Init raw-reader and
1308 // set the input in case of raw data
1310 AliCodeTimerAuto("",0);
1312 if (input) fRawInput = input;
1313 fRawReader = AliRawReader::Create(fRawInput.Data());
1315 if (fRawInput.IsNull()) {
1316 AliInfo("Reconstruction will run over digits");
1319 AliFatal("Can not create raw-data reader ! Exiting...");
1323 if (!fEquipIdMap.IsNull() && fRawReader)
1324 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1326 if (!fUseHLTData.IsNull()) {
1327 // create the RawReaderHLT which performs redirection of HLT input data for
1328 // the specified detectors
1329 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1331 fParentRawReader=fRawReader;
1332 fRawReader=pRawReader;
1334 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1337 AliSysInfo::AddStamp("CreateRawReader");
1340 //_____________________________________________________________________________
1341 void AliReconstruction::InitRun(const char* input)
1343 // Initialization of raw-reader,
1344 // run number, CDB etc.
1345 AliCodeTimerAuto("",0);
1346 AliSysInfo::AddStamp("Start");
1348 // Initialize raw-reader if any
1349 InitRawReader(input);
1351 // Initialize the CDB storage
1354 // Set run number in CDBManager (if it is not already set by the user)
1355 if (!SetRunNumberFromData()) {
1356 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1360 // Set CDB lock: from now on it is forbidden to reset the run number
1361 // or the default storage or to activate any further storage!
1366 //_____________________________________________________________________________
1367 void AliReconstruction::Begin(TTree *)
1369 // Initialize AlReconstruction before
1370 // going into the event loop
1371 // Should follow the TSelector convention
1372 // i.e. initialize only the object on the client side
1373 AliCodeTimerAuto("",0);
1375 AliReconstruction *reco = NULL;
1377 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1380 AliSysInfo::AddStamp("ReadInputInBegin");
1383 // Import ideal TGeo geometry and apply misalignment
1385 TString geom(gSystem->DirName(fGAliceFileName));
1386 geom += "/geometry.root";
1387 AliGeomManager::LoadGeometry(geom.Data());
1389 Abort("LoadGeometry", TSelector::kAbortProcess);
1392 AliSysInfo::AddStamp("LoadGeom");
1393 TString detsToCheck=fRunLocalReconstruction;
1394 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1395 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1398 AliSysInfo::AddStamp("CheckGeom");
1401 if (!MisalignGeometry(fLoadAlignData)) {
1402 Abort("MisalignGeometry", TSelector::kAbortProcess);
1405 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1406 AliSysInfo::AddStamp("MisalignGeom");
1409 Abort("InitGRP", TSelector::kAbortProcess);
1412 AliSysInfo::AddStamp("InitGRP");
1415 Abort("LoadCDB", TSelector::kAbortProcess);
1418 AliSysInfo::AddStamp("LoadCDB");
1420 if (!LoadTriggerScalersCDB()) {
1421 Abort("LoadTriggerScalersCDB", TSelector::kAbortProcess);
1424 AliSysInfo::AddStamp("LoadTriggerScalersCDB");
1426 if (!LoadCTPTimeParamsCDB()) {
1427 Abort("LoadCTPTimeParamsCDB", TSelector::kAbortProcess);
1430 AliSysInfo::AddStamp("LoadCTPTimeParamsCDB");
1432 if (!ReadIntensityInfoCDB()) {
1433 Abort("ReadIntensityInfoCDB", TSelector::kAbortProcess);
1436 AliSysInfo::AddStamp("ReadIntensityInfoCDB");
1438 // Read the reconstruction parameters from OCDB
1439 if (!InitRecoParams()) {
1440 AliWarning("Not all detectors have correct RecoParam objects initialized");
1442 AliSysInfo::AddStamp("InitRecoParams");
1444 if (fInput && gProof) {
1445 if (reco) *reco = *this;
1447 gGeoManager->SetName("Geometry");
1448 gProof->AddInputData(gGeoManager,kTRUE);
1450 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1451 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1452 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1453 magFieldMap->SetName("MagneticFieldMap");
1454 gProof->AddInputData(magFieldMap,kTRUE);
1459 //_____________________________________________________________________________
1460 void AliReconstruction::SlaveBegin(TTree*)
1462 // Initialization related to run-loader,
1463 // vertexer, trackers, recontructors
1464 // In proof mode it is executed on the slave
1465 AliCodeTimerAuto("",0);
1467 TProofOutputFile *outProofFile = NULL;
1469 if (AliDebugLevel() > 0) fInput->Print();
1470 if (AliDebugLevel() > 10) fInput->Dump();
1471 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1474 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1476 AliGeomManager::SetGeometry(tgeo);
1478 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1479 Int_t runNumber = -1;
1480 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1481 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1482 man->SetCacheFlag(kTRUE);
1483 man->SetLock(kTRUE);
1487 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1488 AliMagF *newMap = new AliMagF(*map);
1489 if (!newMap->LoadParameterization()) {
1490 Abort("AliMagF::LoadParameterization", TSelector::kAbortProcess);
1493 TGeoGlobalMagField::Instance()->SetField(newMap);
1494 TGeoGlobalMagField::Instance()->Lock();
1496 if (TNamed *outputFileName = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE"))
1497 fProofOutputFileName = outputFileName->GetTitle();
1498 if (TNamed *outputLocation = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_LOCATION"))
1499 fProofOutputLocation = outputLocation->GetTitle();
1500 if (fInput->FindObject("PROOF_OUTPUTFILE_DATASET"))
1501 fProofOutputDataset = kTRUE;
1502 if (TNamed *archiveList = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_ARCHIVE"))
1503 fProofOutputArchive = archiveList->GetTitle();
1504 if (!fProofOutputFileName.IsNull() &&
1505 !fProofOutputLocation.IsNull() &&
1506 fProofOutputArchive.IsNull()) {
1507 if (!fProofOutputDataset) {
1508 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"M");
1509 outProofFile->SetOutputFileName(Form("%s%s",fProofOutputLocation.Data(),fProofOutputFileName.Data()));
1512 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"DROV",fProofOutputLocation.Data());
1514 if (AliDebugLevel() > 0) outProofFile->Dump();
1515 fOutput->Add(outProofFile);
1517 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1520 // get the run loader
1521 if (!InitRunLoader()) {
1522 Abort("InitRunLoader", TSelector::kAbortProcess);
1525 AliSysInfo::AddStamp("LoadLoader");
1527 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1530 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1531 Abort("CreateTrackers", TSelector::kAbortProcess);
1534 AliSysInfo::AddStamp("CreateTrackers");
1536 // create the ESD output file and tree
1537 if (!outProofFile) {
1538 ffile = TFile::Open("AliESDs.root", "RECREATE");
1539 ffile->SetCompressionLevel(2);
1540 if (!ffile->IsOpen()) {
1541 Abort("OpenESDFile", TSelector::kAbortProcess);
1546 AliInfo(Form("Opening output PROOF file: %s/%s",
1547 outProofFile->GetDir(), outProofFile->GetFileName()));
1548 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1549 Abort(Form("Problems opening output PROOF file: %s/%s",
1550 outProofFile->GetDir(), outProofFile->GetFileName()),
1551 TSelector::kAbortProcess);
1556 ftree = new TTree("esdTree", "Tree with ESD objects");
1557 fesd = new AliESDEvent();
1558 fesd->CreateStdContent();
1559 // add a so far non-std object to the ESD, this will
1560 // become part of the std content
1561 fesd->AddObject(new AliESDHLTDecision);
1563 fesd->WriteToTree(ftree);
1564 if (fWriteESDfriend) {
1565 ffileF = TFile::Open("AliESDfriends.root", "RECREATE");
1566 ftreeF = new TTree("esdFriendTree", "Tree with ESD Friend objects");
1567 fesdf = new AliESDfriend();
1568 ftreeF->Branch("ESDfriend.","AliESDfriend", &fesdf);
1569 fesd->AddObject(fesdf);
1572 ftree->GetUserInfo()->Add(fesd);
1574 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1575 fhltesd = new AliESDEvent();
1576 fhltesd->CreateStdContent();
1577 // read the ESD template from CDB
1578 // HLT is allowed to put non-std content to its ESD, the non-std
1579 // objects need to be created before invocation of WriteToTree in
1580 // order to create all branches. Initialization is done from an
1581 // ESD layout template in CDB
1582 AliCDBManager* man = AliCDBManager::Instance();
1583 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1584 AliCDBEntry* hltESDConfig=NULL;
1585 if (man->GetId(hltESDConfigPath)!=NULL &&
1586 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1587 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1589 // init all internal variables from the list of objects
1590 pESDLayout->GetStdContent();
1592 // copy content and create non-std objects
1593 *fhltesd=*pESDLayout;
1596 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1597 hltESDConfigPath.GetPath().Data()));
1601 fhltesd->WriteToTree(fhlttree);
1602 fhlttree->GetUserInfo()->Add(fhltesd);
1604 ProcInfo_t procInfo;
1605 gSystem->GetProcInfo(&procInfo);
1606 AliInfo(Form("Current memory usage %d %d", procInfo.fMemResident, procInfo.fMemVirtual));
1609 //Initialize the QA and start of cycle
1610 if (fRunQA || fRunGlobalQA)
1613 //Initialize the Plane Efficiency framework
1614 if (fRunPlaneEff && !InitPlaneEff()) {
1615 Abort("InitPlaneEff", TSelector::kAbortProcess);
1619 if (strcmp(gProgName,"alieve") == 0)
1620 fRunAliEVE = InitAliEVE();
1625 //_____________________________________________________________________________
1626 Bool_t AliReconstruction::Process(Long64_t entry)
1628 // run the reconstruction over a single entry
1629 // from the chain with raw data
1630 AliCodeTimerAuto("",0);
1632 TTree *currTree = fChain->GetTree();
1633 AliRawVEvent *event = NULL;
1634 currTree->SetBranchAddress("rawevent",&event);
1635 currTree->GetEntry(entry);
1636 fRawReader = new AliRawReaderRoot(event);
1637 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1645 //_____________________________________________________________________________
1646 void AliReconstruction::Init(TTree *tree)
1648 // Implementation of TSelector::Init()
1651 AliError("The input tree is not found!");
1657 //_____________________________________________________________________________
1658 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1660 // run the reconstruction over a single event
1661 // The event loop is steered in Run method
1664 static Long_t oldMres=0;
1665 static Long_t oldMvir=0;
1666 static Float_t oldCPU=0;
1667 static Long_t aveDMres=0;
1668 static Long_t aveDMvir=0;
1669 static Float_t aveDCPU=0;
1671 AliCodeTimerAuto("",0);
1675 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1676 fRunLoader->SetEventNumber(iEvent);
1677 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1679 fRunLoader->TreeE()->Fill();
1680 if (fRawReader && fRawReader->UseAutoSaveESD())
1681 fRunLoader->TreeE()->AutoSave("SaveSelf");
1684 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1689 fRunLoader->GetEvent(iEvent);
1691 // Fill Event-info object
1693 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo,fListOfCosmicTriggers);
1695 ProcInfo_t procInfo;
1696 if(iEvent==fFirstEvent) {
1697 gSystem->GetProcInfo(&procInfo);
1698 oldMres=procInfo.fMemResident;
1699 oldMvir=procInfo.fMemVirtual;
1700 oldCPU=procInfo.fCpuUser+procInfo.fCpuSys;
1702 AliInfo(Form("================================= Processing event %d of type %-10s ==================================", iEvent,fRecoParam.PrintEventSpecie()));
1704 // Set the reco-params
1706 TString detStr = fLoadCDB;
1707 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1708 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1709 AliReconstructor *reconstructor = GetReconstructor(iDet);
1710 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1711 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1712 reconstructor->SetRecoParam(par);
1713 reconstructor->GetPidSettings(&pid);
1714 reconstructor->SetEventInfo(&fEventInfo);
1716 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1717 if (par) AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(par->GetEventSpecie())) ;
1722 const AliDetectorRecoParam *grppar = fRecoParam.GetDetRecoParam(kNDetectors);
1723 AliQAManager::QAManager()->SetRecoParam(AliQAv1::kGLOBAL, grppar) ;
1724 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(grppar->GetEventSpecie())) ;
1729 if (fRunQA && IsInTasks(AliQAv1::kRAWS)) {
1730 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1731 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1733 // local single event reconstruction
1734 if (!fRunLocalReconstruction.IsNull()) {
1735 TString detectors=fRunLocalReconstruction;
1736 // run HLT event reconstruction first
1737 // ;-( IsSelected changes the string
1738 if (IsSelected("HLT", detectors) &&
1739 !RunLocalEventReconstruction("HLT")) {
1740 if (fStopOnError) {CleanUp(); return kFALSE;}
1742 detectors=fRunLocalReconstruction;
1743 detectors.ReplaceAll("HLT", "");
1744 if (!RunLocalEventReconstruction(detectors)) {
1753 // fill Event header information from the RawEventHeader
1754 if (fRawReader){FillRawEventHeaderESD(fesd);}
1755 if (fRawReader){FillRawEventHeaderESD(fhltesd);}
1757 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1758 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1760 ((AliESDRun*)fesd->GetESDRun())->SetDetectorsInDAQ(fRunInfo->GetDetectorMask());
1761 ((AliESDRun*)fhltesd->GetESDRun())->SetDetectorsInDAQ(fRunInfo->GetDetectorMask());
1762 ((AliESDRun*)fesd->GetESDRun())->SetDetectorsInReco(AliDAQ::DetectorPatternOffline(fFillESD.Data()));
1763 ((AliESDRun*)fhltesd->GetESDRun())->SetDetectorsInReco(AliDAQ::DetectorPatternOffline(fFillESD.Data()));
1765 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1766 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1768 fesd->SetEventSpecie(fRecoParam.GetEventSpecie());
1769 fhltesd->SetEventSpecie(fRecoParam.GetEventSpecie());
1771 // Set magnetic field from the tracker
1772 fesd->SetMagneticField(AliTracker::GetBz());
1773 fhltesd->SetMagneticField(AliTracker::GetBz());
1775 AliESDRun *esdRun,*esdRunH;
1776 esdRun = (AliESDRun*)fesd->GetESDRun();
1777 esdRunH = (AliESDRun*)fhltesd->GetESDRun();
1778 esdRun->SetBeamEnergyIsSqrtSHalfGeV();
1779 esdRunH->SetBeamEnergyIsSqrtSHalfGeV();
1781 for (int ib=2;ib--;) for (int it=2;it--;) {
1782 esdRun->SetMeanIntensity(ib,it, fBeamInt[ib][it]);
1783 esdRunH->SetMeanIntensity(ib,it, fBeamInt[ib][it]);
1786 AliMagF* fld = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1787 if (fld) { // set info needed for field initialization
1788 fesd->SetCurrentL3(fld->GetCurrentSol());
1789 fesd->SetCurrentDip(fld->GetCurrentDip());
1790 fesd->SetBeamEnergy(fld->GetBeamEnergy());
1791 fesd->SetBeamType(fld->GetBeamTypeText());
1792 fesd->SetUniformBMap(fld->IsUniform());
1793 fesd->SetBInfoStored();
1795 fhltesd->SetCurrentL3(fld->GetCurrentSol());
1796 fhltesd->SetCurrentDip(fld->GetCurrentDip());
1797 fhltesd->SetBeamEnergy(fld->GetBeamEnergy());
1798 fhltesd->SetBeamType(fld->GetBeamTypeText());
1799 fhltesd->SetUniformBMap(fld->IsUniform());
1800 fhltesd->SetBInfoStored();
1803 // Set most probable pt, for B=0 tracking
1804 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1805 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1806 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1808 // Fill raw-data error log into the ESD
1809 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1812 if (fRunVertexFinder) {
1813 if (!RunVertexFinder(fesd)) {
1814 if (fStopOnError) {CleanUp(); return kFALSE;}
1818 // For Plane Efficiency: run the SPD trackleter
1819 if (fRunPlaneEff && fSPDTrackleter) {
1820 if (!RunSPDTrackleting(fesd)) {
1821 if (fStopOnError) {CleanUp(); return kFALSE;}
1826 if (!fRunTracking.IsNull()) {
1827 if (fRunMuonTracking) {
1828 if (!RunMuonTracking(fesd)) {
1829 if (fStopOnError) {CleanUp(); return kFALSE;}
1835 if (!fRunTracking.IsNull()) {
1836 if (!RunTracking(fesd,pid)) {
1837 if (fStopOnError) {CleanUp(); return kFALSE;}
1842 if (!fFillESD.IsNull()) {
1843 TString detectors=fFillESD;
1844 // run HLT first and on hltesd
1845 // ;-( IsSelected changes the string
1846 if (IsSelected("HLT", detectors) &&
1847 !FillESD(fhltesd, "HLT")) {
1848 if (fStopOnError) {CleanUp(); return kFALSE;}
1851 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1852 if (detectors.Contains("ALL")) {
1854 for (Int_t idet=0; idet<kNDetectors; ++idet){
1855 detectors += fgkDetectorName[idet];
1859 detectors.ReplaceAll("HLT", "");
1860 if (!FillESD(fesd, detectors)) {
1861 if (fStopOnError) {CleanUp(); return kFALSE;}
1868 if (fFillTriggerESD) {
1869 if (!FillTriggerESD(fesd)) {
1870 if (fStopOnError) {CleanUp(); return kFALSE;}
1873 // Always fill scalers
1874 if (!FillTriggerScalers(fesd)) {
1875 if (fStopOnError) {CleanUp(); return kFALSE;}
1882 // Propagate track to the beam pipe (if not already done by ITS)
1884 const Int_t ntracks = fesd->GetNumberOfTracks();
1885 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1888 UShort_t *selectedIdx=new UShort_t[ntracks];
1890 for (Int_t itrack=0; itrack<ntracks; itrack++){
1891 const Double_t kMaxStep = 1; //max step over the material
1894 AliESDtrack *track = fesd->GetTrack(itrack);
1895 if (!track) continue;
1897 AliExternalTrackParam *tpcTrack =
1898 (AliExternalTrackParam *)track->GetTPCInnerParam();
1902 PropagateTrackToBxByBz(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1905 Int_t n=trkArray.GetEntriesFast();
1906 selectedIdx[n]=track->GetID();
1907 trkArray.AddLast(tpcTrack);
1910 //Tracks refitted by ITS should already be at the SPD vertex
1911 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1914 PropagateTrackToBxByBz(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1915 Double_t x[3]; track->GetXYZ(x);
1916 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1917 track->RelateToVertexBxByBz(fesd->GetPrimaryVertexSPD(), b, kVeryBig);
1922 // Improve the reconstructed primary vertex position using the tracks
1924 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1925 if(fesd->GetPrimaryVertexSPD()) {
1926 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1927 if(vtitle.Contains("cosmics")) {
1928 runVertexFinderTracks=kFALSE;
1932 if (runVertexFinderTracks) {
1933 // TPC + ITS primary vertex
1934 ftVertexer->SetITSMode();
1935 ftVertexer->SetConstraintOff();
1936 // get cuts for vertexer from AliGRPRecoParam
1937 Bool_t constrSPD=kFALSE;
1939 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1940 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1941 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1942 ftVertexer->SetCuts(cutsVertexer);
1943 delete [] cutsVertexer; cutsVertexer = NULL;
1944 if(grpRecoParam->GetVertexerTracksConstraintITS()) {
1945 if(fDiamondProfile && fDiamondProfile->GetXRes()<kRadius){
1946 ftVertexer->SetVtxStart(fDiamondProfile); // apply constraint only if sigmax is smaller than the beam pipe radius
1948 if(fDiamondProfileSPD && fDiamondProfileSPD->GetXRes()<kRadius){
1949 ftVertexer->SetVtxStart(fDiamondProfileSPD);
1955 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1958 TString title=pvtx->GetTitle();
1959 title.Append("SPD");
1960 pvtx->SetTitle(title);
1962 if (pvtx->GetStatus()) {
1963 fesd->SetPrimaryVertexTracks(pvtx);
1964 for (Int_t i=0; i<ntracks; i++) {
1965 AliESDtrack *t = fesd->GetTrack(i);
1966 Double_t x[3]; t->GetXYZ(x);
1967 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1968 t->RelateToVertexBxByBz(pvtx, b, kVeryBig);
1971 delete pvtx; pvtx=NULL;
1974 // TPC-only primary vertex
1975 ftVertexer->SetTPCMode();
1976 ftVertexer->SetConstraintOff();
1977 // get cuts for vertexer from AliGRPRecoParam
1979 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1980 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1981 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1982 ftVertexer->SetCuts(cutsVertexer);
1983 delete [] cutsVertexer; cutsVertexer = NULL;
1984 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC()) {
1985 if(fDiamondProfileTPC->GetXRes()<kRadius) ftVertexer->SetVtxStart(fDiamondProfileTPC); // apply constraint only if sigmax is smaller than the beam pipe radius
1988 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1990 if (pvtx->GetStatus()) {
1991 fesd->SetPrimaryVertexTPC(pvtx);
1992 for (Int_t i=0; i<ntracks; i++) {
1993 AliESDtrack *t = fesd->GetTrack(i);
1994 Double_t x[3]; t->GetXYZ(x);
1995 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1996 t->RelateToVertexTPCBxByBz(pvtx, b, kVeryBig);
1999 delete pvtx; pvtx=NULL;
2003 delete[] selectedIdx;
2005 if(fDiamondProfile && fDiamondProfile->GetXRes()<kRadius) fesd->SetDiamond(fDiamondProfile);
2006 else fesd->SetDiamond(fDiamondProfileSPD);
2010 AliV0vertexer vtxer;
2011 vtxer.Tracks2V0vertices(fesd);
2013 if (fRunCascadeFinder) {
2015 AliCascadeVertexer cvtxer;
2016 cvtxer.V0sTracks2CascadeVertices(fesd);
2020 // RS run updated trackleter: since we want to mark the clusters used by tracks and also mark the
2021 // tracks interpreted as primary, this step should be done in the very end, when full
2022 // ESD info is available (particulalry, V0s)
2024 if (fRunMultFinder) {
2025 if (!RunMultFinder(fesd)) {
2026 if (fStopOnError) {CleanUp(); return kFALSE;}
2031 if (fCleanESD) CleanESD(fesd);
2033 if (fRunQA && IsInTasks(AliQAv1::kESDS)) {
2034 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2035 AliQAManager::QAManager()->RunOneEvent(fesd, fhltesd) ;
2038 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2039 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2040 if (qadm && IsInTasks(AliQAv1::kESDS))
2041 qadm->Exec(AliQAv1::kESDS, fesd);
2044 // copy HLT decision from HLTesd to esd
2045 // the most relevant information is stored in a reduced container in the esd,
2046 // while the full information can be found in the HLTesd
2047 TObject* pHLTSrc=fhltesd->FindListObject(AliESDHLTDecision::Name());
2048 TObject* pHLTTgt=fesd->FindListObject(AliESDHLTDecision::Name());
2049 if (pHLTSrc && pHLTTgt) {
2050 pHLTSrc->Copy(*pHLTTgt);
2053 if (fWriteESDfriend)
2054 fesd->GetESDfriend(fesdf);
2057 if (fWriteESDfriend) {
2059 Double_t rnd = gRandom->Rndm();
2060 if (fFractionFriends < rnd) {
2061 fesdf->~AliESDfriend();
2062 new (fesdf) AliESDfriend(); // Reset...
2063 fesdf->SetSkipBit(kTRUE);
2069 // Auto-save the ESD tree in case of prompt reco @P2
2070 if (fRawReader && fRawReader->UseAutoSaveESD()) {
2071 ftree->AutoSave("SaveSelf");
2072 if (fWriteESDfriend) ftreeF->AutoSave("SaveSelf");
2073 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
2074 if (friendfile) friendfile->Save();
2081 if (fRunAliEVE) RunAliEVE();
2085 if (fWriteESDfriend) {
2086 fesdf->~AliESDfriend();
2087 new (fesdf) AliESDfriend(); // Reset...
2090 gSystem->GetProcInfo(&procInfo);
2091 Long_t dMres=(procInfo.fMemResident-oldMres)/1024;
2092 Long_t dMvir=(procInfo.fMemVirtual-oldMvir)/1024;
2093 Float_t dCPU=procInfo.fCpuUser+procInfo.fCpuSys-oldCPU;
2094 aveDMres+=(dMres-aveDMres)/(iEvent-fFirstEvent+1);
2095 aveDMvir+=(dMvir-aveDMvir)/(iEvent-fFirstEvent+1);
2096 aveDCPU+=(dCPU-aveDCPU)/(iEvent-fFirstEvent+1);
2097 AliInfo(Form("======================= End Event %d: Res %d(%3d <%3d>) Vir %d(%3d <%3d>) CPU %5.2f <%5.2f> ===================",
2098 iEvent, procInfo.fMemResident/1024, dMres, aveDMres, procInfo.fMemVirtual/1024, dMvir, aveDMvir, dCPU, aveDCPU));
2099 oldMres=procInfo.fMemResident;
2100 oldMvir=procInfo.fMemVirtual;
2101 oldCPU=procInfo.fCpuUser+procInfo.fCpuSys;
2104 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2105 if (fReconstructor[iDet]) {
2106 fReconstructor[iDet]->SetRecoParam(NULL);
2107 fReconstructor[iDet]->SetEventInfo(NULL);
2109 if (fTracker[iDet]) fTracker[iDet]->SetEventInfo(NULL);
2112 if (fRunQA || fRunGlobalQA)
2113 AliQAManager::QAManager()->Increment() ;
2118 //_____________________________________________________________________________
2119 void AliReconstruction::SlaveTerminate()
2121 // Finalize the run on the slave side
2122 // Called after the exit
2123 // from the event loop
2124 AliCodeTimerAuto("",0);
2126 if (fIsNewRunLoader) { // galice.root didn't exist
2127 fRunLoader->WriteHeader("OVERWRITE");
2128 fRunLoader->CdGAFile();
2129 fRunLoader->Write(0, TObject::kOverwrite);
2132 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
2133 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
2135 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
2136 cdbMapCopy->SetOwner(1);
2137 cdbMapCopy->SetName("cdbMap");
2138 TIter iter(cdbMap->GetTable());
2141 while((pair = dynamic_cast<TPair*> (iter.Next()))){
2142 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
2143 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
2144 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
2147 TList *cdbListCopy = new TList();
2148 cdbListCopy->SetOwner(1);
2149 cdbListCopy->SetName("cdbList");
2151 TIter iter2(cdbList);
2154 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
2155 cdbListCopy->Add(new TObjString(id->ToString().Data()));
2158 ftree->GetUserInfo()->Add(cdbMapCopy);
2159 ftree->GetUserInfo()->Add(cdbListCopy);
2164 // we want to have only one tree version number
2165 ftree->Write(ftree->GetName(),TObject::kOverwrite);
2166 fhlttree->Write(fhlttree->GetName(),TObject::kOverwrite);
2168 if (fWriteESDfriend) {
2170 ftreeF->Write(ftreeF->GetName(),TObject::kOverwrite);
2173 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
2174 if (fRunPlaneEff && !FinishPlaneEff()) {
2175 AliWarning("Finish PlaneEff evaluation failed");
2178 // End of cycle for the in-loop
2180 if (fRunQA || fRunGlobalQA) {
2181 AliQAManager::QAManager()->EndOfCycle() ;
2183 !fProofOutputLocation.IsNull() &&
2184 fProofOutputArchive.IsNull() &&
2185 !fProofOutputDataset) {
2186 TString qaOutputFile(Form("%sMerged.%s.Data.root",
2187 fProofOutputLocation.Data(),
2188 AliQAv1::GetQADataFileName()));
2189 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",
2190 AliQAv1::GetQADataFileName()));
2191 qaProofFile->SetOutputFileName(qaOutputFile.Data());
2192 if (AliDebugLevel() > 0) qaProofFile->Dump();
2193 fOutput->Add(qaProofFile);
2194 MergeQA(qaProofFile->GetFileName());
2205 if (!fProofOutputFileName.IsNull() &&
2206 !fProofOutputLocation.IsNull() &&
2207 fProofOutputDataset &&
2208 !fProofOutputArchive.IsNull()) {
2209 TProofOutputFile *zipProofFile = new TProofOutputFile(fProofOutputFileName.Data(),
2211 fProofOutputLocation.Data());
2212 if (AliDebugLevel() > 0) zipProofFile->Dump();
2213 fOutput->Add(zipProofFile);
2214 TString fileList(fProofOutputArchive.Data());
2215 fileList.ReplaceAll(","," ");
2217 #if ROOT_SVN_REVISION >= 30174
2218 command.Form("zip -n root %s/%s %s",zipProofFile->GetDir(kTRUE),zipProofFile->GetFileName(),fileList.Data());
2220 command.Form("zip -n root %s/%s %s",zipProofFile->GetDir(),zipProofFile->GetFileName(),fileList.Data());
2222 AliInfo(Form("Executing: %s",command.Data()));
2223 gSystem->Exec(command.Data());
2228 //_____________________________________________________________________________
2229 void AliReconstruction::Terminate()
2231 // Create tags for the events in the ESD tree (the ESD tree is always present)
2232 // In case of empty events the tags will contain dummy values
2233 AliCodeTimerAuto("",0);
2235 // Do not call the ESD tag creator in case of PROOF-based reconstruction
2237 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
2238 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
2239 delete esdtagCreator;
2242 // Cleanup of CDB manager: cache and active storages!
2243 AliCDBManager::Instance()->ClearCache();
2246 //_____________________________________________________________________________
2247 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
2249 // run the local reconstruction
2251 static Int_t eventNr=0;
2252 AliCodeTimerAuto("",0)
2254 TString detStr = detectors;
2255 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2256 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2257 AliReconstructor* reconstructor = GetReconstructor(iDet);
2258 if (!reconstructor) continue;
2259 AliLoader* loader = fLoader[iDet];
2260 // Matthias April 2008: temporary fix to run HLT reconstruction
2261 // although the HLT loader is missing
2262 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2264 reconstructor->Reconstruct(fRawReader, NULL);
2267 reconstructor->Reconstruct(dummy, NULL);
2272 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2275 // conversion of digits
2276 if (fRawReader && reconstructor->HasDigitConversion()) {
2277 AliInfo(Form("converting raw data digits into root objects for %s",
2278 fgkDetectorName[iDet]));
2279 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2280 // fgkDetectorName[iDet]),0);
2281 loader->LoadDigits("update");
2282 loader->CleanDigits();
2283 loader->MakeDigitsContainer();
2284 TTree* digitsTree = loader->TreeD();
2285 reconstructor->ConvertDigits(fRawReader, digitsTree);
2286 loader->WriteDigits("OVERWRITE");
2287 loader->UnloadDigits();
2289 // local reconstruction
2290 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2291 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]),0);
2292 loader->LoadRecPoints("update");
2293 loader->CleanRecPoints();
2294 loader->MakeRecPointsContainer();
2295 TTree* clustersTree = loader->TreeR();
2296 if (fRawReader && !reconstructor->HasDigitConversion()) {
2297 reconstructor->Reconstruct(fRawReader, clustersTree);
2299 loader->LoadDigits("read");
2300 TTree* digitsTree = loader->TreeD();
2302 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2306 reconstructor->Reconstruct(digitsTree, clustersTree);
2307 if (fRunQA && IsInTasks(AliQAv1::kDIGITSR)) {
2308 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2309 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2312 loader->UnloadDigits();
2314 if (fRunQA && IsInTasks(AliQAv1::kRECPOINTS)) {
2315 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2316 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2318 loader->WriteRecPoints("OVERWRITE");
2319 loader->UnloadRecPoints();
2320 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2322 IsSelected("CTP", detStr);
2323 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2324 AliError(Form("the following detectors were not found: %s",
2332 //_____________________________________________________________________________
2333 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2335 // run the SPD trackleting (for SPD efficiency purpouses)
2337 AliCodeTimerAuto("",0)
2339 Double_t vtxPos[3] = {0, 0, 0};
2340 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2342 TArrayF mcVertex(3);
2344 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2345 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2346 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2349 const AliESDVertex *vertex = esd->GetVertex();
2351 AliWarning("Vertex not found");
2354 vertex->GetXYZ(vtxPos);
2355 vertex->GetSigmaXYZ(vtxErr);
2356 if (fSPDTrackleter) {
2357 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2360 fLoader[0]->LoadRecPoints("read");
2361 TTree* tree = fLoader[0]->TreeR();
2363 AliError("Can't get the ITS cluster tree");
2366 fSPDTrackleter->LoadClusters(tree);
2367 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2369 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2370 AliWarning("AliITSTrackleterSPDEff Clusters2Tracks failed");
2371 // fLoader[0]->UnloadRecPoints();
2374 //fSPDTrackleter->UnloadRecPoints();
2376 AliWarning("SPDTrackleter not available");
2382 //_____________________________________________________________________________
2383 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2385 // run the barrel tracking
2387 AliCodeTimerAuto("",0)
2389 AliVertexer *vertexer = CreateVertexer();
2390 if (!vertexer) return kFALSE;
2392 AliInfo(Form("running the ITS vertex finder: %s",vertexer->ClassName()));
2393 AliESDVertex* vertex = NULL;
2395 fLoader[0]->LoadRecPoints();
2396 TTree* cltree = fLoader[0]->TreeR();
2398 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2399 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2402 AliError("Can't get the ITS cluster tree");
2404 fLoader[0]->UnloadRecPoints();
2407 AliError("Can't get the ITS loader");
2410 AliWarning("Vertex not found");
2411 vertex = new AliESDVertex();
2412 vertex->SetName("default");
2415 vertex->SetName("reconstructed");
2420 vertex->GetXYZ(vtxPos);
2421 vertex->GetSigmaXYZ(vtxErr);
2423 esd->SetPrimaryVertexSPD(vertex);
2424 AliESDVertex *vpileup = NULL;
2425 Int_t novertices = 0;
2426 vpileup = vertexer->GetAllVertices(novertices);
2428 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2431 // if SPD multiplicity has been determined, it is stored in the ESD
2432 AliMultiplicity *mult = vertexer->GetMultiplicity();
2433 if(mult)esd->SetMultiplicity(mult);
2435 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2436 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2445 //_____________________________________________________________________________
2446 Bool_t AliReconstruction::RunMultFinder(AliESDEvent*& esd)
2448 // run the trackleter for multiplicity study
2450 AliCodeTimerAuto("",0)
2452 AliTrackleter *trackleter = CreateMultFinder();
2453 if (!trackleter) return kFALSE;
2455 AliInfo(Form("running the ITS multiplicity finder: %s",trackleter->ClassName()));
2458 fLoader[0]->LoadRecPoints();
2459 TTree* cltree = fLoader[0]->TreeR();
2461 trackleter->Reconstruct(esd,cltree);
2462 AliMultiplicity *mult = trackleter->GetMultiplicity();
2463 if(mult) esd->SetMultiplicity(mult);
2466 AliError("Can't get the ITS cluster tree");
2468 fLoader[0]->UnloadRecPoints();
2471 AliError("Can't get the ITS loader");
2479 //_____________________________________________________________________________
2480 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2482 // run the HLT barrel tracking
2484 AliCodeTimerAuto("",0)
2487 AliError("Missing runLoader!");
2491 AliInfo("running HLT tracking");
2493 // Get a pointer to the HLT reconstructor
2494 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2495 if (!reconstructor) return kFALSE;
2498 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2499 TString detName = fgkDetectorName[iDet];
2500 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2501 reconstructor->SetOption(detName.Data());
2502 AliTracker *tracker = reconstructor->CreateTracker();
2504 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2505 if (fStopOnError) return kFALSE;
2509 Double_t vtxErr[3]={0.005,0.005,0.010};
2510 const AliESDVertex *vertex = esd->GetVertex();
2511 vertex->GetXYZ(vtxPos);
2512 tracker->SetVertex(vtxPos,vtxErr);
2514 fLoader[iDet]->LoadRecPoints("read");
2515 TTree* tree = fLoader[iDet]->TreeR();
2517 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2520 tracker->LoadClusters(tree);
2522 if (tracker->Clusters2Tracks(esd) != 0) {
2523 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2527 tracker->UnloadClusters();
2535 //_____________________________________________________________________________
2536 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2538 // run the muon spectrometer tracking
2540 AliCodeTimerAuto("",0)
2543 AliError("Missing runLoader!");
2546 Int_t iDet = 7; // for MUON
2548 AliInfo("is running...");
2550 // Get a pointer to the MUON reconstructor
2551 AliReconstructor *reconstructor = GetReconstructor(iDet);
2552 if (!reconstructor) return kFALSE;
2555 TString detName = fgkDetectorName[iDet];
2556 AliDebug(1, Form("%s tracking", detName.Data()));
2557 AliTracker *tracker = reconstructor->CreateTracker();
2559 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2564 fLoader[iDet]->LoadRecPoints("read");
2566 tracker->LoadClusters(fLoader[iDet]->TreeR());
2568 Int_t rv = tracker->Clusters2Tracks(esd);
2570 fLoader[iDet]->UnloadRecPoints();
2572 tracker->UnloadClusters();
2578 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2586 //_____________________________________________________________________________
2587 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd,AliESDpid &PID)
2589 // run the barrel tracking
2590 static Int_t eventNr=0;
2591 AliCodeTimerAuto("",0)
2593 AliInfo("running tracking");
2595 // Set the event info which is used
2596 // by the trackers in order to obtain
2597 // information about read-out detectors,
2599 AliDebug(1, "Setting event info");
2600 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2601 if (!fTracker[iDet]) continue;
2602 fTracker[iDet]->SetEventInfo(&fEventInfo);
2605 //Fill the ESD with the T0 info (will be used by the TOF)
2606 if (fReconstructor[11] && fLoader[11]) {
2607 fLoader[11]->LoadRecPoints("READ");
2608 TTree *treeR = fLoader[11]->TreeR();
2610 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2614 // pass 1: TPC + ITS inwards
2615 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2616 if (!fTracker[iDet]) continue;
2617 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2620 fLoader[iDet]->LoadRecPoints("read");
2621 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2622 TTree* tree = fLoader[iDet]->TreeR();
2624 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2627 fTracker[iDet]->LoadClusters(tree);
2628 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2630 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2631 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2634 // preliminary PID in TPC needed by the ITS tracker
2636 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2637 PID.MakePID(esd,kTRUE);
2639 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2642 // pass 2: ALL backwards
2644 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2645 if (!fTracker[iDet]) continue;
2646 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2649 if (iDet > 1) { // all except ITS, TPC
2651 fLoader[iDet]->LoadRecPoints("read");
2652 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2653 tree = fLoader[iDet]->TreeR();
2655 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2658 fTracker[iDet]->LoadClusters(tree);
2659 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2663 if (iDet>1) // start filling residuals for the "outer" detectors
2665 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2666 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2668 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2669 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2670 if ( elem && (! elem->At(0)) ) {
2671 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2672 if (qadm) qadm->InitRecPointsForTracker() ;
2676 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2677 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2682 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2683 fTracker[iDet]->UnloadClusters();
2684 fLoader[iDet]->UnloadRecPoints();
2686 // updated PID in TPC needed by the ITS tracker -MI
2688 //GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2689 //AliESDpid::MakePID(esd);
2690 PID.MakePID(esd,kTRUE);
2692 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2694 //stop filling residuals for the "outer" detectors
2695 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2697 // pass 3: TRD + TPC + ITS refit inwards
2699 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2700 if (!fTracker[iDet]) continue;
2701 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2704 if (iDet<2) // start filling residuals for TPC and ITS
2706 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2707 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2709 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2710 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2711 if ( elem && (! elem->At(0)) ) {
2712 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2713 if (qadm) qadm->InitRecPointsForTracker() ;
2718 if (fTracker[iDet]->RefitInward(esd) != 0) {
2719 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2722 // run postprocessing
2723 if (fTracker[iDet]->PostProcess(esd) != 0) {
2724 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2727 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2730 // write space-points to the ESD in case alignment data output
2732 if (fWriteAlignmentData)
2733 WriteAlignmentData(esd);
2735 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2736 if (!fTracker[iDet]) continue;
2738 fTracker[iDet]->UnloadClusters();
2739 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2740 fLoader[iDet]->UnloadRecPoints();
2741 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2743 // stop filling residuals for TPC and ITS
2744 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2750 //_____________________________________________________________________________
2751 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2753 // Remove the data which are not needed for the physics analysis.
2756 Int_t nTracks=esd->GetNumberOfTracks();
2757 Int_t nV0s=esd->GetNumberOfV0s();
2759 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2761 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2762 Bool_t rc=esd->Clean(cleanPars);
2764 nTracks=esd->GetNumberOfTracks();
2765 nV0s=esd->GetNumberOfV0s();
2767 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2772 //_____________________________________________________________________________
2773 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2775 // fill the event summary data
2777 AliCodeTimerAuto("",0)
2778 static Int_t eventNr=0;
2779 TString detStr = detectors;
2781 AliSysInfo::AddStamp(Form("FillESDb%d",eventNr), -19,-19, eventNr);
2782 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2783 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2784 AliReconstructor* reconstructor = GetReconstructor(iDet);
2785 if (!reconstructor) continue;
2786 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2787 TTree* clustersTree = NULL;
2788 if (fLoader[iDet]) {
2789 fLoader[iDet]->LoadRecPoints("read");
2790 clustersTree = fLoader[iDet]->TreeR();
2791 if (!clustersTree) {
2792 AliError(Form("Can't get the %s clusters tree",
2793 fgkDetectorName[iDet]));
2794 if (fStopOnError) return kFALSE;
2797 if (fRawReader && !reconstructor->HasDigitConversion()) {
2798 reconstructor->FillESD(fRawReader, clustersTree, esd);
2800 TTree* digitsTree = NULL;
2801 if (fLoader[iDet]) {
2802 fLoader[iDet]->LoadDigits("read");
2803 digitsTree = fLoader[iDet]->TreeD();
2805 AliError(Form("Can't get the %s digits tree",
2806 fgkDetectorName[iDet]));
2807 if (fStopOnError) return kFALSE;
2810 reconstructor->FillESD(digitsTree, clustersTree, esd);
2811 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2813 if (fLoader[iDet]) {
2814 fLoader[iDet]->UnloadRecPoints();
2818 IsSelected("CTP", detStr);
2819 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2820 AliError(Form("the following detectors were not found: %s",
2822 if (fStopOnError) return kFALSE;
2824 AliSysInfo::AddStamp(Form("FillESDe%d",eventNr), -20,-20, eventNr);
2829 //_____________________________________________________________________________
2830 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2832 // Reads the trigger decision which is
2833 // stored in Trigger.root file and fills
2834 // the corresponding esd entries
2836 AliCodeTimerAuto("",0)
2838 AliInfo("Filling trigger information into the ESD");
2841 AliCTPRawStream input(fRawReader);
2842 if (!input.Next()) {
2843 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2846 if (esd->GetTriggerMask() != input.GetClassMask())
2847 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2848 input.GetClassMask(),esd->GetTriggerMask()));
2849 if (esd->GetOrbitNumber() != input.GetOrbitID())
2850 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2851 input.GetOrbitID(),esd->GetOrbitNumber()));
2852 if (esd->GetBunchCrossNumber() != input.GetBCID())
2853 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2854 input.GetBCID(),esd->GetBunchCrossNumber()));
2855 AliESDHeader* esdheader = esd->GetHeader();
2856 esdheader->SetL0TriggerInputs(input.GetL0Inputs());
2857 esdheader->SetL1TriggerInputs(input.GetL1Inputs());
2858 esdheader->SetL2TriggerInputs(input.GetL2Inputs());
2860 UInt_t orbit=input.GetOrbitID();
2861 for(Int_t i=0 ; i<input.GetNIRs() ; i++ )
2862 if(TMath::Abs(Int_t(orbit-(input.GetIR(i))->GetOrbit()))<=1){
2863 esdheader->AddTriggerIR(input.GetIR(i));
2869 //_____________________________________________________________________________
2870 Bool_t AliReconstruction::FillTriggerScalers(AliESDEvent*& esd)
2873 //fRunScalers->Print();
2874 if(fRunScalers && fRunScalers->CheckRunScalers()){
2875 AliTimeStamp* timestamp = new AliTimeStamp(esd->GetOrbitNumber(), esd->GetPeriodNumber(), esd->GetBunchCrossNumber());
2876 //AliTimeStamp* timestamp = new AliTimeStamp(10308000, 0, (ULong64_t)486238);
2877 AliESDHeader* esdheader = fesd->GetHeader();
2878 for(Int_t i=0;i<50;i++){
2879 if((1ull<<i) & esd->GetTriggerMask()){
2880 AliTriggerScalersESD* scalesd = fRunScalers->GetScalersForEventClass( timestamp, i+1);
2881 if(scalesd)esdheader->SetTriggerScalersRecord(scalesd);
2887 //_____________________________________________________________________________
2888 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2891 // Filling information from RawReader Header
2894 if (!fRawReader) return kFALSE;
2896 AliInfo("Filling information from RawReader Header");
2898 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2899 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2900 esd->SetPeriodNumber(fRawReader->GetPeriod());
2902 esd->SetTimeStamp(fRawReader->GetTimestamp());
2903 esd->SetEventType(fRawReader->GetType());
2909 //_____________________________________________________________________________
2910 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2912 // check whether detName is contained in detectors
2913 // if yes, it is removed from detectors
2915 // check if all detectors are selected
2916 if ((detectors.CompareTo("ALL") == 0) ||
2917 detectors.BeginsWith("ALL ") ||
2918 detectors.EndsWith(" ALL") ||
2919 detectors.Contains(" ALL ")) {
2924 // search for the given detector
2925 Bool_t result = kFALSE;
2926 if ((detectors.CompareTo(detName) == 0) ||
2927 detectors.BeginsWith(detName+" ") ||
2928 detectors.EndsWith(" "+detName) ||
2929 detectors.Contains(" "+detName+" ")) {
2930 detectors.ReplaceAll(detName, "");
2934 // clean up the detectors string
2935 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2936 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2937 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2942 //_____________________________________________________________________________
2943 Bool_t AliReconstruction::InitRunLoader()
2945 // get or create the run loader
2947 if (gAlice) delete gAlice;
2950 TFile *gafile = TFile::Open(fGAliceFileName.Data());
2951 // if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2952 if (gafile) { // galice.root exists
2956 // load all base libraries to get the loader classes
2957 TString libs = gSystem->GetLibraries();
2958 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2959 TString detName = fgkDetectorName[iDet];
2960 if (detName == "HLT") continue;
2961 if (libs.Contains("lib" + detName + "base.so")) continue;
2962 gSystem->Load("lib" + detName + "base.so");
2964 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2966 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2971 fRunLoader->CdGAFile();
2972 fRunLoader->LoadgAlice();
2974 //PH This is a temporary fix to give access to the kinematics
2975 //PH that is needed for the labels of ITS clusters
2976 fRunLoader->LoadHeader();
2977 fRunLoader->LoadKinematics();
2979 } else { // galice.root does not exist
2981 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2983 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2984 AliConfig::GetDefaultEventFolderName(),
2987 AliError(Form("could not create run loader in file %s",
2988 fGAliceFileName.Data()));
2992 fIsNewRunLoader = kTRUE;
2993 fRunLoader->MakeTree("E");
2995 if (fNumberOfEventsPerFile > 0)
2996 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2998 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
3004 //_____________________________________________________________________________
3005 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
3007 // get the reconstructor object and the loader for a detector
3009 if (fReconstructor[iDet]) {
3010 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
3011 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
3012 fReconstructor[iDet]->SetRecoParam(par);
3013 fReconstructor[iDet]->SetRunInfo(fRunInfo);
3015 return fReconstructor[iDet];
3018 // load the reconstructor object
3019 TPluginManager* pluginManager = gROOT->GetPluginManager();
3020 TString detName = fgkDetectorName[iDet];
3021 TString recName = "Ali" + detName + "Reconstructor";
3023 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
3025 AliReconstructor* reconstructor = NULL;
3026 // first check if a plugin is defined for the reconstructor
3027 TPluginHandler* pluginHandler =
3028 pluginManager->FindHandler("AliReconstructor", detName);
3029 // if not, add a plugin for it
3030 if (!pluginHandler) {
3031 AliDebug(1, Form("defining plugin for %s", recName.Data()));
3032 TString libs = gSystem->GetLibraries();
3033 if (libs.Contains("lib" + detName + "base.so") ||
3034 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
3035 pluginManager->AddHandler("AliReconstructor", detName,
3036 recName, detName + "rec", recName + "()");
3038 pluginManager->AddHandler("AliReconstructor", detName,
3039 recName, detName, recName + "()");
3041 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
3043 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
3044 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
3046 if (reconstructor) {
3047 TObject* obj = fOptions.FindObject(detName.Data());
3048 if (obj) reconstructor->SetOption(obj->GetTitle());
3049 reconstructor->SetRunInfo(fRunInfo);
3050 reconstructor->Init();
3051 fReconstructor[iDet] = reconstructor;
3054 // get or create the loader
3055 if (detName != "HLT") {
3056 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
3057 if (!fLoader[iDet]) {
3058 AliConfig::Instance()
3059 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
3061 // first check if a plugin is defined for the loader
3063 pluginManager->FindHandler("AliLoader", detName);
3064 // if not, add a plugin for it
3065 if (!pluginHandler) {
3066 TString loaderName = "Ali" + detName + "Loader";
3067 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
3068 pluginManager->AddHandler("AliLoader", detName,
3069 loaderName, detName + "base",
3070 loaderName + "(const char*, TFolder*)");
3071 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
3073 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
3075 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
3076 fRunLoader->GetEventFolder());
3078 if (!fLoader[iDet]) { // use default loader
3079 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
3081 if (!fLoader[iDet]) {
3082 AliWarning(Form("couldn't get loader for %s", detName.Data()));
3083 if (fStopOnError) return NULL;
3085 fRunLoader->AddLoader(fLoader[iDet]);
3086 fRunLoader->CdGAFile();
3087 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
3088 fRunLoader->Write(0, TObject::kOverwrite);
3093 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
3094 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
3095 reconstructor->SetRecoParam(par);
3096 reconstructor->SetRunInfo(fRunInfo);
3098 return reconstructor;
3101 //_____________________________________________________________________________
3102 AliVertexer* AliReconstruction::CreateVertexer()
3104 // create the vertexer
3105 // Please note that the caller is the owner of the
3108 AliVertexer* vertexer = NULL;
3109 AliReconstructor* itsReconstructor = GetReconstructor(0);
3110 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) || fRunTracking.Contains("ITS"))) {
3111 vertexer = itsReconstructor->CreateVertexer();
3114 AliWarning("couldn't create a vertexer for ITS");
3120 //_____________________________________________________________________________
3121 AliTrackleter* AliReconstruction::CreateMultFinder()
3123 // create the ITS trackleter for mult. estimation
3124 // Please note that the caller is the owner of the
3127 AliTrackleter* trackleter = NULL;
3128 AliReconstructor* itsReconstructor = GetReconstructor(0);
3129 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) || fRunTracking.Contains("ITS"))) {
3130 trackleter = itsReconstructor->CreateMultFinder();
3133 AliWarning("couldn't create a trackleter for ITS");
3139 //_____________________________________________________________________________
3140 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
3142 // create the trackers
3143 AliInfo("Creating trackers");
3145 TString detStr = detectors;
3146 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3147 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3148 AliReconstructor* reconstructor = GetReconstructor(iDet);
3149 if (!reconstructor) continue;
3150 TString detName = fgkDetectorName[iDet];
3151 if (detName == "HLT") {
3152 fRunHLTTracking = kTRUE;
3155 if (detName == "MUON") {
3156 fRunMuonTracking = kTRUE;
3160 fTracker[iDet] = reconstructor->CreateTracker();
3161 if (!fTracker[iDet] && (iDet < 7)) {
3162 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
3163 if (fStopOnError) return kFALSE;
3165 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
3171 //_____________________________________________________________________________
3172 void AliReconstruction::CleanUp()
3174 // delete trackers and the run loader and close and delete the file
3176 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3177 delete fReconstructor[iDet];
3178 fReconstructor[iDet] = NULL;
3179 fLoader[iDet] = NULL;
3180 delete fTracker[iDet];
3181 fTracker[iDet] = NULL;
3186 delete fSPDTrackleter;
3187 fSPDTrackleter = NULL;
3196 delete fParentRawReader;
3197 fParentRawReader=NULL;
3205 if (AliQAManager::QAManager())
3206 AliQAManager::QAManager()->ShowQA() ;
3207 AliQAManager::Destroy() ;
3211 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
3213 // Write space-points which are then used in the alignment procedures
3214 // For the moment only ITS, TPC, TRD and TOF
3216 Int_t ntracks = esd->GetNumberOfTracks();
3217 for (Int_t itrack = 0; itrack < ntracks; itrack++)
3219 AliESDtrack *track = esd->GetTrack(itrack);
3222 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
3223 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
3224 nsp += (iDet==GetDetIndex("TRD")) ? track->GetTRDntracklets():track->GetNcls(iDet);
3226 if (iDet==GetDetIndex("ITS")) { // ITS "extra" clusters
3227 track->GetClusters(iDet,idx);
3228 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
3233 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
3234 track->SetTrackPointArray(sp);
3236 for (Int_t iDet = 5; iDet >= 0; iDet--) {
3237 AliTracker *tracker = fTracker[iDet];
3238 if (!tracker) continue;
3239 Int_t nspdet = (iDet==GetDetIndex("TRD")) ? track->GetTRDtracklets(idx):track->GetClusters(iDet,idx);
3241 if (iDet==GetDetIndex("ITS")) // ITS "extra" clusters
3242 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
3244 if (nspdet <= 0) continue;
3248 while (isp2 < nspdet) {
3249 Bool_t isvalid=kTRUE;
3251 Int_t index=idx[isp++];
3252 if (index < 0) continue;
3254 TString dets = fgkDetectorName[iDet];
3255 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
3256 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
3257 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
3258 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
3259 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
3261 isvalid = tracker->GetTrackPoint(index,p);
3264 if (!isvalid) continue;
3265 if (iDet==GetDetIndex("ITS") && (isp-1)>=6) p.SetExtra();
3266 sp->AddPoint(isptrack,&p); isptrack++;
3273 //_____________________________________________________________________________
3274 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
3276 // The method reads the raw-data error log
3277 // accumulated within the rawReader.
3278 // It extracts the raw-data errors related to
3279 // the current event and stores them into
3280 // a TClonesArray inside the esd object.
3282 if (!fRawReader) return;
3284 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
3286 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
3288 if (iEvent != log->GetEventNumber()) continue;
3290 esd->AddRawDataErrorLog(log);
3295 //_____________________________________________________________________________
3296 // void AliReconstruction::CheckQA()
3298 // check the QA of SIM for this run and remove the detectors
3299 // with status Fatal
3301 // TString newRunLocalReconstruction ;
3302 // TString newRunTracking ;
3303 // TString newFillESD ;
3305 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
3306 // TString detName(AliQAv1::GetDetName(iDet)) ;
3307 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
3308 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
3309 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
3310 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
3312 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
3313 // fRunLocalReconstruction.Contains("ALL") ) {
3314 // newRunLocalReconstruction += detName ;
3315 // newRunLocalReconstruction += " " ;
3317 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
3318 // fRunTracking.Contains("ALL") ) {
3319 // newRunTracking += detName ;
3320 // newRunTracking += " " ;
3322 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
3323 // fFillESD.Contains("ALL") ) {
3324 // newFillESD += detName ;
3325 // newFillESD += " " ;
3329 // fRunLocalReconstruction = newRunLocalReconstruction ;
3330 // fRunTracking = newRunTracking ;
3331 // fFillESD = newFillESD ;
3334 //_____________________________________________________________________________
3335 Int_t AliReconstruction::GetDetIndex(const char* detector)
3337 // return the detector index corresponding to detector
3339 for (index = 0; index < kNDetectors ; index++) {
3340 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3345 //_____________________________________________________________________________
3346 Bool_t AliReconstruction::FinishPlaneEff() {
3348 // Here execute all the necessary operationis, at the end of the tracking phase,
3349 // in case that evaluation of PlaneEfficiencies was required for some detector.
3350 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3352 // This Preliminary version works only FOR ITS !!!!!
3353 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3356 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3359 TString detStr = fLoadCDB;
3360 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3361 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3362 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3363 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
3364 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3365 TString name=planeeff->GetName();
3367 TFile* pefile = TFile::Open(name, "RECREATE");
3368 ret=(Bool_t)planeeff->Write();
3370 if(planeeff->GetCreateHistos()) {
3371 TString hname=planeeff->GetName();
3372 hname+="Histo.root";
3373 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3376 if(fSPDTrackleter) {
3377 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3378 TString name="AliITSPlaneEffSPDtracklet.root";
3379 TFile* pefile = TFile::Open(name, "RECREATE");
3380 ret=(Bool_t)planeeff->Write();
3382 AliESDEvent *dummy=NULL;
3383 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3388 //_____________________________________________________________________________
3389 Bool_t AliReconstruction::InitPlaneEff() {
3391 // Here execute all the necessary operations, before of the tracking phase,
3392 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3393 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3394 // which should be updated/recalculated.
3396 // This Preliminary version will work only FOR ITS !!!!!
3397 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3400 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3403 fSPDTrackleter = NULL;
3404 TString detStr = fLoadCDB;
3405 if (IsSelected(fgkDetectorName[0], detStr)) {
3406 AliReconstructor* itsReconstructor = GetReconstructor(0);
3407 if (itsReconstructor) {
3408 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3410 if (fSPDTrackleter) {
3411 AliInfo("Trackleter for SPD has been created");
3417 //_____________________________________________________________________________
3418 Bool_t AliReconstruction::InitAliEVE()
3420 // This method should be called only in case
3421 // AliReconstruction is run
3422 // within the alieve environment.
3423 // It will initialize AliEVE in a way
3424 // so that it can visualize event processed
3425 // by AliReconstruction.
3426 // The return flag shows whenever the
3427 // AliEVE initialization was successful or not.
3430 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3431 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3432 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3434 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3435 gROOT->ProcessLine("alieve_online_init()");
3440 //_____________________________________________________________________________
3441 void AliReconstruction::RunAliEVE()
3443 // Runs AliEVE visualisation of
3444 // the current event.
3445 // Should be executed only after
3446 // successful initialization of AliEVE.
3448 AliInfo("Running AliEVE...");
3449 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)0x%lx,(AliRawReader*)0x%lx,(AliESDEvent*)0x%lx,(AliESDfriend*)0x%lx);",fRunLoader,fRawReader,fesd,fesdf));
3453 //_____________________________________________________________________________
3454 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3456 // Allows to run QA for a selected set of detectors
3457 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3458 // all selected detectors run the same selected tasks
3460 if (!detAndAction.Contains(":")) {
3461 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3465 Int_t colon = detAndAction.Index(":") ;
3466 fQADetectors = detAndAction(0, colon) ;
3467 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3468 if (fQATasks.Contains("ALL") ) {
3469 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3471 fQATasks.ToUpper() ;
3473 if ( fQATasks.Contains("RAW") )
3474 tempo = Form("%d ", AliQAv1::kRAWS) ;
3475 if ( fQATasks.Contains("DIGIT") )
3476 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3477 if ( fQATasks.Contains("RECPOINT") )
3478 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3479 if ( fQATasks.Contains("ESD") )
3480 tempo += Form("%d ", AliQAv1::kESDS) ;
3482 if (fQATasks.IsNull()) {
3483 AliInfo("No QA requested\n") ;
3488 TString tempo(fQATasks) ;
3489 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3490 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3491 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3492 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3493 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3498 //_____________________________________________________________________________
3499 Bool_t AliReconstruction::InitRecoParams()
3501 // The method accesses OCDB and retrieves all
3502 // the available reco-param objects from there.
3504 Bool_t isOK = kTRUE;
3506 if (fRecoParam.GetDetRecoParamArray(kNDetectors)) {
3507 AliInfo("Using custom GRP reconstruction parameters");
3510 AliInfo("Loading GRP reconstruction parameter objects");
3512 AliCDBPath path("GRP","Calib","RecoParam");
3513 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3515 AliWarning("Couldn't find GRP RecoParam entry in OCDB");
3519 TObject *recoParamObj = entry->GetObject();
3520 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3521 // GRP has a normal TobjArray of AliDetectorRecoParam objects
3522 // Registering them in AliRecoParam
3523 fRecoParam.AddDetRecoParamArray(kNDetectors,dynamic_cast<TObjArray*>(recoParamObj));
3525 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3526 // GRP has only onse set of reco parameters
3527 // Registering it in AliRecoParam
3528 AliInfo("Single set of GRP reconstruction parameters found");
3529 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3530 fRecoParam.AddDetRecoParam(kNDetectors,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3533 AliError("No valid GRP RecoParam object found in the OCDB");
3540 TString detStr = fLoadCDB;
3541 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3543 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3545 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3546 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3550 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3552 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3553 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3555 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3559 TObject *recoParamObj = entry->GetObject();
3560 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3561 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3562 // Registering them in AliRecoParam
3563 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3565 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3566 // The detector has only onse set of reco parameters
3567 // Registering it in AliRecoParam
3568 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3569 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3570 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3573 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3577 // FIX ME: We have to disable the unloading of reco-param CDB
3578 // entries because QA framework is using them. Has to be fix in
3579 // a way that the QA takes the objects already constructed in
3581 // AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3585 if (AliDebugLevel() > 0) fRecoParam.Print();
3590 //_____________________________________________________________________________
3591 Bool_t AliReconstruction::GetEventInfo()
3593 // Fill the event info object
3595 AliCodeTimerAuto("",0)
3597 AliCentralTrigger *aCTP = NULL;
3599 fEventInfo.SetEventType(fRawReader->GetType());
3601 ULong64_t mask = fRawReader->GetClassMask();
3602 fEventInfo.SetTriggerMask(mask);
3603 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3604 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3606 aCTP = new AliCentralTrigger();
3607 TString configstr("");
3608 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3609 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3613 aCTP->SetClassMask(mask);
3614 aCTP->SetClusterMask(clmask);
3617 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3619 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3620 aCTP = fRunLoader->GetTrigger();
3621 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3622 // get inputs from actp - just get
3623 AliESDHeader* esdheader = fesd->GetHeader();
3624 esdheader->SetL0TriggerInputs(aCTP->GetL0TriggerInputs());
3625 esdheader->SetL1TriggerInputs(aCTP->GetL1TriggerInputs());
3626 esdheader->SetL2TriggerInputs(aCTP->GetL2TriggerInputs());
3627 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3630 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3635 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3637 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3638 if (fRawReader) delete aCTP;
3642 UChar_t clustmask = 0;
3644 ULong64_t trmask = fEventInfo.GetTriggerMask();
3645 const TObjArray& classesArray = config->GetClasses();
3646 Int_t nclasses = classesArray.GetEntriesFast();
3647 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3648 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3649 if (trclass && trclass->GetMask()>0) {
3650 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3651 fesd->SetTriggerClass(trclass->GetName(),trindex);
3652 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3653 if (trmask & (1ull << trindex)) {
3655 trclasses += trclass->GetName();
3657 clustmask |= trclass->GetCluster()->GetClusterMask();
3661 fEventInfo.SetTriggerClasses(trclasses);
3663 // Write names of active trigger inputs in ESD Header
3664 const TObjArray& inputsArray = config->GetInputs();
3665 Int_t ninputs = inputsArray.GetEntriesFast();
3666 for( Int_t iinput=0; iinput < ninputs; iinput++ ) {
3667 AliTriggerInput* trginput = (AliTriggerInput*)inputsArray.At(iinput);
3668 if (trginput && trginput->GetMask()>0) {
3669 Int_t inputIndex = (Int_t)TMath::Nint(TMath::Log2(trginput->GetMask()));
3670 AliESDHeader* headeresd = fesd->GetHeader();
3671 Int_t trglevel = (Int_t)trginput->GetLevel();
3672 if (trglevel == 0) headeresd->SetActiveTriggerInputs(trginput->GetInputName(), inputIndex);
3673 if (trglevel == 1) headeresd->SetActiveTriggerInputs(trginput->GetInputName(), inputIndex+24);
3674 if (trglevel == 2) headeresd->SetActiveTriggerInputs(trginput->GetInputName(), inputIndex+48);
3678 // Set the information in ESD
3679 fesd->SetTriggerMask(trmask);
3680 fesd->SetTriggerCluster(clustmask);
3682 if (!aCTP->CheckTriggeredDetectors()) {
3683 if (fRawReader) delete aCTP;
3687 if (fRawReader) delete aCTP;
3689 // We have to fill also the HLT decision here!!
3695 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3697 // Match the detector list found in the rec.C or the default 'ALL'
3698 // to the list found in the GRP (stored there by the shuttle PP which
3699 // gets the information from ECS)
3700 static TString resultList;
3701 TString detList = detectorList;
3705 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3706 if ((detectorMask >> iDet) & 0x1) {
3707 TString det = AliDAQ::OfflineModuleName(iDet);
3708 if ((detList.CompareTo("ALL") == 0) ||
3709 ((detList.BeginsWith("ALL ") ||
3710 detList.EndsWith(" ALL") ||
3711 detList.Contains(" ALL ")) &&
3712 !(detList.BeginsWith("-"+det+" ") ||
3713 detList.EndsWith(" -"+det) ||
3714 detList.Contains(" -"+det+" "))) ||
3715 (detList.CompareTo(det) == 0) ||
3716 detList.BeginsWith(det+" ") ||
3717 detList.EndsWith(" "+det) ||
3718 detList.Contains( " "+det+" " )) {
3719 if (!resultList.EndsWith(det + " ")) {
3728 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3729 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3730 if ((detList.CompareTo("ALL") == 0) ||
3731 ((detList.BeginsWith("ALL ") ||
3732 detList.EndsWith(" ALL") ||
3733 detList.Contains(" ALL ")) &&
3734 !(detList.BeginsWith("-"+hltDet+" ") ||
3735 detList.EndsWith(" -"+hltDet) ||
3736 detList.Contains(" -"+hltDet+" "))) ||
3737 (detList.CompareTo(hltDet) == 0) ||
3738 detList.BeginsWith(hltDet+" ") ||
3739 detList.EndsWith(" "+hltDet) ||
3740 detList.Contains( " "+hltDet+" " )) {
3741 resultList += hltDet;
3745 return resultList.Data();
3749 //______________________________________________________________________________
3750 void AliReconstruction::Abort(const char *method, EAbort what)
3752 // Abort processing. If what = kAbortProcess, the Process() loop will be
3753 // aborted. If what = kAbortFile, the current file in a chain will be
3754 // aborted and the processing will continue with the next file, if there
3755 // is no next file then Process() will be aborted. Abort() can also be
3756 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3757 // the SlaveTerminate() and Terminate() are always called. The abort flag
3758 // can be checked in these methods using GetAbort().
3760 // The method is overwritten in AliReconstruction for better handling of
3761 // reco specific errors
3763 if (!fStopOnError) return;
3767 TString whyMess = method;
3768 whyMess += " failed! Aborting...";
3770 AliError(whyMess.Data());
3773 TString mess = "Abort";
3774 if (fAbort == kAbortProcess)
3775 mess = "AbortProcess";
3776 else if (fAbort == kAbortFile)
3779 Info(mess, whyMess.Data());
3782 //______________________________________________________________________________
3783 Bool_t AliReconstruction::ProcessEvent(void* event)
3785 // Method that is used in case the event loop
3786 // is steered from outside, for example by AMORE
3787 // 'event' is a pointer to the DATE event in the memory
3789 if (fRawReader) delete fRawReader;
3790 fRawReader = new AliRawReaderDate(event);
3791 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
3798 //______________________________________________________________________________
3799 Bool_t AliReconstruction::ParseOutput()
3801 // The method parses the output file
3802 // location string in order to steer
3803 // properly the selector
3805 TPMERegexp re1("(\\w+\\.zip#\\w+\\.root):([,*\\w+\\.root,*]+)@dataset://(\\w++)");
3806 TPMERegexp re2("(\\w+\\.root)?@?dataset://(\\w++)");
3808 if (re1.Match(fESDOutput) == 4) {
3809 // root archive with output files stored and regustered
3811 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",re1[1].Data()));
3812 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re1[3].Data()));
3813 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3814 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_ARCHIVE",re1[2].Data()));
3815 AliInfo(Form("%s files will be stored within %s in dataset %s",
3820 else if (re2.Match(fESDOutput) == 3) {
3821 // output file stored and registered
3823 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",(re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data()));
3824 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re2[2].Data()));
3825 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3826 AliInfo(Form("%s will be stored in dataset %s",
3827 (re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data(),
3831 if (fESDOutput.IsNull()) {
3832 // Output location not given.
3833 // Assuming xrootd has been already started and
3834 // the output file has to be sent back
3835 // to the client machine
3836 TString esdUrl(Form("root://%s/%s/",
3837 TUrl(gSystem->HostName()).GetHostFQDN(),
3839 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE","AliESDs.root"));
3840 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",esdUrl.Data()));
3841 AliInfo(Form("AliESDs.root will be stored in %s",
3845 // User specified an output location.
3846 // Ones has just to parse it here
3847 TUrl outputUrl(fESDOutput.Data());
3848 TString outputFile(gSystem->BaseName(outputUrl.GetFile()));
3849 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.IsNull() ? "AliESDs.root" : outputFile.Data()));
3850 TString outputLocation(outputUrl.GetUrl());
3851 outputLocation.ReplaceAll(outputFile.Data(),"");
3852 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",outputLocation.Data()));
3853 AliInfo(Form("%s will be stored in %s",
3854 outputFile.IsNull() ? "AliESDs.root" : outputFile.Data(),
3855 outputLocation.Data()));