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),
255 fParentRawReader(NULL),
259 fSPDTrackleter(NULL),
261 fDiamondProfileSPD(NULL),
262 fDiamondProfile(NULL),
263 fDiamondProfileTPC(NULL),
264 fListOfCosmicTriggers(NULL),
268 fAlignObjArray(NULL),
272 fInitCDBCalled(kFALSE),
273 fSetRunNumberFromDataCalled(kFALSE),
278 fSameQACycle(kFALSE),
279 fInitQACalled(kFALSE),
280 fWriteQAExpertData(kTRUE),
281 fRunPlaneEff(kFALSE),
292 fIsNewRunLoader(kFALSE),
301 // create reconstruction object with default parameters
304 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
305 fReconstructor[iDet] = NULL;
306 fLoader[iDet] = NULL;
307 fTracker[iDet] = NULL;
309 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
310 fQACycles[iDet] = 999999 ;
311 fQAWriteExpert[iDet] = kFALSE ;
313 fBeamInt[0][0]=fBeamInt[0][1]=fBeamInt[1][0]=fBeamInt[1][1] = -1;
318 //_____________________________________________________________________________
319 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
321 fRunVertexFinder(rec.fRunVertexFinder),
322 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
323 fRunHLTTracking(rec.fRunHLTTracking),
324 fRunMuonTracking(rec.fRunMuonTracking),
325 fRunV0Finder(rec.fRunV0Finder),
326 fRunCascadeFinder(rec.fRunCascadeFinder),
327 fRunMultFinder(rec.fRunMultFinder),
328 fStopOnError(rec.fStopOnError),
329 fWriteAlignmentData(rec.fWriteAlignmentData),
330 fWriteESDfriend(rec.fWriteESDfriend),
331 fFillTriggerESD(rec.fFillTriggerESD),
333 fCleanESD(rec.fCleanESD),
334 fV0DCAmax(rec.fV0DCAmax),
335 fV0CsPmin(rec.fV0CsPmin),
339 fRunLocalReconstruction(rec.fRunLocalReconstruction),
340 fRunTracking(rec.fRunTracking),
341 fFillESD(rec.fFillESD),
342 fLoadCDB(rec.fLoadCDB),
343 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
344 fGAliceFileName(rec.fGAliceFileName),
345 fRawInput(rec.fRawInput),
346 fESDOutput(rec.fESDOutput),
347 fProofOutputFileName(rec.fProofOutputFileName),
348 fProofOutputLocation(rec.fProofOutputLocation),
349 fProofOutputDataset(rec.fProofOutputDataset),
350 fProofOutputArchive(rec.fProofOutputArchive),
351 fEquipIdMap(rec.fEquipIdMap),
352 fFirstEvent(rec.fFirstEvent),
353 fLastEvent(rec.fLastEvent),
354 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
355 fFractionFriends(rec.fFractionFriends),
357 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
358 fLoadAlignData(rec.fLoadAlignData),
359 fUseHLTData(rec.fUseHLTData),
363 fCTPTimeParams(NULL),
368 fParentRawReader(NULL),
370 fRecoParam(rec.fRecoParam),
372 fSPDTrackleter(NULL),
374 fDiamondProfileSPD(rec.fDiamondProfileSPD),
375 fDiamondProfile(rec.fDiamondProfile),
376 fDiamondProfileTPC(rec.fDiamondProfileTPC),
377 fListOfCosmicTriggers(NULL),
381 fAlignObjArray(rec.fAlignObjArray),
382 fCDBUri(rec.fCDBUri),
383 fQARefUri(rec.fQARefUri),
385 fInitCDBCalled(rec.fInitCDBCalled),
386 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
387 fQADetectors(rec.fQADetectors),
388 fQATasks(rec.fQATasks),
390 fRunGlobalQA(rec.fRunGlobalQA),
391 fSameQACycle(rec.fSameQACycle),
392 fInitQACalled(rec.fInitQACalled),
393 fWriteQAExpertData(rec.fWriteQAExpertData),
394 fRunPlaneEff(rec.fRunPlaneEff),
405 fIsNewRunLoader(rec.fIsNewRunLoader),
416 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
417 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
419 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
420 fReconstructor[iDet] = NULL;
421 fLoader[iDet] = NULL;
422 fTracker[iDet] = NULL;
425 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
426 fQACycles[iDet] = rec.fQACycles[iDet];
427 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
430 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
431 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
434 for (int i=2;i--;) for (int j=2;j--;) fBeamInt[i][j] = rec.fBeamInt[i][j];
438 //_____________________________________________________________________________
439 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
441 // assignment operator
442 // Used in PROOF mode
443 // Be very careful while modifing it!
444 // Simple rules to follow:
445 // for persistent data members - use their assignment operators
446 // for non-persistent ones - do nothing or take the default values from constructor
447 // TSelector members should not be touched
448 if(&rec == this) return *this;
450 fRunVertexFinder = rec.fRunVertexFinder;
451 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
452 fRunHLTTracking = rec.fRunHLTTracking;
453 fRunMuonTracking = rec.fRunMuonTracking;
454 fRunV0Finder = rec.fRunV0Finder;
455 fRunCascadeFinder = rec.fRunCascadeFinder;
456 fRunMultFinder = rec.fRunMultFinder;
457 fStopOnError = rec.fStopOnError;
458 fWriteAlignmentData = rec.fWriteAlignmentData;
459 fWriteESDfriend = rec.fWriteESDfriend;
460 fFillTriggerESD = rec.fFillTriggerESD;
462 fCleanESD = rec.fCleanESD;
463 fV0DCAmax = rec.fV0DCAmax;
464 fV0CsPmin = rec.fV0CsPmin;
468 fRunLocalReconstruction = rec.fRunLocalReconstruction;
469 fRunTracking = rec.fRunTracking;
470 fFillESD = rec.fFillESD;
471 fLoadCDB = rec.fLoadCDB;
472 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
473 fGAliceFileName = rec.fGAliceFileName;
474 fRawInput = rec.fRawInput;
475 fESDOutput = rec.fESDOutput;
476 fProofOutputFileName = rec.fProofOutputFileName;
477 fProofOutputLocation = rec.fProofOutputLocation;
478 fProofOutputDataset = rec.fProofOutputDataset;
479 fProofOutputArchive = rec.fProofOutputArchive;
480 fEquipIdMap = rec.fEquipIdMap;
481 fFirstEvent = rec.fFirstEvent;
482 fLastEvent = rec.fLastEvent;
483 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
484 fFractionFriends = rec.fFractionFriends;
486 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
487 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
490 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
491 fLoadAlignData = rec.fLoadAlignData;
492 fUseHLTData = rec.fUseHLTData;
494 delete fRunInfo; fRunInfo = NULL;
495 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
497 fEventInfo = rec.fEventInfo;
499 delete fRunScalers; fRunScalers = NULL;
500 if (rec.fRunScalers) fRunScalers = new AliTriggerRunScalers(*rec.fRunScalers);
502 delete fCTPTimeParams; fCTPTimeParams = NULL;
503 if (rec.fCTPTimeParams) fCTPTimeParams = new AliCTPTimeParams(*rec.fCTPTimeParams);
504 delete fCTPTimeAlign; fCTPTimeAlign = NULL;
505 if (rec.fCTPTimeAlign) fCTPTimeAlign = new AliCTPTimeParams(*rec.fCTPTimeAlign);
509 fParentRawReader = NULL;
511 fRecoParam = rec.fRecoParam;
513 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
514 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
515 delete fLoader[iDet]; fLoader[iDet] = NULL;
516 delete fTracker[iDet]; fTracker[iDet] = NULL;
519 for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
520 fQACycles[iDet] = rec.fQACycles[iDet];
521 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
524 delete fSPDTrackleter; fSPDTrackleter = NULL;
526 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
527 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
528 delete fDiamondProfile; fDiamondProfile = NULL;
529 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
530 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
531 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
533 delete fListOfCosmicTriggers; fListOfCosmicTriggers = NULL;
534 if (rec.fListOfCosmicTriggers) fListOfCosmicTriggers = (THashTable*)((rec.fListOfCosmicTriggers)->Clone());
536 delete fGRPData; fGRPData = NULL;
537 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
538 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
540 delete fAlignObjArray; fAlignObjArray = NULL;
543 fQARefUri = rec.fQARefUri;
544 fSpecCDBUri.Delete();
545 fInitCDBCalled = rec.fInitCDBCalled;
546 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
547 fQADetectors = rec.fQADetectors;
548 fQATasks = rec.fQATasks;
550 fRunGlobalQA = rec.fRunGlobalQA;
551 fSameQACycle = rec.fSameQACycle;
552 fInitQACalled = rec.fInitQACalled;
553 fWriteQAExpertData = rec.fWriteQAExpertData;
554 fRunPlaneEff = rec.fRunPlaneEff;
555 for (int i=2;i--;) for (int j=2;j--;) fBeamInt[i][j] = rec.fBeamInt[i][j];
565 fIsNewRunLoader = rec.fIsNewRunLoader;
577 //_____________________________________________________________________________
578 AliReconstruction::~AliReconstruction()
583 if (fListOfCosmicTriggers) {
584 fListOfCosmicTriggers->Delete();
585 delete fListOfCosmicTriggers;
589 delete fCTPTimeParams;
590 delete fCTPTimeAlign;
592 if (fAlignObjArray) {
593 fAlignObjArray->Delete();
594 delete fAlignObjArray;
596 fSpecCDBUri.Delete();
598 AliCodeTimer::Instance()->Print();
601 //_____________________________________________________________________________
602 void AliReconstruction::InitQA()
604 //Initialize the QA and start of cycle
605 AliCodeTimerAuto("",0);
607 if (fInitQACalled) return;
608 fInitQACalled = kTRUE;
610 AliQAManager * qam = AliQAManager::QAManager(AliQAv1::kRECMODE) ;
611 if (fWriteQAExpertData)
612 qam->SetWriteExpert() ;
614 if (qam->IsDefaultStorageSet()) {
615 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
616 AliWarning("Default QA reference storage has been already set !");
617 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fQARefUri.Data()));
618 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
619 fQARefUri = qam->GetDefaultStorage()->GetURI();
621 if (fQARefUri.Length() > 0) {
622 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
623 AliDebug(2, Form("Default QA reference storage is set to: %s", fQARefUri.Data()));
624 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
626 fQARefUri="local://$ALICE_ROOT/QAref";
627 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
628 AliWarning("Default QA refeference storage not yet set !!!!");
629 AliWarning(Form("Setting it now to: %s", fQARefUri.Data()));
630 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
633 qam->SetDefaultStorage(fQARefUri);
637 qam->SetActiveDetectors(fQADetectors) ;
638 for (Int_t det = 0 ; det < AliQAv1::kNDET ; det++) {
639 qam->SetCycleLength(AliQAv1::DETECTORINDEX_t(det), fQACycles[det]) ;
640 qam->SetWriteExpert(AliQAv1::DETECTORINDEX_t(det)) ;
642 if (!fRawReader && !fInput && IsInTasks(AliQAv1::kRAWS))
643 fQATasks.ReplaceAll(Form("%d",AliQAv1::kRAWS), "") ;
644 qam->SetTasks(fQATasks) ;
645 qam->InitQADataMaker(AliCDBManager::Instance()->GetRun()) ;
648 Bool_t sameCycle = kFALSE ;
649 AliQADataMaker *qadm = qam->GetQADataMaker(AliQAv1::kGLOBAL);
650 AliInfo(Form("Initializing the global QA data maker"));
651 if (IsInTasks(AliQAv1::kRECPOINTS)) {
652 qadm->StartOfCycle(AliQAv1::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
653 TObjArray **arr=qadm->Init(AliQAv1::kRECPOINTS);
654 AliTracker::SetResidualsArray(arr);
657 if (IsInTasks(AliQAv1::kESDS)) {
658 qadm->StartOfCycle(AliQAv1::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
659 qadm->Init(AliQAv1::kESDS);
662 AliSysInfo::AddStamp("InitQA") ;
665 //_____________________________________________________________________________
666 void AliReconstruction::MergeQA(const char *fileName)
668 //Initialize the QA and start of cycle
669 AliCodeTimerAuto("",0) ;
670 AliQAManager::QAManager()->Merge(AliCDBManager::Instance()->GetRun(),fileName) ;
671 AliSysInfo::AddStamp("MergeQA") ;
674 //_____________________________________________________________________________
675 void AliReconstruction::InitCDB()
677 // activate a default CDB storage
678 // First check if we have any CDB storage set, because it is used
679 // to retrieve the calibration and alignment constants
680 AliCodeTimerAuto("",0);
682 if (fInitCDBCalled) return;
683 fInitCDBCalled = kTRUE;
685 AliCDBManager* man = AliCDBManager::Instance();
686 if (man->IsDefaultStorageSet())
688 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
689 AliWarning("Default CDB storage has been already set !");
690 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
691 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
692 fCDBUri = man->GetDefaultStorage()->GetURI();
695 if (fCDBUri.Length() > 0)
697 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
698 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
699 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
700 man->SetDefaultStorage(fCDBUri);
702 else if (!man->GetRaw()){
703 fCDBUri="local://$ALICE_ROOT/OCDB";
704 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
705 AliWarning("Default CDB storage not yet set !!!!");
706 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
707 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
708 man->SetDefaultStorage(fCDBUri);
711 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
712 AliWarning("Default storage will be set after setting the Run Number!!!");
713 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
717 // Now activate the detector specific CDB storage locations
718 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
719 TObject* obj = fSpecCDBUri[i];
721 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
722 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
723 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
724 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
726 AliSysInfo::AddStamp("InitCDB");
729 //_____________________________________________________________________________
730 void AliReconstruction::SetDefaultStorage(const char* uri) {
731 // Store the desired default CDB storage location
732 // Activate it later within the Run() method
738 //_____________________________________________________________________________
739 void AliReconstruction::SetQARefDefaultStorage(const char* uri) {
740 // Store the desired default CDB storage location
741 // Activate it later within the Run() method
744 AliQAv1::SetQARefStorage(fQARefUri.Data()) ;
747 //_____________________________________________________________________________
748 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
749 // Store a detector-specific CDB storage location
750 // Activate it later within the Run() method
752 AliCDBPath aPath(calibType);
753 if(!aPath.IsValid()){
754 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
755 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
756 if(!strcmp(calibType, fgkDetectorName[iDet])) {
757 aPath.SetPath(Form("%s/*", calibType));
758 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
762 if(!aPath.IsValid()){
763 AliError(Form("Not a valid path or detector: %s", calibType));
768 // // check that calibType refers to a "valid" detector name
769 // Bool_t isDetector = kFALSE;
770 // for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
771 // TString detName = fgkDetectorName[iDet];
772 // if(aPath.GetLevel0() == detName) {
773 // isDetector = kTRUE;
779 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
783 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
784 if (obj) fSpecCDBUri.Remove(obj);
785 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
789 //_____________________________________________________________________________
790 Bool_t AliReconstruction::SetRunNumberFromData()
792 // The method is called in Run() in order
793 // to set a correct run number.
794 // In case of raw data reconstruction the
795 // run number is taken from the raw data header
797 if (fSetRunNumberFromDataCalled) return kTRUE;
798 fSetRunNumberFromDataCalled = kTRUE;
800 AliCDBManager* man = AliCDBManager::Instance();
803 if(fRawReader->NextEvent()) {
804 if(man->GetRun() > 0) {
805 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
807 man->SetRun(fRawReader->GetRunNumber());
808 fRawReader->RewindEvents();
811 if(man->GetRun() > 0) {
812 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
815 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
821 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
823 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
828 // read run number from gAlice
829 if(rl->GetHeader()) {
830 man->SetRun(rl->GetHeader()->GetRun());
835 AliError("Neither run-loader header nor RawReader objects are found !");
847 //_____________________________________________________________________________
848 void AliReconstruction::SetCDBLock() {
849 // Set CDB lock: from now on it is forbidden to reset the run number
850 // or the default storage or to activate any further storage!
852 AliCDBManager::Instance()->SetLock(1);
855 //_____________________________________________________________________________
856 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
858 // Read the alignment objects from CDB.
859 // Each detector is supposed to have the
860 // alignment objects in DET/Align/Data CDB path.
861 // All the detector objects are then collected,
862 // sorted by geometry level (starting from ALIC) and
863 // then applied to the TGeo geometry.
864 // Finally an overlaps check is performed.
866 // Load alignment data from CDB and fill fAlignObjArray
867 if(fLoadAlignFromCDB){
869 TString detStr = detectors;
870 TString loadAlObjsListOfDets = "";
872 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
873 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
874 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
876 if(AliGeomManager::GetNalignable(fgkDetectorName[iDet]) != 0)
878 loadAlObjsListOfDets += fgkDetectorName[iDet];
879 loadAlObjsListOfDets += " ";
881 } // end loop over detectors
883 if(AliGeomManager::GetNalignable("GRP") != 0)
884 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
885 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
886 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
888 // Check if the array with alignment objects was
889 // provided by the user. If yes, apply the objects
890 // to the present TGeo geometry
891 if (fAlignObjArray) {
892 if (gGeoManager && gGeoManager->IsClosed()) {
893 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
894 AliError("The misalignment of one or more volumes failed!"
895 "Compare the list of simulated detectors and the list of detector alignment data!");
900 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
906 if (fAlignObjArray) {
907 fAlignObjArray->Delete();
908 delete fAlignObjArray; fAlignObjArray=NULL;
914 //_____________________________________________________________________________
915 void AliReconstruction::SetGAliceFile(const char* fileName)
917 // set the name of the galice file
919 fGAliceFileName = fileName;
922 //_____________________________________________________________________________
923 void AliReconstruction::SetInput(const char* input)
925 // In case the input string starts with 'mem://', we run in an online mode
926 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
927 // file is assumed. One can give as an input:
928 // mem://: - events taken from DAQ monitoring libs online
930 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
931 if (input) fRawInput = input;
934 //_____________________________________________________________________________
935 void AliReconstruction::SetOutput(const char* output)
937 // Set the output ESD filename
938 // 'output' is a normalt ROOT url
939 // The method is used in case of raw-data reco with PROOF
940 if (output) fESDOutput = output;
943 //_____________________________________________________________________________
944 void AliReconstruction::SetOption(const char* detector, const char* option)
946 // set options for the reconstruction of a detector
948 TObject* obj = fOptions.FindObject(detector);
949 if (obj) fOptions.Remove(obj);
950 fOptions.Add(new TNamed(detector, option));
953 //_____________________________________________________________________________
954 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
956 // Set custom reconstruction parameters for a given detector
957 // Single set of parameters for all the events
959 // First check if the reco-params are global
960 if(!strcmp(detector, "GRP")) {
962 fRecoParam.AddDetRecoParam(kNDetectors,par);
966 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
967 if(!strcmp(detector, fgkDetectorName[iDet])) {
969 fRecoParam.AddDetRecoParam(iDet,par);
976 //_____________________________________________________________________________
977 Bool_t AliReconstruction::InitGRP() {
978 //------------------------------------
979 // Initialization of the GRP entry
980 //------------------------------------
981 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
985 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
988 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
990 fGRPData = new AliGRPObject();
991 fGRPData->ReadValuesFromMap(m);
995 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
996 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
1000 // FIX ME: The unloading of GRP entry is temporarily disabled
1001 // because ZDC and VZERO are using it in order to initialize
1002 // their reconstructor objects. In the future one has to think
1003 // of propagating AliRunInfo to the reconstructors.
1004 // AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
1008 AliError("No GRP entry found in OCDB!");
1012 TString lhcState = fGRPData->GetLHCState();
1013 if (lhcState==AliGRPObject::GetInvalidString()) {
1014 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
1015 lhcState = "UNKNOWN";
1018 TString beamType = fGRPData->GetBeamType();
1019 if (beamType==AliGRPObject::GetInvalidString()) {
1020 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
1021 beamType = "UNKNOWN";
1024 Float_t beamEnergy = fGRPData->GetBeamEnergy();
1025 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
1026 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
1030 TString runType = fGRPData->GetRunType();
1031 if (runType==AliGRPObject::GetInvalidString()) {
1032 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
1033 runType = "UNKNOWN";
1036 Int_t activeDetectors = fGRPData->GetDetectorMask();
1037 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
1038 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
1039 activeDetectors = 1074790399;
1042 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
1046 // Process the list of active detectors
1047 if (activeDetectors) {
1048 UInt_t detMask = activeDetectors;
1049 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
1050 fRunTracking = MatchDetectorList(fRunTracking,detMask);
1051 fFillESD = MatchDetectorList(fFillESD,detMask);
1052 fQADetectors = MatchDetectorList(fQADetectors,detMask);
1053 fLoadCDB.Form("%s %s %s %s",
1054 fRunLocalReconstruction.Data(),
1055 fRunTracking.Data(),
1057 fQADetectors.Data());
1058 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
1059 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1) &&
1060 !((detMask >> AliDAQ::DetectorID("ITSSDD")) & 0x1) &&
1061 !((detMask >> AliDAQ::DetectorID("ITSSSD")) & 0x1) ) {
1062 // switch off the vertexer
1063 AliInfo("SPD,SDD,SSD is not in the list of active detectors. Vertexer and Trackleter are switched off.");
1064 fRunVertexFinder = kFALSE;
1065 fRunMultFinder = kFALSE;
1067 if (!((detMask >> AliDAQ::DetectorID("TRG")) & 0x1)) {
1068 // switch off the reading of CTP raw-data payload
1069 if (fFillTriggerESD) {
1070 AliInfo("CTP is not in the list of active detectors. CTP data reading switched off.");
1071 fFillTriggerESD = kFALSE;
1076 AliInfo("===================================================================================");
1077 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
1078 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
1079 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
1080 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
1081 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
1082 AliInfo("===================================================================================");
1084 //*** Dealing with the magnetic field map
1085 if ( TGeoGlobalMagField::Instance()->IsLocked() ) {
1086 if (TGeoGlobalMagField::Instance()->GetField()->TestBit(AliMagF::kOverrideGRP)) {
1087 AliInfo("ExpertMode!!! GRP information will be ignored !");
1088 AliInfo("ExpertMode!!! Running with the externally locked B field !");
1091 AliInfo("Destroying existing B field instance!");
1092 delete TGeoGlobalMagField::Instance();
1095 if ( !TGeoGlobalMagField::Instance()->IsLocked() ) {
1096 // Construct the field map out of the information retrieved from GRP.
1099 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
1100 if (l3Current == AliGRPObject::GetInvalidFloat()) {
1101 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1105 Char_t l3Polarity = fGRPData->GetL3Polarity();
1106 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
1107 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1112 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
1113 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
1114 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1118 Char_t diPolarity = fGRPData->GetDipolePolarity();
1119 if (diPolarity == AliGRPObject::GetInvalidChar()) {
1120 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1124 // read special bits for the polarity convention and map type
1125 Int_t polConvention = fGRPData->IsPolarityConventionLHC() ? AliMagF::kConvLHC : AliMagF::kConvDCS2008;
1126 Bool_t uniformB = fGRPData->IsUniformBMap();
1129 AliMagF* fld = AliMagF::CreateFieldMap(TMath::Abs(l3Current) * (l3Polarity ? -1:1),
1130 TMath::Abs(diCurrent) * (diPolarity ? -1:1),
1131 polConvention,uniformB,beamEnergy, beamType.Data());
1133 TGeoGlobalMagField::Instance()->SetField( fld );
1134 TGeoGlobalMagField::Instance()->Lock();
1135 AliInfo("Running with the B field constructed out of GRP !");
1137 else AliFatal("Failed to create a B field map !");
1139 else AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1142 //*** Get the diamond profiles from OCDB
1143 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1145 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1147 AliError("No SPD diamond profile found in OCDB!");
1150 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1152 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1154 AliError("No diamond profile found in OCDB!");
1157 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1159 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1161 AliError("No TPC diamond profile found in OCDB!");
1164 entry = AliCDBManager::Instance()->Get("GRP/Calib/CosmicTriggers");
1166 fListOfCosmicTriggers = dynamic_cast<THashTable*>(entry->GetObject());
1168 AliCDBManager::Instance()->UnloadFromCache("GRP/Calib/CosmicTriggers");
1171 if (!fListOfCosmicTriggers) {
1172 AliWarning("Can not get list of cosmic triggers from OCDB! Cosmic event specie will be effectively disabled!");
1178 //_____________________________________________________________________________
1179 Bool_t AliReconstruction::LoadCDB()
1181 // Load CDB entries for all active detectors.
1182 // By default we load all the entries in <det>/Calib
1185 AliCodeTimerAuto("",0);
1187 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1189 AliCDBManager::Instance()->Get("GRP/Calib/LHCClockPhase");
1191 TString detStr = fLoadCDB;
1192 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1193 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1194 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1197 // Temporary fix - one has to define the correct policy in order
1198 // to load the trigger OCDB entries only for the detectors that
1199 // in the trigger or that are needed in order to put correct
1200 // information in ESD
1201 AliCDBManager::Instance()->GetAll("TRIGGER/*/*");
1205 //_____________________________________________________________________________
1206 Bool_t AliReconstruction::LoadTriggerScalersCDB()
1208 // Load CTP scalers from OCDB.
1209 // The scalers are checked for consistency.
1211 AliCodeTimerAuto("",0);
1213 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/Scalers");
1217 AliInfo("Found an AliTriggerRunScalers in GRP/CTP/Scalers, reading it");
1218 fRunScalers = dynamic_cast<AliTriggerRunScalers*> (entry->GetObject());
1220 if (fRunScalers && (fRunScalers->CorrectScalersOverflow() == 0)) AliInfo("32bit Trigger counters corrected for overflow");
1225 //_____________________________________________________________________________
1226 Bool_t AliReconstruction::LoadCTPTimeParamsCDB()
1228 // Load CTP timing information (alignment)
1231 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
1232 if (!entry) return kFALSE;
1234 AliInfo("Found an AliCTPTimeParams in GRP/CTP/CTPtiming, reading it");
1235 fCTPTimeParams = dynamic_cast<AliCTPTimeParams*> (entry->GetObject());
1238 AliCDBEntry* entry2 = AliCDBManager::Instance()->Get("GRP/CTP/TimeAlign");
1239 if (!entry2) return kFALSE;
1241 AliInfo("Found an AliCTPTimeParams in GRP/CTP/TimeAlign, reading it");
1242 fCTPTimeAlign = dynamic_cast<AliCTPTimeParams*> (entry2->GetObject());
1243 entry2->SetOwner(0);
1248 //_____________________________________________________________________________
1249 Bool_t AliReconstruction::ReadIntensityInfoCDB()
1251 // Load LHC DIP data
1252 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/LHCData");
1255 AliInfo("Found an AliLHCData in GRP/GRP/LHCData, reading it");
1256 AliLHCData* dipData = dynamic_cast<AliLHCData*> (entry->GetObject());
1257 for (int ib=2;ib--;) {
1259 if (dipData && (dipData->GetMeanIntensity(ib,intI,intNI)>=0)) {
1260 fBeamInt[ib][0] = intI;
1261 fBeamInt[ib][1] = intNI;
1270 //_____________________________________________________________________________
1271 Bool_t AliReconstruction::Run(const char* input)
1274 AliCodeTimerAuto("",0);
1277 if (GetAbort() != TSelector::kContinue) return kFALSE;
1279 TChain *chain = NULL;
1280 if (fRawReader && (chain = fRawReader->GetChain())) {
1281 Long64_t nEntries = (fLastEvent < 0) ? (TChain::kBigNumber) : (fLastEvent - fFirstEvent + 1);
1284 // Temporary fix for long raw-data runs (until socket timeout handling in PROOF is revised)
1285 gProof->Exec("gEnv->SetValue(\"Proof.SocketActivityTimeout\",-1)", kTRUE);
1288 gProof->Exec("TGrid::Connect(\"alien://\")",kTRUE);
1290 TMessage::EnableSchemaEvolutionForAll(kTRUE);
1291 gProof->Exec("TMessage::EnableSchemaEvolutionForAll(kTRUE)",kTRUE);
1293 gProof->AddInput(this);
1295 if (!ParseOutput()) return kFALSE;
1297 gProof->SetParameter("PROOF_MaxSlavesPerNode", 9999);
1299 chain->Process("AliReconstruction","",nEntries,fFirstEvent);
1302 chain->Process(this,"",nEntries,fFirstEvent);
1307 if (GetAbort() != TSelector::kContinue) return kFALSE;
1309 if (GetAbort() != TSelector::kContinue) return kFALSE;
1310 //******* The loop over events
1311 AliInfo("Starting looping over events");
1313 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1314 (fRawReader && fRawReader->NextEvent())) {
1315 if (!ProcessEvent(iEvent)) {
1316 Abort("ProcessEvent",TSelector::kAbortFile);
1322 if (GetAbort() != TSelector::kContinue) return kFALSE;
1324 if (GetAbort() != TSelector::kContinue) return kFALSE;
1330 //_____________________________________________________________________________
1331 void AliReconstruction::InitRawReader(const char* input)
1333 // Init raw-reader and
1334 // set the input in case of raw data
1336 AliCodeTimerAuto("",0);
1338 if (input) fRawInput = input;
1339 fRawReader = AliRawReader::Create(fRawInput.Data());
1341 if (fRawInput.IsNull()) {
1342 AliInfo("Reconstruction will run over digits");
1345 AliFatal("Can not create raw-data reader ! Exiting...");
1349 if (!fEquipIdMap.IsNull() && fRawReader)
1350 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1352 if (!fUseHLTData.IsNull()) {
1353 // create the RawReaderHLT which performs redirection of HLT input data for
1354 // the specified detectors
1355 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1357 fParentRawReader=fRawReader;
1358 fRawReader=pRawReader;
1360 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1363 AliSysInfo::AddStamp("CreateRawReader");
1366 //_____________________________________________________________________________
1367 void AliReconstruction::InitRun(const char* input)
1369 // Initialization of raw-reader,
1370 // run number, CDB etc.
1371 AliCodeTimerAuto("",0);
1372 AliSysInfo::AddStamp("Start");
1374 // Initialize raw-reader if any
1375 InitRawReader(input);
1377 // Initialize the CDB storage
1380 // Set run number in CDBManager (if it is not already set by the user)
1381 if (!SetRunNumberFromData()) {
1382 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1386 // Set CDB lock: from now on it is forbidden to reset the run number
1387 // or the default storage or to activate any further storage!
1392 //_____________________________________________________________________________
1393 void AliReconstruction::Begin(TTree *)
1395 // Initialize AlReconstruction before
1396 // going into the event loop
1397 // Should follow the TSelector convention
1398 // i.e. initialize only the object on the client side
1399 AliCodeTimerAuto("",0);
1401 AliReconstruction *reco = NULL;
1403 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1406 AliSysInfo::AddStamp("ReadInputInBegin");
1409 // Import ideal TGeo geometry and apply misalignment
1411 TString geom(gSystem->DirName(fGAliceFileName));
1412 geom += "/geometry.root";
1413 AliGeomManager::LoadGeometry(geom.Data());
1415 Abort("LoadGeometry", TSelector::kAbortProcess);
1418 AliSysInfo::AddStamp("LoadGeom");
1419 TString detsToCheck=fRunLocalReconstruction;
1420 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1421 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1424 AliSysInfo::AddStamp("CheckGeom");
1427 if (!MisalignGeometry(fLoadAlignData)) {
1428 Abort("MisalignGeometry", TSelector::kAbortProcess);
1431 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1432 AliSysInfo::AddStamp("MisalignGeom");
1435 Abort("InitGRP", TSelector::kAbortProcess);
1438 AliSysInfo::AddStamp("InitGRP");
1441 Abort("LoadCDB", TSelector::kAbortProcess);
1444 AliSysInfo::AddStamp("LoadCDB");
1446 if (!LoadTriggerScalersCDB()) {
1447 Abort("LoadTriggerScalersCDB", TSelector::kAbortProcess);
1450 AliSysInfo::AddStamp("LoadTriggerScalersCDB");
1452 if (!LoadCTPTimeParamsCDB()) {
1453 Abort("LoadCTPTimeParamsCDB", TSelector::kAbortProcess);
1456 AliSysInfo::AddStamp("LoadCTPTimeParamsCDB");
1458 if (!ReadIntensityInfoCDB()) {
1459 Abort("ReadIntensityInfoCDB", TSelector::kAbortProcess);
1462 AliSysInfo::AddStamp("ReadIntensityInfoCDB");
1464 // Read the reconstruction parameters from OCDB
1465 if (!InitRecoParams()) {
1466 AliWarning("Not all detectors have correct RecoParam objects initialized");
1468 AliSysInfo::AddStamp("InitRecoParams");
1470 if (fInput && gProof) {
1471 if (reco) *reco = *this;
1473 gGeoManager->SetName("Geometry");
1474 gProof->AddInputData(gGeoManager,kTRUE);
1476 gProof->AddInputData(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()),kTRUE);
1477 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1478 AliMagF *magFieldMap = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1479 magFieldMap->SetName("MagneticFieldMap");
1480 gProof->AddInputData(magFieldMap,kTRUE);
1485 //_____________________________________________________________________________
1486 void AliReconstruction::SlaveBegin(TTree*)
1488 // Initialization related to run-loader,
1489 // vertexer, trackers, recontructors
1490 // In proof mode it is executed on the slave
1491 AliCodeTimerAuto("",0);
1493 TProofOutputFile *outProofFile = NULL;
1495 if (AliDebugLevel() > 0) fInput->Print();
1496 if (AliDebugLevel() > 10) fInput->Dump();
1497 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1500 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1502 AliGeomManager::SetGeometry(tgeo);
1504 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1505 Int_t runNumber = -1;
1506 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1507 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1508 man->SetCacheFlag(kTRUE);
1509 man->SetLock(kTRUE);
1513 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1514 AliMagF *newMap = new AliMagF(*map);
1515 if (!newMap->LoadParameterization()) {
1516 Abort("AliMagF::LoadParameterization", TSelector::kAbortProcess);
1519 TGeoGlobalMagField::Instance()->SetField(newMap);
1520 TGeoGlobalMagField::Instance()->Lock();
1522 if (TNamed *outputFileName = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE"))
1523 fProofOutputFileName = outputFileName->GetTitle();
1524 if (TNamed *outputLocation = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_LOCATION"))
1525 fProofOutputLocation = outputLocation->GetTitle();
1526 if (fInput->FindObject("PROOF_OUTPUTFILE_DATASET"))
1527 fProofOutputDataset = kTRUE;
1528 if (TNamed *archiveList = (TNamed*)fInput->FindObject("PROOF_OUTPUTFILE_ARCHIVE"))
1529 fProofOutputArchive = archiveList->GetTitle();
1530 if (!fProofOutputFileName.IsNull() &&
1531 !fProofOutputLocation.IsNull() &&
1532 fProofOutputArchive.IsNull()) {
1533 if (!fProofOutputDataset) {
1534 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"M");
1535 outProofFile->SetOutputFileName(Form("%s%s",fProofOutputLocation.Data(),fProofOutputFileName.Data()));
1538 outProofFile = new TProofOutputFile(fProofOutputFileName.Data(),"DROV",fProofOutputLocation.Data());
1540 if (AliDebugLevel() > 0) outProofFile->Dump();
1541 fOutput->Add(outProofFile);
1543 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1546 // get the run loader
1547 if (!InitRunLoader()) {
1548 Abort("InitRunLoader", TSelector::kAbortProcess);
1551 AliSysInfo::AddStamp("LoadLoader");
1553 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1556 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1557 Abort("CreateTrackers", TSelector::kAbortProcess);
1560 AliSysInfo::AddStamp("CreateTrackers");
1562 // create the ESD output file and tree
1563 if (!outProofFile) {
1564 ffile = TFile::Open("AliESDs.root", "RECREATE");
1565 ffile->SetCompressionLevel(2);
1566 if (!ffile->IsOpen()) {
1567 Abort("OpenESDFile", TSelector::kAbortProcess);
1572 AliInfo(Form("Opening output PROOF file: %s/%s",
1573 outProofFile->GetDir(), outProofFile->GetFileName()));
1574 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1575 Abort(Form("Problems opening output PROOF file: %s/%s",
1576 outProofFile->GetDir(), outProofFile->GetFileName()),
1577 TSelector::kAbortProcess);
1582 ftree = new TTree("esdTree", "Tree with ESD objects");
1583 fesd = new AliESDEvent();
1584 fesd->CreateStdContent();
1585 // add a so far non-std object to the ESD, this will
1586 // become part of the std content
1587 fesd->AddObject(new AliESDHLTDecision);
1589 fesd->WriteToTree(ftree);
1590 if (fWriteESDfriend) {
1591 ffileF = TFile::Open("AliESDfriends.root", "RECREATE");
1592 ftreeF = new TTree("esdFriendTree", "Tree with ESD Friend objects");
1593 fesdf = new AliESDfriend();
1594 ftreeF->Branch("ESDfriend.","AliESDfriend", &fesdf);
1595 fesd->AddObject(fesdf);
1598 ftree->GetUserInfo()->Add(fesd);
1600 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1601 fhltesd = new AliESDEvent();
1602 fhltesd->CreateStdContent();
1603 // read the ESD template from CDB
1604 // HLT is allowed to put non-std content to its ESD, the non-std
1605 // objects need to be created before invocation of WriteToTree in
1606 // order to create all branches. Initialization is done from an
1607 // ESD layout template in CDB
1608 AliCDBManager* man = AliCDBManager::Instance();
1609 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1610 AliCDBEntry* hltESDConfig=NULL;
1611 if (man->GetId(hltESDConfigPath)!=NULL &&
1612 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1613 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1615 // init all internal variables from the list of objects
1616 pESDLayout->GetStdContent();
1618 // copy content and create non-std objects
1619 *fhltesd=*pESDLayout;
1622 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1623 hltESDConfigPath.GetPath().Data()));
1627 fhltesd->WriteToTree(fhlttree);
1628 fhlttree->GetUserInfo()->Add(fhltesd);
1630 ProcInfo_t procInfo;
1631 gSystem->GetProcInfo(&procInfo);
1632 AliInfo(Form("Current memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
1635 //Initialize the QA and start of cycle
1636 if (fRunQA || fRunGlobalQA)
1639 //Initialize the Plane Efficiency framework
1640 if (fRunPlaneEff && !InitPlaneEff()) {
1641 Abort("InitPlaneEff", TSelector::kAbortProcess);
1645 if (strcmp(gProgName,"alieve") == 0)
1646 fRunAliEVE = InitAliEVE();
1651 //_____________________________________________________________________________
1652 Bool_t AliReconstruction::Process(Long64_t entry)
1654 // run the reconstruction over a single entry
1655 // from the chain with raw data
1656 AliCodeTimerAuto("",0);
1658 TTree *currTree = fChain->GetTree();
1659 AliRawVEvent *event = NULL;
1660 currTree->SetBranchAddress("rawevent",&event);
1661 currTree->GetEntry(entry);
1662 fRawReader = new AliRawReaderRoot(event);
1663 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1671 //_____________________________________________________________________________
1672 void AliReconstruction::Init(TTree *tree)
1674 // Implementation of TSelector::Init()
1677 AliError("The input tree is not found!");
1683 //_____________________________________________________________________________
1684 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1686 // run the reconstruction over a single event
1687 // The event loop is steered in Run method
1690 static Long_t oldMres=0;
1691 static Long_t oldMvir=0;
1692 static Float_t oldCPU=0;
1693 static Long_t aveDMres=0;
1694 static Long_t aveDMvir=0;
1695 static Float_t aveDCPU=0;
1697 AliCodeTimerAuto("",0);
1701 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1702 fRunLoader->SetEventNumber(iEvent);
1704 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1706 fRunLoader->TreeE()->Fill();
1707 if (fRawReader && fRawReader->UseAutoSaveESD())
1708 fRunLoader->TreeE()->AutoSave("SaveSelf");
1711 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1716 fRunLoader->GetEvent(iEvent);
1718 // Fill Event-info object
1720 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo,fListOfCosmicTriggers);
1722 ProcInfo_t procInfo;
1723 if(iEvent==fFirstEvent) {
1724 gSystem->GetProcInfo(&procInfo);
1725 oldMres=procInfo.fMemResident;
1726 oldMvir=procInfo.fMemVirtual;
1727 oldCPU=procInfo.fCpuUser+procInfo.fCpuSys;
1729 AliInfo(Form("================================= Processing event %d of type %-10s ==================================", iEvent,fRecoParam.PrintEventSpecie()));
1731 // Set the reco-params
1733 TString detStr = fLoadCDB;
1734 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
1735 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1736 AliReconstructor *reconstructor = GetReconstructor(iDet);
1737 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1738 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1739 reconstructor->SetRecoParam(par);
1740 reconstructor->GetPidSettings(&pid);
1741 reconstructor->SetEventInfo(&fEventInfo);
1743 AliQAManager::QAManager()->SetRecoParam(iDet, par) ;
1744 if (par) AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(par->GetEventSpecie())) ;
1749 const AliDetectorRecoParam *grppar = fRecoParam.GetDetRecoParam(kNDetectors);
1750 AliQAManager::QAManager()->SetRecoParam(AliQAv1::kGLOBAL, grppar) ;
1751 AliQAManager::QAManager()->SetEventSpecie(AliRecoParam::Convert(grppar->GetEventSpecie())) ;
1756 if (fRunQA && IsInTasks(AliQAv1::kRAWS)) {
1757 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
1758 AliQAManager::QAManager()->RunOneEvent(fRawReader) ;
1760 // local single event reconstruction
1761 if (!fRunLocalReconstruction.IsNull()) {
1762 TString detectors=fRunLocalReconstruction;
1763 // run HLT event reconstruction first
1764 // ;-( IsSelected changes the string
1765 if (IsSelected("HLT", detectors) &&
1766 !RunLocalEventReconstruction("HLT")) {
1767 if (fStopOnError) {CleanUp(); return kFALSE;}
1769 detectors=fRunLocalReconstruction;
1770 detectors.ReplaceAll("HLT", "");
1771 if (!RunLocalEventReconstruction(detectors)) {
1780 // fill Event header information from the RawEventHeader
1781 if (fRawReader){FillRawEventHeaderESD(fesd);}
1782 if (fRawReader){FillRawEventHeaderESD(fhltesd);}
1784 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1785 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1787 ((AliESDRun*)fesd->GetESDRun())->SetDetectorsInDAQ(fRunInfo->GetDetectorMask());
1788 ((AliESDRun*)fhltesd->GetESDRun())->SetDetectorsInDAQ(fRunInfo->GetDetectorMask());
1789 ((AliESDRun*)fesd->GetESDRun())->SetDetectorsInReco(AliDAQ::DetectorPatternOffline(fFillESD.Data()));
1790 ((AliESDRun*)fhltesd->GetESDRun())->SetDetectorsInReco(AliDAQ::DetectorPatternOffline(fFillESD.Data()));
1792 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1793 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1795 fesd->SetEventSpecie(fRecoParam.GetEventSpecie());
1796 fhltesd->SetEventSpecie(fRecoParam.GetEventSpecie());
1798 // Set magnetic field from the tracker
1799 fesd->SetMagneticField(AliTracker::GetBz());
1800 fhltesd->SetMagneticField(AliTracker::GetBz());
1802 AliESDRun *esdRun,*esdRunH;
1803 esdRun = (AliESDRun*)fesd->GetESDRun();
1804 esdRunH = (AliESDRun*)fhltesd->GetESDRun();
1805 esdRun->SetBeamEnergyIsSqrtSHalfGeV();
1806 esdRunH->SetBeamEnergyIsSqrtSHalfGeV();
1808 for (int ib=2;ib--;) for (int it=2;it--;) {
1809 esdRun->SetMeanIntensity(ib,it, fBeamInt[ib][it]);
1810 esdRunH->SetMeanIntensity(ib,it, fBeamInt[ib][it]);
1813 AliMagF* fld = (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
1814 if (fld) { // set info needed for field initialization
1815 fesd->SetCurrentL3(fld->GetCurrentSol());
1816 fesd->SetCurrentDip(fld->GetCurrentDip());
1817 fesd->SetBeamEnergy(fld->GetBeamEnergy());
1818 fesd->SetBeamType(fld->GetBeamTypeText());
1819 fesd->SetUniformBMap(fld->IsUniform());
1820 fesd->SetBInfoStored();
1822 fhltesd->SetCurrentL3(fld->GetCurrentSol());
1823 fhltesd->SetCurrentDip(fld->GetCurrentDip());
1824 fhltesd->SetBeamEnergy(fld->GetBeamEnergy());
1825 fhltesd->SetBeamType(fld->GetBeamTypeText());
1826 fhltesd->SetUniformBMap(fld->IsUniform());
1827 fhltesd->SetBInfoStored();
1830 // Set most probable pt, for B=0 tracking
1831 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1832 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(kNDetectors));
1833 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1835 // Fill raw-data error log into the ESD
1836 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1839 if (fRunVertexFinder) {
1840 if (!RunVertexFinder(fesd)) {
1841 if (fStopOnError) {CleanUp(); return kFALSE;}
1845 // For Plane Efficiency: run the SPD trackleter
1846 if (fRunPlaneEff && fSPDTrackleter) {
1847 if (!RunSPDTrackleting(fesd)) {
1848 if (fStopOnError) {CleanUp(); return kFALSE;}
1853 if (!fRunTracking.IsNull()) {
1854 if (fRunMuonTracking) {
1855 if (!RunMuonTracking(fesd)) {
1856 if (fStopOnError) {CleanUp(); return kFALSE;}
1862 if (!fRunTracking.IsNull()) {
1863 if (!RunTracking(fesd,pid)) {
1864 if (fStopOnError) {CleanUp(); return kFALSE;}
1869 if (!fFillESD.IsNull()) {
1870 TString detectors=fFillESD;
1871 // run HLT first and on hltesd
1872 // ;-( IsSelected changes the string
1873 if (IsSelected("HLT", detectors) &&
1874 !FillESD(fhltesd, "HLT")) {
1875 if (fStopOnError) {CleanUp(); return kFALSE;}
1878 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1879 if (detectors.Contains("ALL")) {
1881 for (Int_t idet=0; idet<kNDetectors; ++idet){
1882 detectors += fgkDetectorName[idet];
1886 detectors.ReplaceAll("HLT", "");
1887 if (!FillESD(fesd, detectors)) {
1888 if (fStopOnError) {CleanUp(); return kFALSE;}
1893 if (fReconstructor[3])
1894 GetReconstructor(3)->FillEventTimeWithTOF(fesd,&pid);
1899 if (fFillTriggerESD) {
1900 if (!FillTriggerESD(fesd)) {
1901 if (fStopOnError) {CleanUp(); return kFALSE;}
1904 // Always fill scalers
1905 if (!FillTriggerScalers(fesd)) {
1906 if (fStopOnError) {CleanUp(); return kFALSE;}
1913 // Propagate track to the beam pipe (if not already done by ITS)
1915 const Int_t ntracks = fesd->GetNumberOfTracks();
1916 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1919 UShort_t *selectedIdx=new UShort_t[ntracks];
1921 for (Int_t itrack=0; itrack<ntracks; itrack++){
1922 const Double_t kMaxStep = 1; //max step over the material
1925 AliESDtrack *track = fesd->GetTrack(itrack);
1926 if (!track) continue;
1928 AliExternalTrackParam *tpcTrack =
1929 (AliExternalTrackParam *)track->GetTPCInnerParam();
1933 PropagateTrackToBxByBz(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1936 Int_t n=trkArray.GetEntriesFast();
1937 selectedIdx[n]=track->GetID();
1938 trkArray.AddLast(tpcTrack);
1941 //Tracks refitted by ITS should already be at the SPD vertex
1942 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1945 PropagateTrackToBxByBz(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1946 Double_t x[3]; track->GetXYZ(x);
1947 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1948 track->RelateToVertexBxByBz(fesd->GetPrimaryVertexSPD(), b, kVeryBig);
1953 // Improve the reconstructed primary vertex position using the tracks
1955 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1956 if(fesd->GetPrimaryVertexSPD()) {
1957 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1958 if(vtitle.Contains("cosmics")) {
1959 runVertexFinderTracks=kFALSE;
1963 if (runVertexFinderTracks) {
1964 // TPC + ITS primary vertex
1965 ftVertexer->SetITSMode();
1966 ftVertexer->SetConstraintOff();
1967 // get cuts for vertexer from AliGRPRecoParam
1968 Bool_t constrSPD=kFALSE;
1970 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1971 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1972 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1973 ftVertexer->SetCuts(cutsVertexer);
1974 delete [] cutsVertexer; cutsVertexer = NULL;
1975 if(grpRecoParam->GetVertexerTracksConstraintITS()) {
1976 if(fDiamondProfile && fDiamondProfile->GetXRes()<kRadius){
1977 ftVertexer->SetVtxStart(fDiamondProfile); // apply constraint only if sigmax is smaller than the beam pipe radius
1979 if(fDiamondProfileSPD && fDiamondProfileSPD->GetXRes()<kRadius){
1980 ftVertexer->SetVtxStart(fDiamondProfileSPD);
1986 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1989 TString title=pvtx->GetTitle();
1990 title.Append("SPD");
1991 pvtx->SetTitle(title);
1993 if (pvtx->GetStatus()) {
1994 fesd->SetPrimaryVertexTracks(pvtx);
1995 for (Int_t i=0; i<ntracks; i++) {
1996 AliESDtrack *t = fesd->GetTrack(i);
1997 Double_t x[3]; t->GetXYZ(x);
1998 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1999 t->RelateToVertexBxByBz(pvtx, b, kVeryBig);
2002 delete pvtx; pvtx=NULL;
2005 // TPC-only primary vertex
2006 ftVertexer->SetTPCMode();
2007 ftVertexer->SetConstraintOff();
2008 // get cuts for vertexer from AliGRPRecoParam
2010 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
2011 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
2012 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
2013 ftVertexer->SetCuts(cutsVertexer);
2014 delete [] cutsVertexer; cutsVertexer = NULL;
2015 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC()) {
2016 if(fDiamondProfileTPC->GetXRes()<kRadius) ftVertexer->SetVtxStart(fDiamondProfileTPC); // apply constraint only if sigmax is smaller than the beam pipe radius
2019 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
2021 if (pvtx->GetStatus()) {
2022 fesd->SetPrimaryVertexTPC(pvtx);
2023 for (Int_t i=0; i<ntracks; i++) {
2024 AliESDtrack *t = fesd->GetTrack(i);
2025 Double_t x[3]; t->GetXYZ(x);
2026 Double_t b[3]; AliTracker::GetBxByBz(x,b);
2027 t->RelateToVertexTPCBxByBz(pvtx, b, kVeryBig);
2030 delete pvtx; pvtx=NULL;
2034 delete[] selectedIdx;
2036 if(fDiamondProfile && fDiamondProfile->GetXRes()<kRadius) fesd->SetDiamond(fDiamondProfile);
2037 else fesd->SetDiamond(fDiamondProfileSPD);
2041 AliV0vertexer vtxer;
2042 // get cuts for V0vertexer from AliGRPRecoParam
2044 Int_t nCutsV0vertexer = grpRecoParam->GetVertexerV0NCuts();
2045 Double_t *cutsV0vertexer = new Double_t[nCutsV0vertexer];
2046 grpRecoParam->GetVertexerV0Cuts(cutsV0vertexer);
2047 vtxer.SetCuts(cutsV0vertexer);
2048 delete [] cutsV0vertexer; cutsV0vertexer = NULL;
2050 vtxer.Tracks2V0vertices(fesd);
2052 if (fRunCascadeFinder) {
2054 AliCascadeVertexer cvtxer;
2055 // get cuts for CascadeVertexer from AliGRPRecoParam
2057 Int_t nCutsCascadeVertexer = grpRecoParam->GetVertexerCascadeNCuts();
2058 Double_t *cutsCascadeVertexer = new Double_t[nCutsCascadeVertexer];
2059 grpRecoParam->GetVertexerCascadeCuts(cutsCascadeVertexer);
2060 cvtxer.SetCuts(cutsCascadeVertexer);
2061 delete [] cutsCascadeVertexer; cutsCascadeVertexer = NULL;
2063 cvtxer.V0sTracks2CascadeVertices(fesd);
2068 if (fCleanESD) CleanESD(fesd);
2070 // RS run updated trackleter: since we want to mark the clusters used by tracks and also mark the
2071 // tracks interpreted as primary, this step should be done in the very end, when full
2072 // ESD info is available (particulalry, V0s)
2074 if (fRunMultFinder) {
2075 if (!RunMultFinder(fesd)) {
2076 if (fStopOnError) {CleanUp(); return kFALSE;}
2080 if (fRunQA && IsInTasks(AliQAv1::kESDS)) {
2081 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2082 AliQAManager::QAManager()->RunOneEvent(fesd, fhltesd) ;
2085 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2087 qadm->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2088 if (qadm && IsInTasks(AliQAv1::kESDS))
2089 qadm->Exec(AliQAv1::kESDS, fesd);
2092 // copy HLT decision from HLTesd to esd
2093 // the most relevant information is stored in a reduced container in the esd,
2094 // while the full information can be found in the HLTesd
2095 TObject* pHLTSrc=fhltesd->FindListObject(AliESDHLTDecision::Name());
2096 TObject* pHLTTgt=fesd->FindListObject(AliESDHLTDecision::Name());
2097 if (pHLTSrc && pHLTTgt) {
2098 pHLTSrc->Copy(*pHLTTgt);
2101 if (fWriteESDfriend)
2102 fesd->GetESDfriend(fesdf);
2105 if (fWriteESDfriend) {
2109 // Auto-save the ESD tree in case of prompt reco @P2
2110 if (fRawReader && fRawReader->UseAutoSaveESD()) {
2111 ftree->AutoSave("SaveSelf");
2112 if (fWriteESDfriend) ftreeF->AutoSave("SaveSelf");
2119 if (fRunAliEVE) RunAliEVE();
2123 if (fWriteESDfriend) {
2124 fesdf->~AliESDfriend();
2125 new (fesdf) AliESDfriend(); // Reset...
2128 gSystem->GetProcInfo(&procInfo);
2129 Long_t dMres=(procInfo.fMemResident-oldMres)/1024;
2130 Long_t dMvir=(procInfo.fMemVirtual-oldMvir)/1024;
2131 Float_t dCPU=procInfo.fCpuUser+procInfo.fCpuSys-oldCPU;
2132 aveDMres+=(dMres-aveDMres)/(iEvent-fFirstEvent+1);
2133 aveDMvir+=(dMvir-aveDMvir)/(iEvent-fFirstEvent+1);
2134 aveDCPU+=(dCPU-aveDCPU)/(iEvent-fFirstEvent+1);
2135 AliInfo(Form("======================= End Event %d: Res %ld(%3ld <%3ld>) Vir %ld(%3ld <%3ld>) CPU %5.2f <%5.2f> ===================",
2136 iEvent, procInfo.fMemResident/1024, dMres, aveDMres, procInfo.fMemVirtual/1024, dMvir, aveDMvir, dCPU, aveDCPU));
2137 oldMres=procInfo.fMemResident;
2138 oldMvir=procInfo.fMemVirtual;
2139 oldCPU=procInfo.fCpuUser+procInfo.fCpuSys;
2142 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2143 if (fReconstructor[iDet]) {
2144 fReconstructor[iDet]->SetRecoParam(NULL);
2145 fReconstructor[iDet]->SetEventInfo(NULL);
2147 if (fTracker[iDet]) fTracker[iDet]->SetEventInfo(NULL);
2150 if (fRunQA || fRunGlobalQA)
2151 AliQAManager::QAManager()->Increment() ;
2156 //_____________________________________________________________________________
2157 void AliReconstruction::SlaveTerminate()
2159 // Finalize the run on the slave side
2160 // Called after the exit
2161 // from the event loop
2162 AliCodeTimerAuto("",0);
2164 if (fIsNewRunLoader) { // galice.root didn't exist
2165 fRunLoader->WriteHeader("OVERWRITE");
2166 fRunLoader->CdGAFile();
2167 fRunLoader->Write(0, TObject::kOverwrite);
2170 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
2171 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
2173 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
2174 cdbMapCopy->SetOwner(1);
2175 cdbMapCopy->SetName("cdbMap");
2176 TIter iter(cdbMap->GetTable());
2179 while((pair = dynamic_cast<TPair*> (iter.Next()))){
2180 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
2181 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
2182 if (keyStr && valStr)
2183 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
2186 TList *cdbListCopy = new TList();
2187 cdbListCopy->SetOwner(1);
2188 cdbListCopy->SetName("cdbList");
2190 TIter iter2(cdbList);
2193 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
2194 cdbListCopy->Add(new TObjString(id->ToString().Data()));
2197 ftree->GetUserInfo()->Add(cdbMapCopy);
2198 ftree->GetUserInfo()->Add(cdbListCopy);
2203 // we want to have only one tree version number
2204 ftree->Write(ftree->GetName(),TObject::kOverwrite);
2205 fhlttree->Write(fhlttree->GetName(),TObject::kOverwrite);
2207 if (fWriteESDfriend) {
2209 ftreeF->Write(ftreeF->GetName(),TObject::kOverwrite);
2212 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
2213 if (fRunPlaneEff && !FinishPlaneEff()) {
2214 AliWarning("Finish PlaneEff evaluation failed");
2217 // End of cycle for the in-loop
2219 if (fRunQA || fRunGlobalQA) {
2220 AliQAManager::QAManager()->EndOfCycle() ;
2222 !fProofOutputLocation.IsNull() &&
2223 fProofOutputArchive.IsNull() &&
2224 !fProofOutputDataset) {
2225 TString qaOutputFile(Form("%sMerged.%s.Data.root",
2226 fProofOutputLocation.Data(),
2227 AliQAv1::GetQADataFileName()));
2228 TProofOutputFile *qaProofFile = new TProofOutputFile(Form("Merged.%s.Data.root",
2229 AliQAv1::GetQADataFileName()));
2230 qaProofFile->SetOutputFileName(qaOutputFile.Data());
2231 if (AliDebugLevel() > 0) qaProofFile->Dump();
2232 fOutput->Add(qaProofFile);
2233 MergeQA(qaProofFile->GetFileName());
2244 if (!fProofOutputFileName.IsNull() &&
2245 !fProofOutputLocation.IsNull() &&
2246 fProofOutputDataset &&
2247 !fProofOutputArchive.IsNull()) {
2248 TProofOutputFile *zipProofFile = new TProofOutputFile(fProofOutputFileName.Data(),
2250 fProofOutputLocation.Data());
2251 if (AliDebugLevel() > 0) zipProofFile->Dump();
2252 fOutput->Add(zipProofFile);
2253 TString fileList(fProofOutputArchive.Data());
2254 fileList.ReplaceAll(","," ");
2256 #if ROOT_SVN_REVISION >= 30174
2257 command.Form("zip -n root %s/%s %s",zipProofFile->GetDir(kTRUE),zipProofFile->GetFileName(),fileList.Data());
2259 command.Form("zip -n root %s/%s %s",zipProofFile->GetDir(),zipProofFile->GetFileName(),fileList.Data());
2261 AliInfo(Form("Executing: %s",command.Data()));
2262 gSystem->Exec(command.Data());
2267 //_____________________________________________________________________________
2268 void AliReconstruction::Terminate()
2270 // Create tags for the events in the ESD tree (the ESD tree is always present)
2271 // In case of empty events the tags will contain dummy values
2272 AliCodeTimerAuto("",0);
2274 // Do not call the ESD tag creator in case of PROOF-based reconstruction
2276 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
2277 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData, AliQAv1::Instance()->GetQA(), AliQAv1::Instance()->GetEventSpecies(), AliQAv1::kNDET, AliRecoParam::kNSpecies);
2278 delete esdtagCreator;
2281 // Cleanup of CDB manager: cache and active storages!
2282 AliCDBManager::Instance()->ClearCache();
2285 //_____________________________________________________________________________
2286 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
2288 // run the local reconstruction
2290 static Int_t eventNr=0;
2291 AliCodeTimerAuto("",0)
2293 TString detStr = detectors;
2294 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2295 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2296 AliReconstructor* reconstructor = GetReconstructor(iDet);
2297 if (!reconstructor) continue;
2298 AliLoader* loader = fLoader[iDet];
2299 // Matthias April 2008: temporary fix to run HLT reconstruction
2300 // although the HLT loader is missing
2301 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
2303 reconstructor->Reconstruct(fRawReader, NULL);
2306 reconstructor->Reconstruct(dummy, NULL);
2311 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
2314 // conversion of digits
2315 if (fRawReader && reconstructor->HasDigitConversion()) {
2316 AliInfo(Form("converting raw data digits into root objects for %s",
2317 fgkDetectorName[iDet]));
2318 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
2319 // fgkDetectorName[iDet]),0);
2320 loader->LoadDigits("update");
2321 loader->CleanDigits();
2322 loader->MakeDigitsContainer();
2323 TTree* digitsTree = loader->TreeD();
2324 reconstructor->ConvertDigits(fRawReader, digitsTree);
2325 loader->WriteDigits("OVERWRITE");
2326 loader->UnloadDigits();
2328 // local reconstruction
2329 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
2330 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]),0);
2331 loader->LoadRecPoints("update");
2332 loader->CleanRecPoints();
2333 loader->MakeRecPointsContainer();
2334 TTree* clustersTree = loader->TreeR();
2335 if (fRawReader && !reconstructor->HasDigitConversion()) {
2336 reconstructor->Reconstruct(fRawReader, clustersTree);
2338 loader->LoadDigits("read");
2339 TTree* digitsTree = loader->TreeD();
2341 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
2345 reconstructor->Reconstruct(digitsTree, clustersTree);
2346 if (fRunQA && IsInTasks(AliQAv1::kDIGITSR)) {
2347 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2348 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, digitsTree) ;
2351 loader->UnloadDigits();
2353 if (fRunQA && IsInTasks(AliQAv1::kRECPOINTS)) {
2354 AliQAManager::QAManager()->SetEventSpecie(fRecoParam.GetEventSpecie()) ;
2355 AliQAManager::QAManager()->RunOneEventInOneDetector(iDet, clustersTree) ;
2357 loader->WriteRecPoints("OVERWRITE");
2358 loader->UnloadRecPoints();
2359 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
2361 if (!IsSelected("CTP", detStr)) AliDebug(10,"No CTP");
2362 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2363 AliError(Form("the following detectors were not found: %s",
2371 //_____________________________________________________________________________
2372 Bool_t AliReconstruction::RunSPDTrackleting(AliESDEvent*& esd)
2374 // run the SPD trackleting (for SPD efficiency purpouses)
2376 AliCodeTimerAuto("",0)
2378 Double_t vtxPos[3] = {0, 0, 0};
2379 Double_t vtxErr[3] = {0.0, 0.0, 0.0};
2381 TArrayF mcVertex(3);
2383 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
2384 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
2385 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
2388 const AliESDVertex *vertex = esd->GetVertex();
2390 AliWarning("Vertex not found");
2393 vertex->GetXYZ(vtxPos);
2394 vertex->GetSigmaXYZ(vtxErr);
2395 if (fSPDTrackleter) {
2396 AliInfo("running the SPD Trackleter for Plane Efficiency Evaluation");
2399 fLoader[0]->LoadRecPoints("read");
2400 TTree* tree = fLoader[0]->TreeR();
2402 AliError("Can't get the ITS cluster tree");
2405 fSPDTrackleter->LoadClusters(tree);
2406 fSPDTrackleter->SetVertex(vtxPos, vtxErr);
2408 if (fSPDTrackleter->Clusters2Tracks(esd) != 0) {
2409 AliWarning("AliITSTrackleterSPDEff Clusters2Tracks failed");
2410 // fLoader[0]->UnloadRecPoints();
2413 //fSPDTrackleter->UnloadRecPoints();
2415 AliWarning("SPDTrackleter not available");
2421 //_____________________________________________________________________________
2422 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
2424 // run the barrel tracking
2426 AliCodeTimerAuto("",0)
2428 AliVertexer *vertexer = CreateVertexer();
2429 if (!vertexer) return kFALSE;
2431 AliInfo(Form("running the ITS vertex finder: %s",vertexer->ClassName()));
2432 AliESDVertex* vertex = NULL;
2434 fLoader[0]->LoadRecPoints();
2435 TTree* cltree = fLoader[0]->TreeR();
2437 if(fDiamondProfileSPD) vertexer->SetVtxStart(fDiamondProfileSPD);
2438 vertex = vertexer->FindVertexForCurrentEvent(cltree);
2441 AliError("Can't get the ITS cluster tree");
2443 fLoader[0]->UnloadRecPoints();
2446 AliError("Can't get the ITS loader");
2449 AliWarning("Vertex not found");
2450 vertex = new AliESDVertex();
2451 vertex->SetName("default");
2454 vertex->SetName("reconstructed");
2459 vertex->GetXYZ(vtxPos);
2460 vertex->GetSigmaXYZ(vtxErr);
2462 esd->SetPrimaryVertexSPD(vertex);
2463 AliESDVertex *vpileup = NULL;
2464 Int_t novertices = 0;
2465 vpileup = vertexer->GetAllVertices(novertices);
2467 for (Int_t kk=1; kk<novertices; kk++)esd->AddPileupVertexSPD(&vpileup[kk]);
2470 // if SPD multiplicity has been determined, it is stored in the ESD
2471 AliMultiplicity *mult = vertexer->GetMultiplicity();
2472 if(mult)esd->SetMultiplicity(mult);
2474 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2475 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2484 //_____________________________________________________________________________
2485 Bool_t AliReconstruction::RunMultFinder(AliESDEvent*& esd)
2487 // run the trackleter for multiplicity study
2489 AliCodeTimerAuto("",0)
2491 AliTrackleter *trackleter = CreateMultFinder();
2492 if (!trackleter) return kFALSE;
2494 AliInfo(Form("running the ITS multiplicity finder: %s",trackleter->ClassName()));
2497 fLoader[0]->LoadRecPoints();
2498 TTree* cltree = fLoader[0]->TreeR();
2500 trackleter->Reconstruct(esd,cltree);
2501 AliMultiplicity *mult = trackleter->GetMultiplicity();
2502 if(mult) esd->SetMultiplicity(mult);
2505 AliError("Can't get the ITS cluster tree");
2507 fLoader[0]->UnloadRecPoints();
2510 AliError("Can't get the ITS loader");
2518 //_____________________________________________________________________________
2519 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2521 // run the HLT barrel tracking
2523 AliCodeTimerAuto("",0)
2526 AliError("Missing runLoader!");
2530 AliInfo("running HLT tracking");
2532 // Get a pointer to the HLT reconstructor
2533 AliReconstructor *reconstructor = GetReconstructor(kNDetectors-1);
2534 if (!reconstructor) return kFALSE;
2537 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2538 TString detName = fgkDetectorName[iDet];
2539 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2540 reconstructor->SetOption(detName.Data());
2541 AliTracker *tracker = reconstructor->CreateTracker();
2543 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2544 if (fStopOnError) return kFALSE;
2548 Double_t vtxErr[3]={0.005,0.005,0.010};
2549 const AliESDVertex *vertex = esd->GetVertex();
2550 vertex->GetXYZ(vtxPos);
2551 tracker->SetVertex(vtxPos,vtxErr);
2553 fLoader[iDet]->LoadRecPoints("read");
2554 TTree* tree = fLoader[iDet]->TreeR();
2556 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2559 tracker->LoadClusters(tree);
2561 if (tracker->Clusters2Tracks(esd) != 0) {
2562 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2566 tracker->UnloadClusters();
2574 //_____________________________________________________________________________
2575 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2577 // run the muon spectrometer tracking
2579 AliCodeTimerAuto("",0)
2582 AliError("Missing runLoader!");
2585 Int_t iDet = 7; // for MUON
2587 AliInfo("is running...");
2589 // Get a pointer to the MUON reconstructor
2590 AliReconstructor *reconstructor = GetReconstructor(iDet);
2591 if (!reconstructor) return kFALSE;
2594 TString detName = fgkDetectorName[iDet];
2595 AliDebug(1, Form("%s tracking", detName.Data()));
2596 AliTracker *tracker = reconstructor->CreateTracker();
2598 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2603 fLoader[iDet]->LoadRecPoints("read");
2605 tracker->LoadClusters(fLoader[iDet]->TreeR());
2607 Int_t rv = tracker->Clusters2Tracks(esd);
2609 fLoader[iDet]->UnloadRecPoints();
2611 tracker->UnloadClusters();
2617 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2625 //_____________________________________________________________________________
2626 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd,AliESDpid &PID)
2628 // run the barrel tracking
2629 static Int_t eventNr=0;
2630 AliCodeTimerAuto("",0)
2632 AliInfo("running tracking");
2634 // Set the event info which is used
2635 // by the trackers in order to obtain
2636 // information about read-out detectors,
2638 AliDebug(1, "Setting event info");
2639 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2640 if (!fTracker[iDet]) continue;
2641 fTracker[iDet]->SetEventInfo(&fEventInfo);
2644 //Fill the ESD with the T0 info (will be used by the TOF)
2645 if (fReconstructor[11] && fLoader[11]) {
2646 fLoader[11]->LoadRecPoints("READ");
2647 TTree *treeR = fLoader[11]->TreeR();
2649 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2653 // pass 1: TPC + ITS inwards
2654 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2655 if (!fTracker[iDet]) continue;
2656 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2659 fLoader[iDet]->LoadRecPoints("read");
2660 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2661 TTree* tree = fLoader[iDet]->TreeR();
2663 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2666 fTracker[iDet]->LoadClusters(tree);
2667 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2669 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2670 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2673 // preliminary PID in TPC needed by the ITS tracker
2675 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2676 PID.MakePID(esd,kTRUE);
2678 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2681 // pass 2: ALL backwards
2683 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2684 if (!fTracker[iDet]) continue;
2685 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2688 if (iDet > 1) { // all except ITS, TPC
2690 fLoader[iDet]->LoadRecPoints("read");
2691 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2692 tree = fLoader[iDet]->TreeR();
2694 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2697 fTracker[iDet]->LoadClusters(tree);
2698 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2702 if (iDet>1) // start filling residuals for the "outer" detectors
2704 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2705 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2707 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2708 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2709 if ( elem && (! elem->At(0)) ) {
2710 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2711 if (qadm) qadm->InitRecPointsForTracker() ;
2715 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2716 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2721 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2722 fTracker[iDet]->UnloadClusters();
2723 fLoader[iDet]->UnloadRecPoints();
2725 // updated PID in TPC needed by the ITS tracker -MI
2727 //GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2728 //AliESDpid::MakePID(esd);
2729 PID.MakePID(esd,kTRUE);
2731 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2733 //stop filling residuals for the "outer" detectors
2734 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2736 // pass 3: TRD + TPC + ITS refit inwards
2738 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2739 if (!fTracker[iDet]) continue;
2740 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2743 if (iDet<2) // start filling residuals for TPC and ITS
2745 AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kTRUE);
2746 TObjArray ** arr = AliTracker::GetResidualsArray() ;
2748 AliRecoParam::EventSpecie_t es=fRecoParam.GetEventSpecie();
2749 TObjArray * elem = arr[AliRecoParam::AConvert(es)];
2750 if ( elem && (! elem->At(0)) ) {
2751 AliQADataMaker *qadm = AliQAManager::QAManager()->GetQADataMaker(AliQAv1::kGLOBAL);
2752 if (qadm) qadm->InitRecPointsForTracker() ;
2757 if (fTracker[iDet]->RefitInward(esd) != 0) {
2758 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2761 // run postprocessing
2762 if (fTracker[iDet]->PostProcess(esd) != 0) {
2763 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2766 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2769 // write space-points to the ESD in case alignment data output
2771 if (fWriteAlignmentData)
2772 WriteAlignmentData(esd);
2774 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2775 if (!fTracker[iDet]) continue;
2777 fTracker[iDet]->UnloadClusters();
2778 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2779 fLoader[iDet]->UnloadRecPoints();
2780 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2782 // stop filling residuals for TPC and ITS
2783 if (fRunGlobalQA) AliTracker::SetFillResiduals(fRecoParam.GetEventSpecie(), kFALSE);
2789 //_____________________________________________________________________________
2790 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2792 // Remove the data which are not needed for the physics analysis.
2795 Int_t nTracks=esd->GetNumberOfTracks();
2796 Int_t nV0s=esd->GetNumberOfV0s();
2798 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2800 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2801 Bool_t rc=esd->Clean(cleanPars);
2803 nTracks=esd->GetNumberOfTracks();
2804 nV0s=esd->GetNumberOfV0s();
2806 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2811 //_____________________________________________________________________________
2812 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2814 // fill the event summary data
2816 AliCodeTimerAuto("",0)
2817 static Int_t eventNr=0;
2818 TString detStr = detectors;
2820 AliSysInfo::AddStamp(Form("FillESDb%d",eventNr), -19,-19, eventNr);
2821 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2822 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2823 AliReconstructor* reconstructor = GetReconstructor(iDet);
2824 if (!reconstructor) continue;
2825 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2826 TTree* clustersTree = NULL;
2827 if (fLoader[iDet]) {
2828 fLoader[iDet]->LoadRecPoints("read");
2829 clustersTree = fLoader[iDet]->TreeR();
2830 if (!clustersTree) {
2831 AliError(Form("Can't get the %s clusters tree",
2832 fgkDetectorName[iDet]));
2833 if (fStopOnError) return kFALSE;
2836 if (fRawReader && !reconstructor->HasDigitConversion()) {
2837 reconstructor->FillESD(fRawReader, clustersTree, esd);
2839 TTree* digitsTree = NULL;
2840 if (fLoader[iDet]) {
2841 fLoader[iDet]->LoadDigits("read");
2842 digitsTree = fLoader[iDet]->TreeD();
2844 AliError(Form("Can't get the %s digits tree",
2845 fgkDetectorName[iDet]));
2846 if (fStopOnError) return kFALSE;
2849 reconstructor->FillESD(digitsTree, clustersTree, esd);
2850 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2852 if (fLoader[iDet]) {
2853 fLoader[iDet]->UnloadRecPoints();
2857 if (!IsSelected("CTP", detStr)) AliDebug(10,"No CTP");
2858 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2859 AliError(Form("the following detectors were not found: %s",
2861 if (fStopOnError) return kFALSE;
2863 AliSysInfo::AddStamp(Form("FillESDe%d",eventNr), -20,-20, eventNr);
2868 //_____________________________________________________________________________
2869 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2871 // Reads the trigger decision which is
2872 // stored in Trigger.root file and fills
2873 // the corresponding esd entries
2875 AliCodeTimerAuto("",0)
2877 AliInfo("Filling trigger information into the ESD");
2880 AliCTPRawStream input(fRawReader);
2881 if (!input.Next()) {
2882 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2885 if (esd->GetTriggerMask() != input.GetClassMask())
2886 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2887 input.GetClassMask(),esd->GetTriggerMask()));
2888 if (esd->GetOrbitNumber() != input.GetOrbitID())
2889 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2890 input.GetOrbitID(),esd->GetOrbitNumber()));
2891 if (esd->GetBunchCrossNumber() != input.GetBCID())
2892 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2893 input.GetBCID(),esd->GetBunchCrossNumber()));
2894 AliESDHeader* esdheader = esd->GetHeader();
2895 esdheader->SetL0TriggerInputs(input.GetL0Inputs());
2896 esdheader->SetL1TriggerInputs(input.GetL1Inputs());
2897 esdheader->SetL2TriggerInputs(input.GetL2Inputs());
2899 UInt_t orbit=input.GetOrbitID();
2900 for(Int_t i=0 ; i<input.GetNIRs() ; i++ )
2901 if(TMath::Abs(Int_t(orbit-(input.GetIR(i))->GetOrbit()))<=1){
2902 esdheader->AddTriggerIR(input.GetIR(i));
2908 //_____________________________________________________________________________
2909 Bool_t AliReconstruction::FillTriggerScalers(AliESDEvent*& esd)
2912 //fRunScalers->Print();
2913 if(fRunScalers && fRunScalers->CheckRunScalers()){
2914 AliTimeStamp* timestamp = new AliTimeStamp(esd->GetOrbitNumber(), esd->GetPeriodNumber(), esd->GetBunchCrossNumber());
2915 //AliTimeStamp* timestamp = new AliTimeStamp(10308000, 0, (ULong64_t)486238);
2916 AliESDHeader* esdheader = fesd->GetHeader();
2917 for(Int_t i=0;i<50;i++){
2918 if((1ull<<i) & esd->GetTriggerMask()){
2919 AliTriggerScalersESD* scalesd = fRunScalers->GetScalersForEventClass( timestamp, i+1);
2920 if(scalesd)esdheader->SetTriggerScalersRecord(scalesd);
2926 //_____________________________________________________________________________
2927 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2930 // Filling information from RawReader Header
2933 if (!fRawReader) return kFALSE;
2935 AliInfo("Filling information from RawReader Header");
2937 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2938 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2939 esd->SetPeriodNumber(fRawReader->GetPeriod());
2941 esd->SetTimeStamp(fRawReader->GetTimestamp());
2942 esd->SetEventType(fRawReader->GetType());
2948 //_____________________________________________________________________________
2949 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2951 // check whether detName is contained in detectors
2952 // if yes, it is removed from detectors
2954 // check if all detectors are selected
2955 if ((detectors.CompareTo("ALL") == 0) ||
2956 detectors.BeginsWith("ALL ") ||
2957 detectors.EndsWith(" ALL") ||
2958 detectors.Contains(" ALL ")) {
2963 // search for the given detector
2964 Bool_t result = kFALSE;
2965 if ((detectors.CompareTo(detName) == 0) ||
2966 detectors.BeginsWith(detName+" ") ||
2967 detectors.EndsWith(" "+detName) ||
2968 detectors.Contains(" "+detName+" ")) {
2969 detectors.ReplaceAll(detName, "");
2973 // clean up the detectors string
2974 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2975 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2976 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2981 //_____________________________________________________________________________
2982 Bool_t AliReconstruction::InitRunLoader()
2984 // get or create the run loader
2986 if (gAlice) delete gAlice;
2989 TFile *gafile = TFile::Open(fGAliceFileName.Data());
2990 // if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2991 if (gafile) { // galice.root exists
2995 // load all base libraries to get the loader classes
2996 TString libs = gSystem->GetLibraries();
2997 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
2998 TString detName = fgkDetectorName[iDet];
2999 if (detName == "HLT") continue;
3000 if (libs.Contains("lib" + detName + "base.so")) continue;
3001 gSystem->Load("lib" + detName + "base.so");
3003 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
3005 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
3010 fRunLoader->CdGAFile();
3011 fRunLoader->LoadgAlice();
3013 //PH This is a temporary fix to give access to the kinematics
3014 //PH that is needed for the labels of ITS clusters
3015 fRunLoader->LoadHeader();
3016 fRunLoader->LoadKinematics();
3018 } else { // galice.root does not exist
3020 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
3022 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
3023 AliConfig::GetDefaultEventFolderName(),
3026 AliError(Form("could not create run loader in file %s",
3027 fGAliceFileName.Data()));
3031 fIsNewRunLoader = kTRUE;
3032 fRunLoader->MakeTree("E");
3034 if (fNumberOfEventsPerFile > 0)
3035 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
3037 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
3043 //_____________________________________________________________________________
3044 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
3046 // get the reconstructor object and the loader for a detector
3048 if (fReconstructor[iDet]) {
3049 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
3050 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
3051 fReconstructor[iDet]->SetRecoParam(par);
3052 fReconstructor[iDet]->SetRunInfo(fRunInfo);
3054 return fReconstructor[iDet];
3057 // load the reconstructor object
3058 TPluginManager* pluginManager = gROOT->GetPluginManager();
3059 TString detName = fgkDetectorName[iDet];
3060 TString recName = "Ali" + detName + "Reconstructor";
3062 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
3064 AliReconstructor* reconstructor = NULL;
3065 // first check if a plugin is defined for the reconstructor
3066 TPluginHandler* pluginHandler =
3067 pluginManager->FindHandler("AliReconstructor", detName);
3068 // if not, add a plugin for it
3069 if (!pluginHandler) {
3070 AliDebug(1, Form("defining plugin for %s", recName.Data()));
3071 TString libs = gSystem->GetLibraries();
3072 if (libs.Contains("lib" + detName + "base.so") ||
3073 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
3074 pluginManager->AddHandler("AliReconstructor", detName,
3075 recName, detName + "rec", recName + "()");
3077 pluginManager->AddHandler("AliReconstructor", detName,
3078 recName, detName, recName + "()");
3080 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
3082 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
3083 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
3085 if (reconstructor) {
3086 TObject* obj = fOptions.FindObject(detName.Data());
3087 if (obj) reconstructor->SetOption(obj->GetTitle());
3088 reconstructor->SetRunInfo(fRunInfo);
3089 reconstructor->Init();
3090 fReconstructor[iDet] = reconstructor;
3093 // get or create the loader
3094 if (detName != "HLT") {
3095 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
3096 if (!fLoader[iDet]) {
3097 AliConfig::Instance()
3098 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
3100 // first check if a plugin is defined for the loader
3102 pluginManager->FindHandler("AliLoader", detName);
3103 // if not, add a plugin for it
3104 if (!pluginHandler) {
3105 TString loaderName = "Ali" + detName + "Loader";
3106 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
3107 pluginManager->AddHandler("AliLoader", detName,
3108 loaderName, detName + "base",
3109 loaderName + "(const char*, TFolder*)");
3110 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
3112 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
3114 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
3115 fRunLoader->GetEventFolder());
3117 if (!fLoader[iDet]) { // use default loader
3118 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
3120 if (!fLoader[iDet]) {
3121 AliWarning(Form("couldn't get loader for %s", detName.Data()));
3122 if (fStopOnError) return NULL;
3124 fRunLoader->AddLoader(fLoader[iDet]);
3125 fRunLoader->CdGAFile();
3126 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
3127 fRunLoader->Write(0, TObject::kOverwrite);
3132 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
3133 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
3134 if (reconstructor) {
3135 reconstructor->SetRecoParam(par);
3136 reconstructor->SetRunInfo(fRunInfo);
3139 return reconstructor;
3142 //_____________________________________________________________________________
3143 AliVertexer* AliReconstruction::CreateVertexer()
3145 // create the vertexer
3146 // Please note that the caller is the owner of the
3149 AliVertexer* vertexer = NULL;
3150 AliReconstructor* itsReconstructor = GetReconstructor(0);
3151 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) ||
3152 fRunTracking.Contains("ITS") || fFillESD.Contains("ITS") )) {
3153 vertexer = itsReconstructor->CreateVertexer();
3156 AliWarning("couldn't create a vertexer for ITS");
3162 //_____________________________________________________________________________
3163 AliTrackleter* AliReconstruction::CreateMultFinder()
3165 // create the ITS trackleter for mult. estimation
3166 // Please note that the caller is the owner of the
3169 AliTrackleter* trackleter = NULL;
3170 AliReconstructor* itsReconstructor = GetReconstructor(0);
3171 if (itsReconstructor && ((fRunLocalReconstruction.Contains("ITS")) ||
3172 fRunTracking.Contains("ITS") || fFillESD.Contains("ITS") )) {
3173 trackleter = itsReconstructor->CreateMultFinder();
3176 AliWarning("ITS is not in reconstruction, switching off RunMultFinder");
3177 fRunMultFinder = kFALSE;
3183 //_____________________________________________________________________________
3184 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
3186 // create the trackers
3187 AliInfo("Creating trackers");
3189 TString detStr = detectors;
3190 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3191 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3192 AliReconstructor* reconstructor = GetReconstructor(iDet);
3193 if (!reconstructor) continue;
3194 TString detName = fgkDetectorName[iDet];
3195 if (detName == "HLT") {
3196 fRunHLTTracking = kTRUE;
3199 if (detName == "MUON") {
3200 fRunMuonTracking = kTRUE;
3204 fTracker[iDet] = reconstructor->CreateTracker();
3205 if (!fTracker[iDet] && (iDet < 7)) {
3206 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
3207 if (fStopOnError) return kFALSE;
3209 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
3215 //_____________________________________________________________________________
3216 void AliReconstruction::CleanUp()
3218 // delete trackers and the run loader and close and delete the file
3220 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3221 delete fReconstructor[iDet];
3222 fReconstructor[iDet] = NULL;
3223 fLoader[iDet] = NULL;
3224 delete fTracker[iDet];
3225 fTracker[iDet] = NULL;
3230 delete fSPDTrackleter;
3231 fSPDTrackleter = NULL;
3240 delete fParentRawReader;
3241 fParentRawReader=NULL;
3249 if (AliQAManager::QAManager())
3250 AliQAManager::QAManager()->ShowQA() ;
3251 AliQAManager::Destroy() ;
3255 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
3257 // Write space-points which are then used in the alignment procedures
3258 // For the moment only ITS, TPC, TRD and TOF
3260 Int_t ntracks = esd->GetNumberOfTracks();
3261 for (Int_t itrack = 0; itrack < ntracks; itrack++)
3263 AliESDtrack *track = esd->GetTrack(itrack);
3266 for (Int_t i=0; i<200; ++i) idx[i] = -1; //PH avoid uninitialized values
3267 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
3268 nsp += (iDet==GetDetIndex("TRD")) ? track->GetTRDntracklets():track->GetNcls(iDet);
3270 if (iDet==GetDetIndex("ITS")) { // ITS "extra" clusters
3271 track->GetClusters(iDet,idx);
3272 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
3277 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
3278 track->SetTrackPointArray(sp);
3280 for (Int_t iDet = 5; iDet >= 0; iDet--) {
3281 AliTracker *tracker = fTracker[iDet];
3282 if (!tracker) continue;
3283 Int_t nspdet = (iDet==GetDetIndex("TRD")) ? track->GetTRDtracklets(idx):track->GetClusters(iDet,idx);
3285 if (iDet==GetDetIndex("ITS")) // ITS "extra" clusters
3286 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
3288 if (nspdet <= 0) continue;
3292 while (isp2 < nspdet) {
3293 Bool_t isvalid=kTRUE;
3295 Int_t index=idx[isp++];
3296 if (index < 0) continue;
3298 TString dets = fgkDetectorName[iDet];
3299 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
3300 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
3301 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
3302 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
3303 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
3305 isvalid = tracker->GetTrackPoint(index,p);
3308 if (!isvalid) continue;
3309 if (iDet==GetDetIndex("ITS") && (isp-1)>=6) p.SetExtra();
3310 sp->AddPoint(isptrack,&p); isptrack++;
3317 //_____________________________________________________________________________
3318 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
3320 // The method reads the raw-data error log
3321 // accumulated within the rawReader.
3322 // It extracts the raw-data errors related to
3323 // the current event and stores them into
3324 // a TClonesArray inside the esd object.
3326 if (!fRawReader) return;
3328 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
3330 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
3332 if (iEvent != log->GetEventNumber()) continue;
3334 esd->AddRawDataErrorLog(log);
3339 //_____________________________________________________________________________
3340 // void AliReconstruction::CheckQA()
3342 // check the QA of SIM for this run and remove the detectors
3343 // with status Fatal
3345 // TString newRunLocalReconstruction ;
3346 // TString newRunTracking ;
3347 // TString newFillESD ;
3349 // for (Int_t iDet = 0; iDet < AliQAv1::kNDET; iDet++) {
3350 // TString detName(AliQAv1::GetDetName(iDet)) ;
3351 // AliQAv1 * qa = AliQAv1::Instance(AliQAv1::DETECTORINDEX_t(iDet)) ;
3352 // if ( qa->IsSet(AliQAv1::DETECTORINDEX_t(iDet), AliQAv1::kSIM, specie, AliQAv1::kFATAL)) {
3353 // AliInfo(Form("QA status for %s %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed",
3354 // detName.Data(), AliRecoParam::GetEventSpecieName(es))) ;
3356 // if ( fRunLocalReconstruction.Contains(AliQAv1::GetDetName(iDet)) ||
3357 // fRunLocalReconstruction.Contains("ALL") ) {
3358 // newRunLocalReconstruction += detName ;
3359 // newRunLocalReconstruction += " " ;
3361 // if ( fRunTracking.Contains(AliQAv1::GetDetName(iDet)) ||
3362 // fRunTracking.Contains("ALL") ) {
3363 // newRunTracking += detName ;
3364 // newRunTracking += " " ;
3366 // if ( fFillESD.Contains(AliQAv1::GetDetName(iDet)) ||
3367 // fFillESD.Contains("ALL") ) {
3368 // newFillESD += detName ;
3369 // newFillESD += " " ;
3373 // fRunLocalReconstruction = newRunLocalReconstruction ;
3374 // fRunTracking = newRunTracking ;
3375 // fFillESD = newFillESD ;
3378 //_____________________________________________________________________________
3379 Int_t AliReconstruction::GetDetIndex(const char* detector)
3381 // return the detector index corresponding to detector
3383 for (index = 0; index < kNDetectors ; index++) {
3384 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
3389 //_____________________________________________________________________________
3390 Bool_t AliReconstruction::FinishPlaneEff() {
3392 // Here execute all the necessary operationis, at the end of the tracking phase,
3393 // in case that evaluation of PlaneEfficiencies was required for some detector.
3394 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
3396 // This Preliminary version works 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 TString detStr = fLoadCDB;
3404 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
3405 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
3406 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3407 if(fTracker[iDet] && fTracker[iDet]->GetPlaneEff()) {
3408 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
3409 TString name=planeeff->GetName();
3411 TFile* pefile = TFile::Open(name, "RECREATE");
3412 ret=(Bool_t)planeeff->Write();
3414 if(planeeff->GetCreateHistos()) {
3415 TString hname=planeeff->GetName();
3416 hname+="Histo.root";
3417 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
3420 if(fSPDTrackleter) {
3421 AliPlaneEff *planeeff=fSPDTrackleter->GetPlaneEff();
3422 TString name="AliITSPlaneEffSPDtracklet.root";
3423 TFile* pefile = TFile::Open(name, "RECREATE");
3424 ret=(Bool_t)planeeff->Write();
3426 AliESDEvent *dummy=NULL;
3427 ret=(Bool_t)fSPDTrackleter->PostProcess(dummy); // take care of writing other files
3432 //_____________________________________________________________________________
3433 Bool_t AliReconstruction::InitPlaneEff() {
3435 // Here execute all the necessary operations, before of the tracking phase,
3436 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
3437 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
3438 // which should be updated/recalculated.
3440 // This Preliminary version will work only FOR ITS !!!!!
3441 // other detectors (TOF,TRD, etc. have to develop their specific codes)
3444 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
3447 fSPDTrackleter = NULL;
3448 TString detStr = fLoadCDB;
3449 if (IsSelected(fgkDetectorName[0], detStr)) {
3450 AliReconstructor* itsReconstructor = GetReconstructor(0);
3451 if (itsReconstructor) {
3452 fSPDTrackleter = itsReconstructor->CreateTrackleter(); // this is NULL unless required in RecoParam
3454 if (fSPDTrackleter) {
3455 AliInfo("Trackleter for SPD has been created");
3461 //_____________________________________________________________________________
3462 Bool_t AliReconstruction::InitAliEVE()
3464 // This method should be called only in case
3465 // AliReconstruction is run
3466 // within the alieve environment.
3467 // It will initialize AliEVE in a way
3468 // so that it can visualize event processed
3469 // by AliReconstruction.
3470 // The return flag shows whenever the
3471 // AliEVE initialization was successful or not.
3473 TString macroStr(gSystem->Getenv("ALIEVE_ONLINE_MACRO"));
3475 if (macroStr.IsNull())
3476 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
3478 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
3480 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
3482 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
3483 gROOT->ProcessLine("alieve_online_init()");
3488 //_____________________________________________________________________________
3489 void AliReconstruction::RunAliEVE()
3491 // Runs AliEVE visualisation of
3492 // the current event.
3493 // Should be executed only after
3494 // successful initialization of AliEVE.
3496 AliInfo("Running AliEVE...");
3497 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p,(AliESDfriend*)%p);",fRunLoader,fRawReader,fesd,fesdf));
3501 //_____________________________________________________________________________
3502 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
3504 // Allows to run QA for a selected set of detectors
3505 // and a selected set of tasks among RAWS, DIGITSR, RECPOINTS and ESDS
3506 // all selected detectors run the same selected tasks
3508 if (!detAndAction.Contains(":")) {
3509 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
3513 Int_t colon = detAndAction.Index(":") ;
3514 fQADetectors = detAndAction(0, colon) ;
3515 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
3516 if (fQATasks.Contains("ALL") ) {
3517 fQATasks = Form("%d %d %d %d", AliQAv1::kRAWS, AliQAv1::kDIGITSR, AliQAv1::kRECPOINTS, AliQAv1::kESDS) ;
3519 fQATasks.ToUpper() ;
3521 if ( fQATasks.Contains("RAW") )
3522 tempo = Form("%d ", AliQAv1::kRAWS) ;
3523 if ( fQATasks.Contains("DIGIT") )
3524 tempo += Form("%d ", AliQAv1::kDIGITSR) ;
3525 if ( fQATasks.Contains("RECPOINT") )
3526 tempo += Form("%d ", AliQAv1::kRECPOINTS) ;
3527 if ( fQATasks.Contains("ESD") )
3528 tempo += Form("%d ", AliQAv1::kESDS) ;
3530 if (fQATasks.IsNull()) {
3531 AliInfo("No QA requested\n") ;
3536 TString tempo(fQATasks) ;
3537 tempo.ReplaceAll(Form("%d", AliQAv1::kRAWS), AliQAv1::GetTaskName(AliQAv1::kRAWS)) ;
3538 tempo.ReplaceAll(Form("%d", AliQAv1::kDIGITSR), AliQAv1::GetTaskName(AliQAv1::kDIGITSR)) ;
3539 tempo.ReplaceAll(Form("%d", AliQAv1::kRECPOINTS), AliQAv1::GetTaskName(AliQAv1::kRECPOINTS)) ;
3540 tempo.ReplaceAll(Form("%d", AliQAv1::kESDS), AliQAv1::GetTaskName(AliQAv1::kESDS)) ;
3541 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
3546 //_____________________________________________________________________________
3547 Bool_t AliReconstruction::InitRecoParams()
3549 // The method accesses OCDB and retrieves all
3550 // the available reco-param objects from there.
3552 Bool_t isOK = kTRUE;
3554 if (fRecoParam.GetDetRecoParamArray(kNDetectors)) {
3555 AliInfo("Using custom GRP reconstruction parameters");
3558 AliInfo("Loading GRP reconstruction parameter objects");
3560 AliCDBPath path("GRP","Calib","RecoParam");
3561 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3563 AliWarning("Couldn't find GRP RecoParam entry in OCDB");
3567 TObject *recoParamObj = entry->GetObject();
3568 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3569 // GRP has a normal TobjArray of AliDetectorRecoParam objects
3570 // Registering them in AliRecoParam
3571 fRecoParam.AddDetRecoParamArray(kNDetectors,dynamic_cast<TObjArray*>(recoParamObj));
3573 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3574 // GRP has only onse set of reco parameters
3575 // Registering it in AliRecoParam
3576 AliInfo("Single set of GRP reconstruction parameters found");
3577 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3578 fRecoParam.AddDetRecoParam(kNDetectors,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3581 AliError("No valid GRP RecoParam object found in the OCDB");
3588 TString detStr = fLoadCDB;
3589 for (Int_t iDet = 0; iDet < kNDetectors; iDet++) {
3591 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3593 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3594 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3598 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3600 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3601 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3603 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3607 TObject *recoParamObj = entry->GetObject();
3608 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3609 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3610 // Registering them in AliRecoParam
3611 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3613 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3614 // The detector has only onse set of reco parameters
3615 // Registering it in AliRecoParam
3616 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3617 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3618 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3621 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3625 // FIX ME: We have to disable the unloading of reco-param CDB
3626 // entries because QA framework is using them. Has to be fix in
3627 // a way that the QA takes the objects already constructed in
3629 // AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3633 if (AliDebugLevel() > 0) fRecoParam.Print();
3638 //_____________________________________________________________________________
3639 Bool_t AliReconstruction::GetEventInfo()
3641 // Fill the event info object
3643 AliCodeTimerAuto("",0)
3645 AliCentralTrigger *aCTP = NULL;
3647 fEventInfo.SetEventType(fRawReader->GetType());
3649 ULong64_t mask = fRawReader->GetClassMask();
3650 fEventInfo.SetTriggerMask(mask);
3651 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3652 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3654 aCTP = new AliCentralTrigger();
3655 TString configstr("");
3656 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3657 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3661 aCTP->SetClassMask(mask);
3662 aCTP->SetClusterMask(clmask);
3665 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3667 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3668 aCTP = fRunLoader->GetTrigger();
3669 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3670 // get inputs from actp - just get
3671 AliESDHeader* esdheader = fesd->GetHeader();
3672 esdheader->SetL0TriggerInputs(aCTP->GetL0TriggerInputs());
3673 esdheader->SetL1TriggerInputs(aCTP->GetL1TriggerInputs());
3674 esdheader->SetL2TriggerInputs(aCTP->GetL2TriggerInputs());
3675 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3678 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3683 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3685 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3686 if (fRawReader) delete aCTP;
3690 UChar_t clustmask = 0;
3692 ULong64_t trmask = fEventInfo.GetTriggerMask();
3693 const TObjArray& classesArray = config->GetClasses();
3694 Int_t nclasses = classesArray.GetEntriesFast();
3695 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3696 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3697 if (trclass && trclass->GetMask()>0) {
3698 Int_t trindex = TMath::Nint(TMath::Log2(trclass->GetMask()));
3699 fesd->SetTriggerClass(trclass->GetName(),trindex);
3700 if (fRawReader) fRawReader->LoadTriggerClass(trclass->GetName(),trindex);
3701 if (trmask & (1ull << trindex)) {
3703 trclasses += trclass->GetName();
3705 clustmask |= trclass->GetCluster()->GetClusterMask();
3709 fEventInfo.SetTriggerClasses(trclasses);
3711 // Write names of active trigger inputs in ESD Header
3712 const TObjArray& inputsArray = config->GetInputs();
3713 Int_t ninputs = inputsArray.GetEntriesFast();
3714 for( Int_t iinput=0; iinput < ninputs; iinput++ ) {
3715 AliTriggerInput* trginput = (AliTriggerInput*)inputsArray.At(iinput);
3716 if (trginput && trginput->GetMask()>0) {
3717 Int_t inputIndex = (Int_t)TMath::Nint(TMath::Log2(trginput->GetMask()));
3718 AliESDHeader* headeresd = fesd->GetHeader();
3719 Int_t trglevel = (Int_t)trginput->GetLevel();
3720 if (trglevel == 0) headeresd->SetActiveTriggerInputs(trginput->GetInputName(), inputIndex);
3721 if (trglevel == 1) headeresd->SetActiveTriggerInputs(trginput->GetInputName(), inputIndex+24);
3722 if (trglevel == 2) headeresd->SetActiveTriggerInputs(trginput->GetInputName(), inputIndex+48);
3726 // Set the information in ESD
3727 fesd->SetTriggerMask(trmask);
3728 fesd->SetTriggerCluster(clustmask);
3730 if (!aCTP->CheckTriggeredDetectors()) {
3731 if (fRawReader) delete aCTP;
3735 if (fRawReader) delete aCTP;
3737 // We have to fill also the HLT decision here!!
3743 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3745 // Match the detector list found in the rec.C or the default 'ALL'
3746 // to the list found in the GRP (stored there by the shuttle PP which
3747 // gets the information from ECS)
3748 static TString resultList;
3749 TString detList = detectorList;
3753 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3754 if ((detectorMask >> iDet) & 0x1) {
3755 TString det = AliDAQ::OfflineModuleName(iDet);
3756 if ((detList.CompareTo("ALL") == 0) ||
3757 ((detList.BeginsWith("ALL ") ||
3758 detList.EndsWith(" ALL") ||
3759 detList.Contains(" ALL ")) &&
3760 !(detList.BeginsWith("-"+det+" ") ||
3761 detList.EndsWith(" -"+det) ||
3762 detList.Contains(" -"+det+" "))) ||
3763 (detList.CompareTo(det) == 0) ||
3764 detList.BeginsWith(det+" ") ||
3765 detList.EndsWith(" "+det) ||
3766 detList.Contains( " "+det+" " )) {
3767 if (!resultList.EndsWith(det + " ")) {
3776 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3777 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3778 if ((detList.CompareTo("ALL") == 0) ||
3779 ((detList.BeginsWith("ALL ") ||
3780 detList.EndsWith(" ALL") ||
3781 detList.Contains(" ALL ")) &&
3782 !(detList.BeginsWith("-"+hltDet+" ") ||
3783 detList.EndsWith(" -"+hltDet) ||
3784 detList.Contains(" -"+hltDet+" "))) ||
3785 (detList.CompareTo(hltDet) == 0) ||
3786 detList.BeginsWith(hltDet+" ") ||
3787 detList.EndsWith(" "+hltDet) ||
3788 detList.Contains( " "+hltDet+" " )) {
3789 resultList += hltDet;
3793 return resultList.Data();
3797 //______________________________________________________________________________
3798 void AliReconstruction::Abort(const char *method, EAbort what)
3800 // Abort processing. If what = kAbortProcess, the Process() loop will be
3801 // aborted. If what = kAbortFile, the current file in a chain will be
3802 // aborted and the processing will continue with the next file, if there
3803 // is no next file then Process() will be aborted. Abort() can also be
3804 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3805 // the SlaveTerminate() and Terminate() are always called. The abort flag
3806 // can be checked in these methods using GetAbort().
3808 // The method is overwritten in AliReconstruction for better handling of
3809 // reco specific errors
3811 if (!fStopOnError) return;
3815 TString whyMess = method;
3816 whyMess += " failed! Aborting...";
3818 AliError(whyMess.Data());
3821 TString mess = "Abort";
3822 if (fAbort == kAbortProcess)
3823 mess = "AbortProcess";
3824 else if (fAbort == kAbortFile)
3827 Info(mess, whyMess.Data());
3830 //______________________________________________________________________________
3831 Bool_t AliReconstruction::ProcessEvent(void* event)
3833 // Method that is used in case the event loop
3834 // is steered from outside, for example by AMORE
3835 // 'event' is a pointer to the DATE event in the memory
3837 if (fRawReader) delete fRawReader;
3838 fRawReader = new AliRawReaderDate(event);
3839 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
3846 //______________________________________________________________________________
3847 Bool_t AliReconstruction::ParseOutput()
3849 // The method parses the output file
3850 // location string in order to steer
3851 // properly the selector
3853 TPMERegexp re1("(\\w+\\.zip#\\w+\\.root):([,*\\w+\\.root,*]+)@dataset://(\\w++)");
3854 TPMERegexp re2("(\\w+\\.root)?@?dataset://(\\w++)");
3856 if (re1.Match(fESDOutput) == 4) {
3857 // root archive with output files stored and regustered
3859 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",re1[1].Data()));
3860 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re1[3].Data()));
3861 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3862 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_ARCHIVE",re1[2].Data()));
3863 AliInfo(Form("%s files will be stored within %s in dataset %s",
3868 else if (re2.Match(fESDOutput) == 3) {
3869 // output file stored and registered
3871 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",(re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data()));
3872 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",re2[2].Data()));
3873 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_DATASET",""));
3874 AliInfo(Form("%s will be stored in dataset %s",
3875 (re2[1].IsNull()) ? "AliESDs.root" : re2[1].Data(),
3879 if (fESDOutput.IsNull()) {
3880 // Output location not given.
3881 // Assuming xrootd has been already started and
3882 // the output file has to be sent back
3883 // to the client machine
3884 TString esdUrl(Form("root://%s/%s/",
3885 TUrl(gSystem->HostName()).GetHostFQDN(),
3887 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE","AliESDs.root"));
3888 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",esdUrl.Data()));
3889 AliInfo(Form("AliESDs.root will be stored in %s",
3893 // User specified an output location.
3894 // Ones has just to parse it here
3895 TUrl outputUrl(fESDOutput.Data());
3896 TString outputFile(gSystem->BaseName(outputUrl.GetFile()));
3897 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.IsNull() ? "AliESDs.root" : outputFile.Data()));
3898 TString outputLocation(outputUrl.GetUrl());
3899 outputLocation.ReplaceAll(outputFile.Data(),"");
3900 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE_LOCATION",outputLocation.Data()));
3901 AliInfo(Form("%s will be stored in %s",
3902 outputFile.IsNull() ? "AliESDs.root" : outputFile.Data(),
3903 outputLocation.Data()));
3910 //______________________________________________________________________________
3911 Bool_t AliReconstruction::IsHighPt() const {
3912 // Selection of events containing "high" pT tracks
3913 // If at least one track is found within 1.5 and 100 GeV (pT)
3914 // that was reconstructed by both ITS and TPC, the event is accepted
3918 const Double_t pTmin = 1.5;
3919 const Double_t pTmax = 100;
3921 mask |= (AliESDtrack::kITSrefit);
3922 mask |= (AliESDtrack::kTPCrefit);
3924 Bool_t isOK = kFALSE;
3926 if (fesd && fesd->GetEventType()==AliRawEventHeaderBase::kPhysicsEvent) {
3927 // Check if this ia a physics event (code 7)
3928 Int_t ntrk = fesd->GetNumberOfTracks();
3929 for (Int_t itrk=0; itrk<ntrk; ++itrk) {
3931 AliESDtrack * trk = fesd->GetTrack(itrk);
3933 && trk->Pt() > pTmin
3934 && trk->Pt() < pTmax
3935 && (trk->GetStatus() & mask) == mask ) {
3945 //______________________________________________________________________________
3946 Bool_t AliReconstruction::IsCosmicOrCalibSpecie() const {
3947 // Select cosmic or calibration events
3949 Bool_t isOK = kFALSE;
3951 if (fesd && fesd->GetEventType()==AliRawEventHeaderBase::kPhysicsEvent) {
3952 // Check if this ia a physics event (code 7)
3954 UInt_t specie = fesd->GetEventSpecie();
3955 if (specie==AliRecoParam::kCosmic || specie==AliRecoParam::kCalib) {
3962 //______________________________________________________________________________
3963 void AliReconstruction::WriteESDfriend() {
3964 // Fill the ESD friend in the tree. The required fraction of ESD friends is stored
3965 // in fFractionFriends. We select events where we store the ESD friends according
3966 // to the following algorithm:
3967 // 1. Store all Cosmic or Calibration events within the required fraction
3968 // 2. Sample "high Pt" events within the remaining fraction after step 1.
3969 // 3. Sample randomly events if we still have remaining slot
3973 Bool_t isSelected = kFALSE;
3975 if (IsCosmicOrCalibSpecie()) { // Selection of calib or cosmic events
3977 Double_t curentSpecieFraction = ((Double_t)(fNspecie+1))/((Double_t)(fNall+1));
3978 // "Bayesian" estimate supposing that without events all the events are of the required type
3980 Double_t rnd = gRandom->Rndm()*curentSpecieFraction;
3981 if (rnd<fFractionFriends) {
3987 Double_t remainingFraction = fFractionFriends;
3988 remainingFraction -= ((Double_t)(fSspecie)/(Double_t)(fNall));
3990 if (IsHighPt()) { // Selection of "high Pt" events
3992 Double_t curentHighPtFraction = ((Double_t)(fNhighPt+1))/((Double_t)(fNall+1));
3993 // "Bayesian" estimate supposing that without events all the events are of the required type
3996 Double_t rnd = gRandom->Rndm()*curentHighPtFraction;
3997 if (rnd<remainingFraction) {
4003 remainingFraction -= ((Double_t)(fShighPt)/(Double_t)(fNall));
4005 // Random selection to fill the remaining fraction (if any)
4007 Double_t rnd = gRandom->Rndm();
4008 if (rnd<remainingFraction) {
4014 fesdf->~AliESDfriend();
4015 new (fesdf) AliESDfriend(); // Reset...
4016 fesdf->SetSkipBit(kTRUE);