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("..."); //
108 // For debug purposes the method SetCheckPointLevel can be used. If the //
109 // argument is greater than 0, files with ESD events will be written after //
110 // selected steps of the reconstruction for each event: //
111 // level 1: after tracking and after filling of ESD (final) //
112 // level 2: in addition after each tracking step //
113 // level 3: in addition after the filling of ESD for each detector //
114 // If a final check point file exists for an event, this event will be //
115 // skipped in the reconstruction. The tracking and the filling of ESD for //
116 // a detector will be skipped as well, if the corresponding check point //
117 // file exists. The ESD event will then be loaded from the file instead. //
119 ///////////////////////////////////////////////////////////////////////////////
126 #include <TPluginManager.h>
127 #include <TGeoManager.h>
128 #include <TLorentzVector.h>
131 #include <TObjArray.h>
134 #include "AliReconstruction.h"
135 #include "AliCodeTimer.h"
136 #include "AliReconstructor.h"
138 #include "AliRunLoader.h"
140 #include "AliRawReaderFile.h"
141 #include "AliRawReaderDate.h"
142 #include "AliRawReaderRoot.h"
143 #include "AliRawEventHeaderBase.h"
144 #include "AliESDEvent.h"
145 #include "AliESDMuonTrack.h"
146 #include "AliESDfriend.h"
147 #include "AliESDVertex.h"
148 #include "AliESDcascade.h"
149 #include "AliESDkink.h"
150 #include "AliESDtrack.h"
151 #include "AliESDCaloCluster.h"
152 #include "AliESDCaloCells.h"
153 #include "AliMultiplicity.h"
154 #include "AliTracker.h"
155 #include "AliVertexer.h"
156 #include "AliVertexerTracks.h"
157 #include "AliV0vertexer.h"
158 #include "AliCascadeVertexer.h"
159 #include "AliHeader.h"
160 #include "AliGenEventHeader.h"
162 #include "AliESDpid.h"
163 #include "AliESDtrack.h"
164 #include "AliESDPmdTrack.h"
166 #include "AliESDTagCreator.h"
167 #include "AliAODTagCreator.h"
169 #include "AliGeomManager.h"
170 #include "AliTrackPointArray.h"
171 #include "AliCDBManager.h"
172 #include "AliCDBStorage.h"
173 #include "AliCDBEntry.h"
174 #include "AliAlignObj.h"
176 #include "AliCentralTrigger.h"
177 #include "AliTriggerConfiguration.h"
178 #include "AliTriggerClass.h"
179 #include "AliCTPRawStream.h"
181 #include "AliQADataMakerRec.h"
182 #include "AliGlobalQADataMaker.h"
184 #include "AliQADataMakerSteer.h"
186 #include "AliPlaneEff.h"
188 #include "AliSysInfo.h" // memory snapshots
189 #include "AliRawHLTManager.h"
191 #include "AliMagWrapCheb.h"
193 ClassImp(AliReconstruction)
196 //_____________________________________________________________________________
197 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
199 //_____________________________________________________________________________
200 AliReconstruction::AliReconstruction(const char* gAliceFilename,
201 const char* name, const char* title) :
204 fUniformField(kFALSE),
205 fForcedFieldMap(0x0),
206 fRunVertexFinder(kTRUE),
207 fRunVertexFinderTracks(kTRUE),
208 fRunHLTTracking(kFALSE),
209 fRunMuonTracking(kFALSE),
211 fRunCascadeFinder(kTRUE),
212 fStopOnError(kFALSE),
213 fWriteAlignmentData(kFALSE),
214 fWriteESDfriend(kFALSE),
216 fFillTriggerESD(kTRUE),
224 fRunLocalReconstruction("ALL"),
227 fUseTrackingErrorsForAlignment(""),
228 fGAliceFileName(gAliceFilename),
233 fNumberOfEventsPerFile(1),
236 fLoadAlignFromCDB(kTRUE),
237 fLoadAlignData("ALL"),
243 fParentRawReader(NULL),
246 fDiamondProfile(NULL),
247 fDiamondProfileTPC(NULL),
248 fMeanVertexConstraint(kTRUE),
252 fAlignObjArray(NULL),
255 fInitCDBCalled(kFALSE),
256 fSetRunNumberFromDataCalled(kFALSE),
262 fSameQACycle(kFALSE),
264 fRunPlaneEff(kFALSE),
276 fIsNewRunLoader(kFALSE),
279 // create reconstruction object with default parameters
281 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
282 fReconstructor[iDet] = NULL;
283 fLoader[iDet] = NULL;
284 fTracker[iDet] = NULL;
285 fQADataMaker[iDet] = NULL;
286 fQACycles[iDet] = 999999;
288 fQADataMaker[fgkNDetectors]=NULL; //Global QA
292 //_____________________________________________________________________________
293 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
296 fUniformField(rec.fUniformField),
297 fForcedFieldMap(0x0),
298 fRunVertexFinder(rec.fRunVertexFinder),
299 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
300 fRunHLTTracking(rec.fRunHLTTracking),
301 fRunMuonTracking(rec.fRunMuonTracking),
302 fRunV0Finder(rec.fRunV0Finder),
303 fRunCascadeFinder(rec.fRunCascadeFinder),
304 fStopOnError(rec.fStopOnError),
305 fWriteAlignmentData(rec.fWriteAlignmentData),
306 fWriteESDfriend(rec.fWriteESDfriend),
307 fWriteAOD(rec.fWriteAOD),
308 fFillTriggerESD(rec.fFillTriggerESD),
310 fCleanESD(rec.fCleanESD),
311 fV0DCAmax(rec.fV0DCAmax),
312 fV0CsPmin(rec.fV0CsPmin),
316 fRunLocalReconstruction(rec.fRunLocalReconstruction),
317 fRunTracking(rec.fRunTracking),
318 fFillESD(rec.fFillESD),
319 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
320 fGAliceFileName(rec.fGAliceFileName),
322 fEquipIdMap(rec.fEquipIdMap),
323 fFirstEvent(rec.fFirstEvent),
324 fLastEvent(rec.fLastEvent),
325 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
328 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
329 fLoadAlignData(rec.fLoadAlignData),
330 fESDPar(rec.fESDPar),
331 fUseHLTData(rec.fUseHLTData),
335 fParentRawReader(NULL),
338 fDiamondProfile(NULL),
339 fDiamondProfileTPC(NULL),
340 fMeanVertexConstraint(rec.fMeanVertexConstraint),
344 fAlignObjArray(rec.fAlignObjArray),
345 fCDBUri(rec.fCDBUri),
347 fInitCDBCalled(rec.fInitCDBCalled),
348 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
349 fQADetectors(rec.fQADetectors),
350 fQATasks(rec.fQATasks),
352 fRunGlobalQA(rec.fRunGlobalQA),
353 fInLoopQA(rec.fInLoopQA),
354 fSameQACycle(rec.fSameQACycle),
355 fRunPlaneEff(rec.fRunPlaneEff),
367 fIsNewRunLoader(rec.fIsNewRunLoader),
372 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
373 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
375 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
376 fReconstructor[iDet] = NULL;
377 fLoader[iDet] = NULL;
378 fTracker[iDet] = NULL;
379 fQADataMaker[iDet] = NULL;
380 fQACycles[iDet] = rec.fQACycles[iDet];
382 fQADataMaker[fgkNDetectors]=NULL; //Global QA
383 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
384 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
387 fForcedFieldMap=new AliMagWrapCheb(*((AliMagWrapCheb*)rec.fForcedFieldMap));
390 //_____________________________________________________________________________
391 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
393 // assignment operator
395 this->~AliReconstruction();
396 new(this) AliReconstruction(rec);
400 //_____________________________________________________________________________
401 AliReconstruction::~AliReconstruction()
407 fSpecCDBUri.Delete();
408 delete fForcedFieldMap;
410 AliCodeTimer::Instance()->Print();
413 //_____________________________________________________________________________
414 void AliReconstruction::InitCDB()
416 // activate a default CDB storage
417 // First check if we have any CDB storage set, because it is used
418 // to retrieve the calibration and alignment constants
420 if (fInitCDBCalled) return;
421 fInitCDBCalled = kTRUE;
423 AliCDBManager* man = AliCDBManager::Instance();
424 if (man->IsDefaultStorageSet())
426 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
427 AliWarning("Default CDB storage has been already set !");
428 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
429 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
430 fCDBUri = man->GetDefaultStorage()->GetURI();
433 if (fCDBUri.Length() > 0)
435 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
436 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
437 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
439 fCDBUri="local://$ALICE_ROOT";
440 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
441 AliWarning("Default CDB storage not yet set !!!!");
442 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
443 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
446 man->SetDefaultStorage(fCDBUri);
449 // Now activate the detector specific CDB storage locations
450 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
451 TObject* obj = fSpecCDBUri[i];
453 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
454 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
455 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
456 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
461 //_____________________________________________________________________________
462 void AliReconstruction::SetDefaultStorage(const char* uri) {
463 // Store the desired default CDB storage location
464 // Activate it later within the Run() method
470 //_____________________________________________________________________________
471 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
472 // Store a detector-specific CDB storage location
473 // Activate it later within the Run() method
475 AliCDBPath aPath(calibType);
476 if(!aPath.IsValid()){
477 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
478 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
479 if(!strcmp(calibType, fgkDetectorName[iDet])) {
480 aPath.SetPath(Form("%s/*", calibType));
481 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
485 if(!aPath.IsValid()){
486 AliError(Form("Not a valid path or detector: %s", calibType));
491 // // check that calibType refers to a "valid" detector name
492 // Bool_t isDetector = kFALSE;
493 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
494 // TString detName = fgkDetectorName[iDet];
495 // if(aPath.GetLevel0() == detName) {
496 // isDetector = kTRUE;
502 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
506 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
507 if (obj) fSpecCDBUri.Remove(obj);
508 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
512 //_____________________________________________________________________________
513 Bool_t AliReconstruction::SetRunNumberFromData()
515 // The method is called in Run() in order
516 // to set a correct run number.
517 // In case of raw data reconstruction the
518 // run number is taken from the raw data header
520 if (fSetRunNumberFromDataCalled) return kTRUE;
521 fSetRunNumberFromDataCalled = kTRUE;
523 AliCDBManager* man = AliCDBManager::Instance();
525 if(man->GetRun() > 0) {
526 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
530 AliError("No run loader is found !");
533 // read run number from gAlice
534 if(fRunLoader->GetAliRun())
535 AliCDBManager::Instance()->SetRun(fRunLoader->GetHeader()->GetRun());
538 if(fRawReader->NextEvent()) {
539 AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
540 fRawReader->RewindEvents();
543 if(man->GetRun() > 0) {
544 AliWarning("No raw events is found ! Using settings in AliCDBManager !");
549 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
555 AliError("Neither gAlice nor RawReader objects are found !");
565 //_____________________________________________________________________________
566 void AliReconstruction::SetCDBLock() {
567 // Set CDB lock: from now on it is forbidden to reset the run number
568 // or the default storage or to activate any further storage!
570 AliCDBManager::Instance()->SetLock(1);
573 //_____________________________________________________________________________
574 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
576 // Read the alignment objects from CDB.
577 // Each detector is supposed to have the
578 // alignment objects in DET/Align/Data CDB path.
579 // All the detector objects are then collected,
580 // sorted by geometry level (starting from ALIC) and
581 // then applied to the TGeo geometry.
582 // Finally an overlaps check is performed.
584 // Load alignment data from CDB and fill fAlignObjArray
585 if(fLoadAlignFromCDB){
587 TString detStr = detectors;
588 TString loadAlObjsListOfDets = "";
590 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
591 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
592 loadAlObjsListOfDets += fgkDetectorName[iDet];
593 loadAlObjsListOfDets += " ";
594 } // end loop over detectors
595 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
596 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
598 // Check if the array with alignment objects was
599 // provided by the user. If yes, apply the objects
600 // to the present TGeo geometry
601 if (fAlignObjArray) {
602 if (gGeoManager && gGeoManager->IsClosed()) {
603 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
604 AliError("The misalignment of one or more volumes failed!"
605 "Compare the list of simulated detectors and the list of detector alignment data!");
610 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
616 delete fAlignObjArray; fAlignObjArray=0;
621 //_____________________________________________________________________________
622 void AliReconstruction::SetGAliceFile(const char* fileName)
624 // set the name of the galice file
626 fGAliceFileName = fileName;
629 //_____________________________________________________________________________
630 void AliReconstruction::SetInput(const char* input)
632 // In case the input string starts with 'mem://', we run in an online mode
633 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
634 // file is assumed. One can give as an input:
635 // mem://: - events taken from DAQ monitoring libs online
637 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
641 //_____________________________________________________________________________
642 void AliReconstruction::SetOption(const char* detector, const char* option)
644 // set options for the reconstruction of a detector
646 TObject* obj = fOptions.FindObject(detector);
647 if (obj) fOptions.Remove(obj);
648 fOptions.Add(new TNamed(detector, option));
651 //_____________________________________________________________________________
652 Bool_t AliReconstruction::SetFieldMap(Float_t l3Current, Float_t diCurrent, Float_t factor, const char *path) {
653 //------------------------------------------------
654 // The magnetic field map, defined externally...
655 // L3 current 30000 A -> 0.5 T
656 // L3 current 12000 A -> 0.2 T
657 // dipole current 6000 A
658 // The polarities must be the same
659 //------------------------------------------------
660 const Float_t l3NominalCurrent1=30000.; // (A)
661 const Float_t l3NominalCurrent2=12000.; // (A)
662 const Float_t diNominalCurrent =6000. ; // (A)
664 const Float_t tolerance=0.03; // relative current tolerance
665 const Float_t zero=77.; // "zero" current (A)
668 Bool_t dipoleON=kFALSE;
670 TString s=(factor < 0) ? "L3: -" : "L3: +";
672 if (TMath::Abs(l3Current-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
673 map=AliMagWrapCheb::k5kG;
676 if (TMath::Abs(l3Current-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
677 map=AliMagWrapCheb::k2kG;
680 if (TMath::Abs(l3Current) < zero) {
681 map=AliMagWrapCheb::k2kG;
683 factor=0.; // in fact, this is a global factor...
685 AliError("Wrong L3 current !");
689 if (TMath::Abs(diCurrent-diNominalCurrent)/diNominalCurrent < tolerance) {
690 // 3% current tolerance...
694 if (TMath::Abs(diCurrent) < zero) { // some small current..
698 AliError("Wrong dipole current !");
702 delete fForcedFieldMap;
704 new AliMagWrapCheb("B field map ",s,2,factor,10.,map,dipoleON,path);
706 fForcedFieldMap->Print();
708 AliTracker::SetFieldMap(fForcedFieldMap,fUniformField);
714 Bool_t AliReconstruction::InitGRP() {
715 //------------------------------------
716 // Initialization of the GRP entry
717 //------------------------------------
718 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
720 if (entry) fGRPData = dynamic_cast<TMap*>(entry->GetObject());
723 AliError("No GRP entry found in OCDB!");
728 //*** Dealing with the magnetic field map
729 if (AliTracker::GetFieldMap()) {
730 AliInfo("Running with the externally set B field !");
732 // Construct the field map out of the information retrieved from GRP.
737 TObjString *l3Current=
738 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
740 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
743 TObjString *l3Polarity=
744 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
746 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
751 TObjString *diCurrent=
752 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
754 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
757 TObjString *diPolarity=
758 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
760 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
765 Float_t l3Cur=TMath::Abs(atof(l3Current->GetName()));
766 Float_t diCur=TMath::Abs(atof(diCurrent->GetName()));
767 Float_t l3Pol=atof(l3Polarity->GetName());
769 if (l3Pol != 0.) factor=-1.;
772 if (!SetFieldMap(l3Cur, diCur, factor)) {
773 AliFatal("Failed to creat a B field map ! Exiting...");
775 AliInfo("Running with the B field constructed out of GRP !");
778 AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
784 //*** Get the diamond profile from OCDB
785 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
787 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
789 AliError("No diamond profile found in OCDB!");
792 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
794 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
796 AliError("No diamond profile found in OCDB!");
802 //_____________________________________________________________________________
803 Bool_t AliReconstruction::Run(const char* input)
806 AliCodeTimerAuto("");
808 if (!InitRun(input)) return kFALSE;
809 //******* The loop over events
811 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
812 (fRawReader && fRawReader->NextEvent())) {
813 if (!RunEvent(iEvent)) return kFALSE;
817 if (!FinishRun()) return kFALSE;
822 //_____________________________________________________________________________
823 Bool_t AliReconstruction::InitRun(const char* input)
825 // Initialize all the stuff before
826 // going into the event loop
827 // If the second argument is given, the first one is ignored and
828 // the reconstruction works in an online mode
829 AliCodeTimerAuto("");
831 // Overwrite the previous setting
832 if (input) fInput = input;
834 // set the input in case of raw data
835 fRawReader = AliRawReader::Create(fInput.Data());
837 AliInfo("Reconstruction will run over digits");
839 if (!fEquipIdMap.IsNull() && fRawReader)
840 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
842 if (!fUseHLTData.IsNull()) {
843 // create the RawReaderHLT which performs redirection of HLT input data for
844 // the specified detectors
845 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
847 fParentRawReader=fRawReader;
848 fRawReader=pRawReader;
850 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
854 AliSysInfo::AddStamp("Start");
855 // get the run loader
856 if (!InitRunLoader()) return kFALSE;
857 AliSysInfo::AddStamp("LoadLoader");
859 // Initialize the CDB storage
862 AliSysInfo::AddStamp("LoadCDB");
864 // Set run number in CDBManager (if it is not already set by the user)
865 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
867 // Set CDB lock: from now on it is forbidden to reset the run number
868 // or the default storage or to activate any further storage!
871 // Import ideal TGeo geometry and apply misalignment
873 TString geom(gSystem->DirName(fGAliceFileName));
874 geom += "/geometry.root";
875 AliGeomManager::LoadGeometry(geom.Data());
877 TString detsToCheck=fRunLocalReconstruction;
878 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data()))
879 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
880 if (!gGeoManager) if (fStopOnError) return kFALSE;
883 if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
884 AliSysInfo::AddStamp("LoadGeom");
887 if (!InitGRP()) return kFALSE;
890 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
891 if(fDiamondProfile && fMeanVertexConstraint) ftVertexer->SetVtxStart(fDiamondProfile);
894 if (fRunVertexFinder && !CreateVertexer()) {
900 AliSysInfo::AddStamp("Vertexer");
903 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
909 AliSysInfo::AddStamp("LoadTrackers");
911 // get the possibly already existing ESD file and tree
912 fesd = new AliESDEvent(); fhltesd = new AliESDEvent();
913 if (!gSystem->AccessPathName("AliESDs.root")){
914 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
915 ffileOld = TFile::Open("AliESDs.old.root");
916 if (ffileOld && ffileOld->IsOpen()) {
917 ftreeOld = (TTree*) ffileOld->Get("esdTree");
918 if (ftreeOld)fesd->ReadFromTree(ftreeOld);
919 fhlttreeOld = (TTree*) ffileOld->Get("HLTesdTree");
920 if (fhlttreeOld) fhltesd->ReadFromTree(fhlttreeOld);
924 // create the ESD output file and tree
925 ffile = TFile::Open("AliESDs.root", "RECREATE");
926 ffile->SetCompressionLevel(2);
927 if (!ffile->IsOpen()) {
928 AliError("opening AliESDs.root failed");
929 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
932 ftree = new TTree("esdTree", "Tree with ESD objects");
933 fesd = new AliESDEvent();
934 fesd->CreateStdContent();
935 fesd->WriteToTree(ftree);
937 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
938 fhltesd = new AliESDEvent();
939 fhltesd->CreateStdContent();
940 fhltesd->WriteToTree(fhlttree);
943 if (fWriteESDfriend) {
944 fesdf = new AliESDfriend();
945 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
946 br->SetFile("AliESDfriends.root");
947 fesd->AddObject(fesdf);
952 if (fRawReader) fRawReader->RewindEvents();
955 gSystem->GetProcInfo(&ProcInfo);
956 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
959 if (fRunQA && fRawReader && fQATasks.Contains(Form("%d", AliQA::kRAWS))) {
960 AliQADataMakerSteer qas ;
961 qas.Run(fRunLocalReconstruction, fRawReader) ;
962 fSameQACycle = kTRUE ;
964 //Initialize the QA and start of cycle for out-of-cycle QA
966 TString detStr(fQADetectors) ;
967 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
968 if (!IsSelected(fgkDetectorName[iDet], detStr))
970 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
973 AliInfo(Form("Initializing the QA data maker for %s",
974 fgkDetectorName[iDet]));
975 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
976 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
977 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
978 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
980 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
981 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
982 fSameQACycle = kTRUE;
984 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
985 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle);
986 fSameQACycle = kTRUE;
991 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
992 AliInfo(Form("Initializing the global QA data maker"));
993 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
995 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
996 AliTracker::SetResidualsArray(arr);
998 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
999 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
1002 fSameQACycle = kFALSE;
1003 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1004 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
1005 fSameQACycle = kTRUE;
1007 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1008 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle);
1009 fSameQACycle = kTRUE;
1015 //Initialize the Plane Efficiency framework
1016 if (fRunPlaneEff && !InitPlaneEff()) {
1017 if(fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1020 if (strcmp(gProgName,"alieve") == 0)
1021 fRunAliEVE = InitAliEVE();
1026 //_____________________________________________________________________________
1027 Bool_t AliReconstruction::RunEvent(Int_t iEvent)
1029 // run the reconstruction over a single event
1030 // The event loop is steered in Run method
1032 AliCodeTimerAuto("");
1034 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1035 fRunLoader->SetEventNumber(iEvent);
1036 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1038 //?? fRunLoader->MakeTree("H");
1039 fRunLoader->TreeE()->Fill();
1042 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1043 // copy old ESD to the new one
1045 fesd->ReadFromTree(ftreeOld);
1046 ftreeOld->GetEntry(iEvent);
1050 fhltesd->ReadFromTree(fhlttreeOld);
1051 fhlttreeOld->GetEntry(iEvent);
1057 AliInfo(Form("processing event %d", iEvent));
1059 //Start of cycle for the in-loop QA
1062 fSameQACycle = kFALSE ;
1063 TString detStr(fQADetectors);
1064 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1065 if (!IsSelected(fgkDetectorName[iDet], detStr))
1067 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
1070 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1071 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
1072 fSameQACycle = kTRUE;
1074 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1075 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle) ;
1076 fSameQACycle = kTRUE;
1080 fSameQACycle = kFALSE;
1081 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1082 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1083 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
1084 fSameQACycle = kTRUE;
1086 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1087 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle);
1088 fSameQACycle = kTRUE;
1094 fRunLoader->GetEvent(iEvent);
1096 char aFileName[256];
1097 sprintf(aFileName, "ESD_%d.%d_final.root",
1098 fRunLoader->GetHeader()->GetRun(),
1099 fRunLoader->GetHeader()->GetEventNrInRun());
1100 if (!gSystem->AccessPathName(aFileName)) return kTRUE;
1102 // local single event reconstruction
1103 if (!fRunLocalReconstruction.IsNull()) {
1104 TString detectors=fRunLocalReconstruction;
1105 // run HLT event reconstruction first
1106 // ;-( IsSelected changes the string
1107 if (IsSelected("HLT", detectors) &&
1108 !RunLocalEventReconstruction("HLT")) {
1109 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1111 detectors=fRunLocalReconstruction;
1112 detectors.ReplaceAll("HLT", "");
1113 if (!RunLocalEventReconstruction(detectors)) {
1114 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1118 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1119 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1120 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1121 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1123 // Set magnetic field from the tracker
1124 fesd->SetMagneticField(AliTracker::GetBz());
1125 fhltesd->SetMagneticField(AliTracker::GetBz());
1129 // Fill raw-data error log into the ESD
1130 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1133 if (fRunVertexFinder) {
1134 if (!ReadESD(fesd, "vertex")) {
1135 if (!RunVertexFinder(fesd)) {
1136 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1138 if (fCheckPointLevel > 0) WriteESD(fesd, "vertex");
1143 if (!fRunTracking.IsNull()) {
1144 if (fRunMuonTracking) {
1145 if (!RunMuonTracking(fesd)) {
1146 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1152 if (!fRunTracking.IsNull()) {
1153 if (!ReadESD(fesd, "tracking")) {
1154 if (!RunTracking(fesd)) {
1155 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1157 if (fCheckPointLevel > 0) WriteESD(fesd, "tracking");
1162 if (!fFillESD.IsNull()) {
1163 TString detectors=fFillESD;
1164 // run HLT first and on hltesd
1165 // ;-( IsSelected changes the string
1166 if (IsSelected("HLT", detectors) &&
1167 !FillESD(fhltesd, "HLT")) {
1168 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1171 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1172 if (detectors.Contains("ALL")) {
1174 for (Int_t idet=0; idet<fgkNDetectors; ++idet){
1175 detectors += fgkDetectorName[idet];
1179 detectors.ReplaceAll("HLT", "");
1180 if (!FillESD(fesd, detectors)) {
1181 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1185 // fill Event header information from the RawEventHeader
1186 if (fRawReader){FillRawEventHeaderESD(fesd);}
1189 AliESDpid::MakePID(fesd);
1190 if (fCheckPointLevel > 1) WriteESD(fesd, "PID");
1192 if (fFillTriggerESD) {
1193 if (!ReadESD(fesd, "trigger")) {
1194 if (!FillTriggerESD(fesd)) {
1195 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1197 if (fCheckPointLevel > 1) WriteESD(fesd, "trigger");
1204 // Propagate track to the beam pipe (if not already done by ITS)
1206 const Int_t ntracks = fesd->GetNumberOfTracks();
1207 const Double_t kBz = fesd->GetMagneticField();
1208 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1211 UShort_t *selectedIdx=new UShort_t[ntracks];
1213 for (Int_t itrack=0; itrack<ntracks; itrack++){
1214 const Double_t kMaxStep = 5; //max step over the material
1217 AliESDtrack *track = fesd->GetTrack(itrack);
1218 if (!track) continue;
1220 AliExternalTrackParam *tpcTrack =
1221 (AliExternalTrackParam *)track->GetTPCInnerParam();
1225 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
1228 Int_t n=trkArray.GetEntriesFast();
1229 selectedIdx[n]=track->GetID();
1230 trkArray.AddLast(tpcTrack);
1233 //Tracks refitted by ITS should already be at the SPD vertex
1234 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1237 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1238 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1243 // Improve the reconstructed primary vertex position using the tracks
1245 TObject *obj = fOptions.FindObject("ITS");
1247 TString optITS = obj->GetTitle();
1248 if (optITS.Contains("cosmics") || optITS.Contains("COSMICS"))
1249 fRunVertexFinderTracks=kFALSE;
1251 if (fRunVertexFinderTracks) {
1252 // TPC + ITS primary vertex
1253 ftVertexer->SetITSrefitRequired();
1254 if(fDiamondProfile && fMeanVertexConstraint) {
1255 ftVertexer->SetVtxStart(fDiamondProfile);
1257 ftVertexer->SetConstraintOff();
1259 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1261 if (pvtx->GetStatus()) {
1262 fesd->SetPrimaryVertex(pvtx);
1263 for (Int_t i=0; i<ntracks; i++) {
1264 AliESDtrack *t = fesd->GetTrack(i);
1265 t->RelateToVertex(pvtx, kBz, kVeryBig);
1270 // TPC-only primary vertex
1271 ftVertexer->SetITSrefitNotRequired();
1272 if(fDiamondProfileTPC && fMeanVertexConstraint) {
1273 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1275 ftVertexer->SetConstraintOff();
1277 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1279 if (pvtx->GetStatus()) {
1280 fesd->SetPrimaryVertexTPC(pvtx);
1281 for (Int_t i=0; i<ntracks; i++) {
1282 AliESDtrack *t = fesd->GetTrack(i);
1283 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1289 delete[] selectedIdx;
1291 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1296 AliV0vertexer vtxer;
1297 vtxer.Tracks2V0vertices(fesd);
1299 if (fRunCascadeFinder) {
1301 AliCascadeVertexer cvtxer;
1302 cvtxer.V0sTracks2CascadeVertices(fesd);
1307 if (fCleanESD) CleanESD(fesd);
1311 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1312 if (qadm && fQATasks.Contains(Form("%d", AliQA::kESDS)))
1313 qadm->Exec(AliQA::kESDS, fesd);
1317 if (fWriteESDfriend) {
1318 fesdf->~AliESDfriend();
1319 new (fesdf) AliESDfriend(); // Reset...
1320 fesd->GetESDfriend(fesdf);
1328 if (fRunAliEVE) RunAliEVE();
1330 if (fCheckPointLevel > 0) WriteESD(fesd, "final");
1333 if (fWriteESDfriend) {
1334 fesdf->~AliESDfriend();
1335 new (fesdf) AliESDfriend(); // Reset...
1338 ProcInfo_t ProcInfo;
1339 gSystem->GetProcInfo(&ProcInfo);
1340 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1343 // End of cycle for the in-loop
1347 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1348 if (!IsSelected(fgkDetectorName[iDet], fQADetectors))
1350 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1353 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1354 qadm->EndOfCycle(AliQA::kRECPOINTS);
1355 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1356 qadm->EndOfCycle(AliQA::kESDS);
1361 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1363 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1364 qadm->EndOfCycle(AliQA::kRECPOINTS);
1365 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1366 qadm->EndOfCycle(AliQA::kESDS);
1375 //_____________________________________________________________________________
1376 Bool_t AliReconstruction::FinishRun()
1379 // Called after the exit
1380 // from the event loop
1381 AliCodeTimerAuto("");
1383 if (fIsNewRunLoader) { // galice.root didn't exist
1384 fRunLoader->WriteHeader("OVERWRITE");
1385 fRunLoader->CdGAFile();
1386 fRunLoader->Write(0, TObject::kOverwrite);
1389 ftree->GetUserInfo()->Add(fesd);
1390 fhlttree->GetUserInfo()->Add(fhltesd);
1392 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1393 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1395 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1396 cdbMapCopy->SetOwner(1);
1397 cdbMapCopy->SetName("cdbMap");
1398 TIter iter(cdbMap->GetTable());
1401 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1402 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1403 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1404 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1407 TList *cdbListCopy = new TList();
1408 cdbListCopy->SetOwner(1);
1409 cdbListCopy->SetName("cdbList");
1411 TIter iter2(cdbList);
1414 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1415 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1418 ftree->GetUserInfo()->Add(cdbMapCopy);
1419 ftree->GetUserInfo()->Add(cdbListCopy);
1422 if(fESDPar.Contains("ESD.par")){
1423 AliInfo("Attaching ESD.par to Tree");
1424 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1425 ftree->GetUserInfo()->Add(fn);
1431 if (fWriteESDfriend)
1432 ftree->SetBranchStatus("ESDfriend*",0);
1433 // we want to have only one tree version number
1434 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1437 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1438 if (fRunPlaneEff && !FinishPlaneEff()) {
1439 AliWarning("Finish PlaneEff evaluation failed");
1443 CleanUp(ffile, ffileOld);
1446 AliWarning("AOD creation not supported anymore during reconstruction. See ANALYSIS/AliAnalysisTaskESDfilter.cxx instead.");
1449 // Create tags for the events in the ESD tree (the ESD tree is always present)
1450 // In case of empty events the tags will contain dummy values
1451 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1452 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1454 AliWarning("AOD tag creation not supported anymore during reconstruction.");
1457 //Finish QA and end of cycle for out-of-loop QA
1460 AliQADataMakerSteer qas;
1461 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1462 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS, fSameQACycle);
1464 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1465 qas.Run(fRunLocalReconstruction.Data(), AliQA::kESDS, fSameQACycle);
1467 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1469 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1470 qadm->EndOfCycle(AliQA::kRECPOINTS);
1471 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1472 qadm->EndOfCycle(AliQA::kESDS);
1479 // Cleanup of CDB manager: cache and active storages!
1480 AliCDBManager::Instance()->ClearCache();
1486 //_____________________________________________________________________________
1487 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1489 // run the local reconstruction
1490 static Int_t eventNr=0;
1491 AliCodeTimerAuto("")
1493 // AliCDBManager* man = AliCDBManager::Instance();
1494 // Bool_t origCache = man->GetCacheFlag();
1496 // TString detStr = detectors;
1497 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1498 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1499 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1500 // if (!reconstructor) continue;
1501 // if (reconstructor->HasLocalReconstruction()) continue;
1503 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1504 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1506 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1507 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1509 // man->SetCacheFlag(kTRUE);
1510 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1511 // man->GetAll(calibPath); // entries are cached!
1513 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1515 // if (fRawReader) {
1516 // fRawReader->RewindEvents();
1517 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1519 // reconstructor->Reconstruct(fRunLoader);
1522 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1523 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1525 // // unload calibration data
1526 // man->UnloadFromCache(calibPath);
1527 // //man->ClearCache();
1530 // man->SetCacheFlag(origCache);
1532 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1533 // AliError(Form("the following detectors were not found: %s",
1535 // if (fStopOnError) return kFALSE;
1542 //_____________________________________________________________________________
1543 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1545 // run the local reconstruction
1547 static Int_t eventNr=0;
1548 AliCodeTimerAuto("")
1550 TString detStr = detectors;
1551 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1552 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1553 AliReconstructor* reconstructor = GetReconstructor(iDet);
1554 if (!reconstructor) continue;
1555 AliLoader* loader = fLoader[iDet];
1556 // Matthias April 2008: temporary fix to run HLT reconstruction
1557 // although the HLT loader is missing
1558 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1560 reconstructor->Reconstruct(fRawReader, NULL);
1563 reconstructor->Reconstruct(dummy, NULL);
1568 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1571 // conversion of digits
1572 if (fRawReader && reconstructor->HasDigitConversion()) {
1573 AliInfo(Form("converting raw data digits into root objects for %s",
1574 fgkDetectorName[iDet]));
1575 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1576 fgkDetectorName[iDet]));
1577 loader->LoadDigits("update");
1578 loader->CleanDigits();
1579 loader->MakeDigitsContainer();
1580 TTree* digitsTree = loader->TreeD();
1581 reconstructor->ConvertDigits(fRawReader, digitsTree);
1582 loader->WriteDigits("OVERWRITE");
1583 loader->UnloadDigits();
1585 // local reconstruction
1586 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1587 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1588 loader->LoadRecPoints("update");
1589 loader->CleanRecPoints();
1590 loader->MakeRecPointsContainer();
1591 TTree* clustersTree = loader->TreeR();
1592 if (fRawReader && !reconstructor->HasDigitConversion()) {
1593 reconstructor->Reconstruct(fRawReader, clustersTree);
1595 loader->LoadDigits("read");
1596 TTree* digitsTree = loader->TreeD();
1598 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1599 if (fStopOnError) return kFALSE;
1601 reconstructor->Reconstruct(digitsTree, clustersTree);
1603 loader->UnloadDigits();
1606 // In-loop QA for local reconstrucion
1607 if (fRunQA && fInLoopQA) {
1608 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1611 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1613 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1615 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1616 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1618 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1622 loader->WriteRecPoints("OVERWRITE");
1623 loader->UnloadRecPoints();
1624 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1627 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1628 AliError(Form("the following detectors were not found: %s",
1630 if (fStopOnError) return kFALSE;
1636 //_____________________________________________________________________________
1637 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1639 // run the barrel tracking
1641 AliCodeTimerAuto("")
1643 AliESDVertex* vertex = NULL;
1644 Double_t vtxPos[3] = {0, 0, 0};
1645 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1646 TArrayF mcVertex(3);
1647 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1648 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1649 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1653 AliInfo("running the ITS vertex finder");
1655 fLoader[0]->LoadRecPoints();
1656 TTree* cltree = fLoader[0]->TreeR();
1658 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1659 vertex = fVertexer->FindVertexForCurrentEvent(cltree);
1662 AliError("Can't get the ITS cluster tree");
1664 fLoader[0]->UnloadRecPoints();
1667 AliError("Can't get the ITS loader");
1670 AliWarning("Vertex not found");
1671 vertex = new AliESDVertex();
1672 vertex->SetName("default");
1675 vertex->SetName("reconstructed");
1679 AliInfo("getting the primary vertex from MC");
1680 vertex = new AliESDVertex(vtxPos, vtxErr);
1684 vertex->GetXYZ(vtxPos);
1685 vertex->GetSigmaXYZ(vtxErr);
1687 AliWarning("no vertex reconstructed");
1688 vertex = new AliESDVertex(vtxPos, vtxErr);
1690 esd->SetPrimaryVertexSPD(vertex);
1691 // if SPD multiplicity has been determined, it is stored in the ESD
1692 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1693 if(mult)esd->SetMultiplicity(mult);
1695 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1696 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1703 //_____________________________________________________________________________
1704 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1706 // run the HLT barrel tracking
1708 AliCodeTimerAuto("")
1711 AliError("Missing runLoader!");
1715 AliInfo("running HLT tracking");
1717 // Get a pointer to the HLT reconstructor
1718 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1719 if (!reconstructor) return kFALSE;
1722 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1723 TString detName = fgkDetectorName[iDet];
1724 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1725 reconstructor->SetOption(detName.Data());
1726 AliTracker *tracker = reconstructor->CreateTracker();
1728 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1729 if (fStopOnError) return kFALSE;
1733 Double_t vtxErr[3]={0.005,0.005,0.010};
1734 const AliESDVertex *vertex = esd->GetVertex();
1735 vertex->GetXYZ(vtxPos);
1736 tracker->SetVertex(vtxPos,vtxErr);
1738 fLoader[iDet]->LoadRecPoints("read");
1739 TTree* tree = fLoader[iDet]->TreeR();
1741 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1744 tracker->LoadClusters(tree);
1746 if (tracker->Clusters2Tracks(esd) != 0) {
1747 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1751 tracker->UnloadClusters();
1759 //_____________________________________________________________________________
1760 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1762 // run the muon spectrometer tracking
1764 AliCodeTimerAuto("")
1767 AliError("Missing runLoader!");
1770 Int_t iDet = 7; // for MUON
1772 AliInfo("is running...");
1774 // Get a pointer to the MUON reconstructor
1775 AliReconstructor *reconstructor = GetReconstructor(iDet);
1776 if (!reconstructor) return kFALSE;
1779 TString detName = fgkDetectorName[iDet];
1780 AliDebug(1, Form("%s tracking", detName.Data()));
1781 AliTracker *tracker = reconstructor->CreateTracker();
1783 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1788 fLoader[iDet]->LoadRecPoints("read");
1790 tracker->LoadClusters(fLoader[iDet]->TreeR());
1792 Int_t rv = tracker->Clusters2Tracks(esd);
1796 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1800 fLoader[iDet]->UnloadRecPoints();
1802 tracker->UnloadClusters();
1810 //_____________________________________________________________________________
1811 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1813 // run the barrel tracking
1814 static Int_t eventNr=0;
1815 AliCodeTimerAuto("")
1817 AliInfo("running tracking");
1819 //Fill the ESD with the T0 info (will be used by the TOF)
1820 if (fReconstructor[11] && fLoader[11]) {
1821 fLoader[11]->LoadRecPoints("READ");
1822 TTree *treeR = fLoader[11]->TreeR();
1823 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1826 // pass 1: TPC + ITS inwards
1827 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1828 if (!fTracker[iDet]) continue;
1829 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1832 fLoader[iDet]->LoadRecPoints("read");
1833 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1834 TTree* tree = fLoader[iDet]->TreeR();
1836 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1839 fTracker[iDet]->LoadClusters(tree);
1840 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1842 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1843 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1846 if (fCheckPointLevel > 1) {
1847 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1849 // preliminary PID in TPC needed by the ITS tracker
1851 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1852 AliESDpid::MakePID(esd);
1854 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1857 // pass 2: ALL backwards
1859 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1860 if (!fTracker[iDet]) continue;
1861 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1864 if (iDet > 1) { // all except ITS, TPC
1866 fLoader[iDet]->LoadRecPoints("read");
1867 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1868 tree = fLoader[iDet]->TreeR();
1870 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1873 fTracker[iDet]->LoadClusters(tree);
1874 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1878 if (iDet>1) // start filling residuals for the "outer" detectors
1879 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1881 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1882 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1885 if (fCheckPointLevel > 1) {
1886 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1890 if (iDet > 2) { // all except ITS, TPC, TRD
1891 fTracker[iDet]->UnloadClusters();
1892 fLoader[iDet]->UnloadRecPoints();
1894 // updated PID in TPC needed by the ITS tracker -MI
1896 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1897 AliESDpid::MakePID(esd);
1899 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1901 //stop filling residuals for the "outer" detectors
1902 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1904 // write space-points to the ESD in case alignment data output
1906 if (fWriteAlignmentData)
1907 WriteAlignmentData(esd);
1909 // pass 3: TRD + TPC + ITS refit inwards
1911 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1912 if (!fTracker[iDet]) continue;
1913 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1916 if (iDet<2) // start filling residuals for TPC and ITS
1917 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1919 if (fTracker[iDet]->RefitInward(esd) != 0) {
1920 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1923 // run postprocessing
1924 if (fTracker[iDet]->PostProcess(esd) != 0) {
1925 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1928 if (fCheckPointLevel > 1) {
1929 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1931 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1933 fTracker[iDet]->UnloadClusters();
1934 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1935 fLoader[iDet]->UnloadRecPoints();
1936 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1938 // stop filling residuals for TPC and ITS
1939 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1945 //_____________________________________________________________________________
1946 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1948 // Remove the data which are not needed for the physics analysis.
1951 Int_t nTracks=esd->GetNumberOfTracks();
1952 Int_t nV0s=esd->GetNumberOfV0s();
1954 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1956 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1957 Bool_t rc=esd->Clean(cleanPars);
1959 nTracks=esd->GetNumberOfTracks();
1960 nV0s=esd->GetNumberOfV0s();
1962 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1967 //_____________________________________________________________________________
1968 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1970 // fill the event summary data
1972 AliCodeTimerAuto("")
1973 static Int_t eventNr=0;
1974 TString detStr = detectors;
1976 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1977 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1978 AliReconstructor* reconstructor = GetReconstructor(iDet);
1979 if (!reconstructor) continue;
1980 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1981 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1982 TTree* clustersTree = NULL;
1983 if (fLoader[iDet]) {
1984 fLoader[iDet]->LoadRecPoints("read");
1985 clustersTree = fLoader[iDet]->TreeR();
1986 if (!clustersTree) {
1987 AliError(Form("Can't get the %s clusters tree",
1988 fgkDetectorName[iDet]));
1989 if (fStopOnError) return kFALSE;
1992 if (fRawReader && !reconstructor->HasDigitConversion()) {
1993 reconstructor->FillESD(fRawReader, clustersTree, esd);
1995 TTree* digitsTree = NULL;
1996 if (fLoader[iDet]) {
1997 fLoader[iDet]->LoadDigits("read");
1998 digitsTree = fLoader[iDet]->TreeD();
2000 AliError(Form("Can't get the %s digits tree",
2001 fgkDetectorName[iDet]));
2002 if (fStopOnError) return kFALSE;
2005 reconstructor->FillESD(digitsTree, clustersTree, esd);
2006 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2008 if (fLoader[iDet]) {
2009 fLoader[iDet]->UnloadRecPoints();
2012 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
2016 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2017 AliError(Form("the following detectors were not found: %s",
2019 if (fStopOnError) return kFALSE;
2021 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2026 //_____________________________________________________________________________
2027 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2029 // Reads the trigger decision which is
2030 // stored in Trigger.root file and fills
2031 // the corresponding esd entries
2033 AliCodeTimerAuto("")
2035 AliInfo("Filling trigger information into the ESD");
2037 AliCentralTrigger *aCTP = NULL;
2040 AliCTPRawStream input(fRawReader);
2041 if (!input.Next()) {
2042 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger mask will be taken from the event header, trigger cluster mask will be empty !");
2043 ULong64_t mask = (((ULong64_t)fRawReader->GetTriggerPattern()[1]) << 32) +
2044 fRawReader->GetTriggerPattern()[0];
2045 esd->SetTriggerMask(mask);
2046 esd->SetTriggerCluster(0);
2049 esd->SetTriggerMask(input.GetClassMask());
2050 esd->SetTriggerCluster(input.GetClusterMask());
2053 aCTP = new AliCentralTrigger();
2054 TString configstr("");
2055 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
2056 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
2062 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
2064 if (!runloader->LoadTrigger()) {
2065 aCTP = runloader->GetTrigger();
2066 esd->SetTriggerMask(aCTP->GetClassMask());
2067 esd->SetTriggerCluster(aCTP->GetClusterMask());
2070 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
2075 AliError("No run loader is available! The trigger information is not stored in the ESD !");
2080 // Now fill the trigger class names into AliESDRun object
2081 AliTriggerConfiguration *config = aCTP->GetConfiguration();
2083 AliError("No trigger configuration has been found! The trigger classes information will not be stored in ESD!");
2084 if (fRawReader) delete aCTP;
2088 const TObjArray& classesArray = config->GetClasses();
2089 Int_t nclasses = classesArray.GetEntriesFast();
2090 for( Int_t j=0; j<nclasses; j++ ) {
2091 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
2092 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
2093 esd->SetTriggerClass(trclass->GetName(),trindex);
2096 if (fRawReader) delete aCTP;
2104 //_____________________________________________________________________________
2105 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2108 // Filling information from RawReader Header
2111 AliInfo("Filling information from RawReader Header");
2112 esd->SetBunchCrossNumber(0);
2113 esd->SetOrbitNumber(0);
2114 esd->SetPeriodNumber(0);
2115 esd->SetTimeStamp(0);
2116 esd->SetEventType(0);
2117 const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
2120 const UInt_t *id = eventHeader->GetP("Id");
2121 esd->SetBunchCrossNumber((id)[1]&0x00000fff);
2122 esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
2123 esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
2125 esd->SetTimeStamp((eventHeader->Get("Timestamp")));
2126 esd->SetEventType((eventHeader->Get("Type")));
2133 //_____________________________________________________________________________
2134 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2136 // check whether detName is contained in detectors
2137 // if yes, it is removed from detectors
2139 // check if all detectors are selected
2140 if ((detectors.CompareTo("ALL") == 0) ||
2141 detectors.BeginsWith("ALL ") ||
2142 detectors.EndsWith(" ALL") ||
2143 detectors.Contains(" ALL ")) {
2148 // search for the given detector
2149 Bool_t result = kFALSE;
2150 if ((detectors.CompareTo(detName) == 0) ||
2151 detectors.BeginsWith(detName+" ") ||
2152 detectors.EndsWith(" "+detName) ||
2153 detectors.Contains(" "+detName+" ")) {
2154 detectors.ReplaceAll(detName, "");
2158 // clean up the detectors string
2159 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2160 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2161 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2166 //_____________________________________________________________________________
2167 Bool_t AliReconstruction::InitRunLoader()
2169 // get or create the run loader
2171 if (gAlice) delete gAlice;
2174 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2175 // load all base libraries to get the loader classes
2176 TString libs = gSystem->GetLibraries();
2177 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2178 TString detName = fgkDetectorName[iDet];
2179 if (detName == "HLT") continue;
2180 if (libs.Contains("lib" + detName + "base.so")) continue;
2181 gSystem->Load("lib" + detName + "base.so");
2183 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2185 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2190 fRunLoader->CdGAFile();
2191 fRunLoader->LoadgAlice();
2193 //PH This is a temporary fix to give access to the kinematics
2194 //PH that is needed for the labels of ITS clusters
2195 fRunLoader->LoadHeader();
2196 fRunLoader->LoadKinematics();
2198 } else { // galice.root does not exist
2200 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2204 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2205 AliConfig::GetDefaultEventFolderName(),
2208 AliError(Form("could not create run loader in file %s",
2209 fGAliceFileName.Data()));
2213 fIsNewRunLoader = kTRUE;
2214 fRunLoader->MakeTree("E");
2216 if (fNumberOfEventsPerFile > 0)
2217 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2219 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2225 //_____________________________________________________________________________
2226 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2228 // get the reconstructor object and the loader for a detector
2230 if (fReconstructor[iDet]) return fReconstructor[iDet];
2232 // load the reconstructor object
2233 TPluginManager* pluginManager = gROOT->GetPluginManager();
2234 TString detName = fgkDetectorName[iDet];
2235 TString recName = "Ali" + detName + "Reconstructor";
2237 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2239 AliReconstructor* reconstructor = NULL;
2240 // first check if a plugin is defined for the reconstructor
2241 TPluginHandler* pluginHandler =
2242 pluginManager->FindHandler("AliReconstructor", detName);
2243 // if not, add a plugin for it
2244 if (!pluginHandler) {
2245 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2246 TString libs = gSystem->GetLibraries();
2247 if (libs.Contains("lib" + detName + "base.so") ||
2248 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2249 pluginManager->AddHandler("AliReconstructor", detName,
2250 recName, detName + "rec", recName + "()");
2252 pluginManager->AddHandler("AliReconstructor", detName,
2253 recName, detName, recName + "()");
2255 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2257 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2258 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2260 if (reconstructor) {
2261 TObject* obj = fOptions.FindObject(detName.Data());
2262 if (obj) reconstructor->SetOption(obj->GetTitle());
2263 reconstructor->Init();
2264 fReconstructor[iDet] = reconstructor;
2267 // get or create the loader
2268 if (detName != "HLT") {
2269 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2270 if (!fLoader[iDet]) {
2271 AliConfig::Instance()
2272 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2274 // first check if a plugin is defined for the loader
2276 pluginManager->FindHandler("AliLoader", detName);
2277 // if not, add a plugin for it
2278 if (!pluginHandler) {
2279 TString loaderName = "Ali" + detName + "Loader";
2280 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2281 pluginManager->AddHandler("AliLoader", detName,
2282 loaderName, detName + "base",
2283 loaderName + "(const char*, TFolder*)");
2284 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2286 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2288 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2289 fRunLoader->GetEventFolder());
2291 if (!fLoader[iDet]) { // use default loader
2292 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2294 if (!fLoader[iDet]) {
2295 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2296 if (fStopOnError) return NULL;
2298 fRunLoader->AddLoader(fLoader[iDet]);
2299 fRunLoader->CdGAFile();
2300 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2301 fRunLoader->Write(0, TObject::kOverwrite);
2306 return reconstructor;
2309 //_____________________________________________________________________________
2310 Bool_t AliReconstruction::CreateVertexer()
2312 // create the vertexer
2315 AliReconstructor* itsReconstructor = GetReconstructor(0);
2316 if (itsReconstructor) {
2317 fVertexer = itsReconstructor->CreateVertexer();
2320 AliWarning("couldn't create a vertexer for ITS");
2321 if (fStopOnError) return kFALSE;
2327 //_____________________________________________________________________________
2328 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2330 // create the trackers
2332 TString detStr = detectors;
2333 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2334 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2335 AliReconstructor* reconstructor = GetReconstructor(iDet);
2336 if (!reconstructor) continue;
2337 TString detName = fgkDetectorName[iDet];
2338 if (detName == "HLT") {
2339 fRunHLTTracking = kTRUE;
2342 if (detName == "MUON") {
2343 fRunMuonTracking = kTRUE;
2348 fTracker[iDet] = reconstructor->CreateTracker();
2349 if (!fTracker[iDet] && (iDet < 7)) {
2350 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2351 if (fStopOnError) return kFALSE;
2353 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2359 //_____________________________________________________________________________
2360 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2362 // delete trackers and the run loader and close and delete the file
2364 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2365 delete fReconstructor[iDet];
2366 fReconstructor[iDet] = NULL;
2367 fLoader[iDet] = NULL;
2368 delete fTracker[iDet];
2369 fTracker[iDet] = NULL;
2370 // delete fQADataMaker[iDet];
2371 // fQADataMaker[iDet] = NULL;
2376 if (ftVertexer) delete ftVertexer;
2379 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2380 delete fDiamondProfile;
2381 fDiamondProfile = NULL;
2382 delete fDiamondProfileTPC;
2383 fDiamondProfileTPC = NULL;
2393 if (fParentRawReader) delete fParentRawReader;
2394 fParentRawReader=NULL;
2404 gSystem->Unlink("AliESDs.old.root");
2409 //_____________________________________________________________________________
2411 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2413 // read the ESD event from a file
2415 if (!esd) return kFALSE;
2417 sprintf(fileName, "ESD_%d.%d_%s.root",
2418 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2419 if (gSystem->AccessPathName(fileName)) return kFALSE;
2421 AliInfo(Form("reading ESD from file %s", fileName));
2422 AliDebug(1, Form("reading ESD from file %s", fileName));
2423 TFile* file = TFile::Open(fileName);
2424 if (!file || !file->IsOpen()) {
2425 AliError(Form("opening %s failed", fileName));
2432 esd = (AliESDEvent*) file->Get("ESD");
2441 //_____________________________________________________________________________
2442 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2444 // write the ESD event to a file
2448 sprintf(fileName, "ESD_%d.%d_%s.root",
2449 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2451 AliDebug(1, Form("writing ESD to file %s", fileName));
2452 TFile* file = TFile::Open(fileName, "recreate");
2453 if (!file || !file->IsOpen()) {
2454 AliError(Form("opening %s failed", fileName));
2463 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2465 // Write space-points which are then used in the alignment procedures
2466 // For the moment only ITS, TRD and TPC
2468 // Load TOF clusters
2470 fLoader[3]->LoadRecPoints("read");
2471 TTree* tree = fLoader[3]->TreeR();
2473 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2476 fTracker[3]->LoadClusters(tree);
2478 Int_t ntracks = esd->GetNumberOfTracks();
2479 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2481 AliESDtrack *track = esd->GetTrack(itrack);
2484 for (Int_t iDet = 3; iDet >= 0; iDet--)
2485 nsp += track->GetNcls(iDet);
2487 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2488 track->SetTrackPointArray(sp);
2490 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2491 AliTracker *tracker = fTracker[iDet];
2492 if (!tracker) continue;
2493 Int_t nspdet = track->GetNcls(iDet);
2494 if (nspdet <= 0) continue;
2495 track->GetClusters(iDet,idx);
2499 while (isp2 < nspdet) {
2501 TString dets = fgkDetectorName[iDet];
2502 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2503 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2504 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2505 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2506 isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2508 isvalid = tracker->GetTrackPoint(idx[isp2],p);
2511 const Int_t kNTPCmax = 159;
2512 if (iDet==1 && isp2>kNTPCmax) break; // to be fixed
2513 if (!isvalid) continue;
2514 sp->AddPoint(isptrack,&p); isptrack++; isp++;
2520 fTracker[3]->UnloadClusters();
2521 fLoader[3]->UnloadRecPoints();
2525 //_____________________________________________________________________________
2526 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2528 // The method reads the raw-data error log
2529 // accumulated within the rawReader.
2530 // It extracts the raw-data errors related to
2531 // the current event and stores them into
2532 // a TClonesArray inside the esd object.
2534 if (!fRawReader) return;
2536 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2538 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2540 if (iEvent != log->GetEventNumber()) continue;
2542 esd->AddRawDataErrorLog(log);
2547 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2548 // Dump a file content into a char in TNamed
2550 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2551 Int_t kBytes = (Int_t)in.tellg();
2552 printf("Size: %d \n",kBytes);
2555 char* memblock = new char [kBytes];
2556 in.seekg (0, ios::beg);
2557 in.read (memblock, kBytes);
2559 TString fData(memblock,kBytes);
2560 fn = new TNamed(pName,fData);
2561 printf("fData Size: %d \n",fData.Sizeof());
2562 printf("pName Size: %d \n",pName.Sizeof());
2563 printf("fn Size: %d \n",fn->Sizeof());
2567 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2573 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2574 // This is not really needed in AliReconstruction at the moment
2575 // but can serve as a template
2577 TList *fList = fTree->GetUserInfo();
2578 TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2579 printf("fn Size: %d \n",fn->Sizeof());
2581 TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2582 const char* cdata = fn->GetTitle();
2583 printf("fTmp Size %d\n",fTmp.Sizeof());
2585 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2586 printf("calculated size %d\n",size);
2587 ofstream out(pName.Data(),ios::out | ios::binary);
2588 out.write(cdata,size);
2593 //_____________________________________________________________________________
2594 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2596 // get the quality assurance data maker object and the loader for a detector
2598 if (fQADataMaker[iDet])
2599 return fQADataMaker[iDet];
2601 AliQADataMakerRec * qadm = NULL;
2602 if (iDet == fgkNDetectors) { //Global QA
2603 qadm = new AliGlobalQADataMaker();
2604 fQADataMaker[iDet] = qadm;
2608 // load the QA data maker object
2609 TPluginManager* pluginManager = gROOT->GetPluginManager();
2610 TString detName = fgkDetectorName[iDet];
2611 TString qadmName = "Ali" + detName + "QADataMakerRec";
2612 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT"))
2615 // first check if a plugin is defined for the quality assurance data maker
2616 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2617 // if not, add a plugin for it
2618 if (!pluginHandler) {
2619 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
2620 TString libs = gSystem->GetLibraries();
2621 if (libs.Contains("lib" + detName + "base.so") ||
2622 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2623 pluginManager->AddHandler("AliQADataMakerRec", detName,
2624 qadmName, detName + "qadm", qadmName + "()");
2626 pluginManager->AddHandler("AliQADataMakerRec", detName,
2627 qadmName, detName, qadmName + "()");
2629 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2631 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2632 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
2635 fQADataMaker[iDet] = qadm;
2640 //_____________________________________________________________________________
2641 Bool_t AliReconstruction::RunQA(AliESDEvent *& esd)
2643 // run the Quality Assurance data producer
2645 AliCodeTimerAuto("")
2646 TString detStr = fQADetectors ;
2647 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2648 if (!IsSelected(fgkDetectorName[iDet], detStr))
2650 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
2653 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2654 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2656 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
2657 qadm->Exec(AliQA::kESDS, esd) ;
2660 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2662 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2663 AliError(Form("the following detectors were not found: %s",
2673 //_____________________________________________________________________________
2674 void AliReconstruction::CheckQA()
2676 // check the QA of SIM for this run and remove the detectors
2677 // with status Fatal
2679 TString newRunLocalReconstruction ;
2680 TString newRunTracking ;
2681 TString newFillESD ;
2683 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2684 TString detName(AliQA::GetDetName(iDet)) ;
2685 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2686 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2687 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2689 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2690 fRunLocalReconstruction.Contains("ALL") ) {
2691 newRunLocalReconstruction += detName ;
2692 newRunLocalReconstruction += " " ;
2694 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2695 fRunTracking.Contains("ALL") ) {
2696 newRunTracking += detName ;
2697 newRunTracking += " " ;
2699 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2700 fFillESD.Contains("ALL") ) {
2701 newFillESD += detName ;
2706 fRunLocalReconstruction = newRunLocalReconstruction ;
2707 fRunTracking = newRunTracking ;
2708 fFillESD = newFillESD ;
2711 //_____________________________________________________________________________
2712 Int_t AliReconstruction::GetDetIndex(const char* detector)
2714 // return the detector index corresponding to detector
2716 for (index = 0; index < fgkNDetectors ; index++) {
2717 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2722 //_____________________________________________________________________________
2723 Bool_t AliReconstruction::FinishPlaneEff() {
2725 // Here execute all the necessary operationis, at the end of the tracking phase,
2726 // in case that evaluation of PlaneEfficiencies was required for some detector.
2727 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2729 // This Preliminary version works only FOR ITS !!!!!
2730 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2733 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2736 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2737 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2738 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2739 if(fTracker[iDet]) {
2740 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2741 ret=planeeff->WriteIntoCDB();
2742 if(planeeff->GetCreateHistos()) {
2743 TString name="PlaneEffHisto";
2744 name+=fgkDetectorName[iDet];
2746 ret*=planeeff->WriteHistosToFile(name,"RECREATE");
2752 //_____________________________________________________________________________
2753 Bool_t AliReconstruction::InitPlaneEff() {
2755 // Here execute all the necessary operations, before of the tracking phase,
2756 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2757 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2758 // which should be updated/recalculated.
2760 // This Preliminary version will work only FOR ITS !!!!!
2761 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2764 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2766 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2770 //_____________________________________________________________________________
2771 Bool_t AliReconstruction::InitAliEVE()
2773 // This method should be called only in case
2774 // AliReconstruction is run
2775 // within the alieve environment.
2776 // It will initialize AliEVE in a way
2777 // so that it can visualize event processed
2778 // by AliReconstruction.
2779 // The return flag shows whenever the
2780 // AliEVE initialization was successful or not.
2783 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2784 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2785 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2787 gROOT->ProcessLine("if (!gAliEveEvent) {gAliEveEvent = new AliEveEventManager();gAliEveEvent->SetAutoLoad(kTRUE);gAliEveEvent->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(gAliEveEvent);};");
2788 gROOT->ProcessLine("alieve_online_init()");
2793 //_____________________________________________________________________________
2794 void AliReconstruction::RunAliEVE()
2796 // Runs AliEVE visualisation of
2797 // the current event.
2798 // Should be executed only after
2799 // successful initialization of AliEVE.
2801 AliInfo("Running AliEVE...");
2802 gROOT->ProcessLine(Form("gAliEveEvent->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2803 gROOT->ProcessLine("gAliEveEvent->StartStopAutoLoadTimer();");
2807 //_____________________________________________________________________________
2808 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2810 // Allows to run QA for a selected set of detectors
2811 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2812 // all selected detectors run the same selected tasks
2814 if (!detAndAction.Contains(":")) {
2815 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2819 Int_t colon = detAndAction.Index(":") ;
2820 fQADetectors = detAndAction(0, colon) ;
2821 if (fQADetectors.Contains("ALL") )
2822 fQADetectors = fFillESD ;
2823 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2824 if (fQATasks.Contains("ALL") ) {
2825 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2827 fQATasks.ToUpper() ;
2829 if ( fQATasks.Contains("RAW") )
2830 tempo = Form("%d ", AliQA::kRAWS) ;
2831 if ( fQATasks.Contains("RECPOINT") )
2832 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2833 if ( fQATasks.Contains("ESD") )
2834 tempo += Form("%d ", AliQA::kESDS) ;
2836 if (fQATasks.IsNull()) {
2837 AliInfo("No QA requested\n") ;
2842 TString tempo(fQATasks) ;
2843 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2844 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2845 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2846 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;