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 l3Current = TMath::Abs(l3Current);
673 if (TMath::Abs(l3Current-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
674 map=AliMagWrapCheb::k5kG;
677 if (TMath::Abs(l3Current-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
678 map=AliMagWrapCheb::k2kG;
681 if (l3Current < zero) {
682 map=AliMagWrapCheb::k2kG;
684 factor=0.; // in fact, this is a global factor...
685 fUniformField=kTRUE; // track with the uniform (zero) B field
687 AliError(Form("Wrong L3 current (%f A)!",l3Current));
691 diCurrent = TMath::Abs(diCurrent);
692 if (TMath::Abs(diCurrent-diNominalCurrent)/diNominalCurrent < tolerance) {
693 // 3% current tolerance...
697 if (diCurrent < zero) { // some small current..
701 AliError(Form("Wrong dipole current (%f A)!",diCurrent));
705 delete fForcedFieldMap;
707 new AliMagWrapCheb("B field map ",s,2,factor,10.,map,dipoleON,path);
709 fForcedFieldMap->Print();
711 AliTracker::SetFieldMap(fForcedFieldMap,fUniformField);
717 Bool_t AliReconstruction::InitGRP() {
718 //------------------------------------
719 // Initialization of the GRP entry
720 //------------------------------------
721 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
723 if (entry) fGRPData = dynamic_cast<TMap*>(entry->GetObject());
726 AliError("No GRP entry found in OCDB!");
731 //*** Dealing with the magnetic field map
732 if (AliTracker::GetFieldMap()) {
733 AliInfo("Running with the externally set B field !");
735 // Construct the field map out of the information retrieved from GRP.
740 TObjString *l3Current=
741 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
743 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
746 TObjString *l3Polarity=
747 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
749 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
754 TObjString *diCurrent=
755 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
757 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
760 TObjString *diPolarity=
761 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
763 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
768 Float_t l3Cur=TMath::Abs(atof(l3Current->GetName()));
769 Float_t diCur=TMath::Abs(atof(diCurrent->GetName()));
770 Float_t l3Pol=atof(l3Polarity->GetName());
772 if (l3Pol != 0.) factor=-1.;
775 if (!SetFieldMap(l3Cur, diCur, factor)) {
776 AliFatal("Failed to creat a B field map ! Exiting...");
778 AliInfo("Running with the B field constructed out of GRP !");
781 AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
787 //*** Get the diamond profile from OCDB
788 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
790 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
792 AliError("No diamond profile found in OCDB!");
795 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
797 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
799 AliError("No diamond profile found in OCDB!");
805 //_____________________________________________________________________________
806 Bool_t AliReconstruction::Run(const char* input)
809 AliCodeTimerAuto("");
811 if (!InitRun(input)) return kFALSE;
812 //******* The loop over events
814 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
815 (fRawReader && fRawReader->NextEvent())) {
816 if (!RunEvent(iEvent)) return kFALSE;
820 if (!FinishRun()) return kFALSE;
825 //_____________________________________________________________________________
826 Bool_t AliReconstruction::InitRun(const char* input)
828 // Initialize all the stuff before
829 // going into the event loop
830 // If the second argument is given, the first one is ignored and
831 // the reconstruction works in an online mode
832 AliCodeTimerAuto("");
834 // Overwrite the previous setting
835 if (input) fInput = input;
837 // set the input in case of raw data
838 fRawReader = AliRawReader::Create(fInput.Data());
840 AliInfo("Reconstruction will run over digits");
842 if (!fEquipIdMap.IsNull() && fRawReader)
843 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
845 if (!fUseHLTData.IsNull()) {
846 // create the RawReaderHLT which performs redirection of HLT input data for
847 // the specified detectors
848 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
850 fParentRawReader=fRawReader;
851 fRawReader=pRawReader;
853 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
857 AliSysInfo::AddStamp("Start");
858 // get the run loader
859 if (!InitRunLoader()) return kFALSE;
860 AliSysInfo::AddStamp("LoadLoader");
862 // Initialize the CDB storage
865 AliSysInfo::AddStamp("LoadCDB");
867 // Set run number in CDBManager (if it is not already set by the user)
868 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
870 // Set CDB lock: from now on it is forbidden to reset the run number
871 // or the default storage or to activate any further storage!
874 // Import ideal TGeo geometry and apply misalignment
876 TString geom(gSystem->DirName(fGAliceFileName));
877 geom += "/geometry.root";
878 AliGeomManager::LoadGeometry(geom.Data());
880 TString detsToCheck=fRunLocalReconstruction;
881 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data()))
882 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
883 if (!gGeoManager) if (fStopOnError) return kFALSE;
886 if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
887 AliSysInfo::AddStamp("LoadGeom");
890 if (!InitGRP()) return kFALSE;
893 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
894 if(fDiamondProfile && fMeanVertexConstraint) ftVertexer->SetVtxStart(fDiamondProfile);
897 if (fRunVertexFinder && !CreateVertexer()) {
903 AliSysInfo::AddStamp("Vertexer");
906 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
912 AliSysInfo::AddStamp("LoadTrackers");
914 // get the possibly already existing ESD file and tree
915 fesd = new AliESDEvent(); fhltesd = new AliESDEvent();
916 if (!gSystem->AccessPathName("AliESDs.root")){
917 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
918 ffileOld = TFile::Open("AliESDs.old.root");
919 if (ffileOld && ffileOld->IsOpen()) {
920 ftreeOld = (TTree*) ffileOld->Get("esdTree");
921 if (ftreeOld)fesd->ReadFromTree(ftreeOld);
922 fhlttreeOld = (TTree*) ffileOld->Get("HLTesdTree");
923 if (fhlttreeOld) fhltesd->ReadFromTree(fhlttreeOld);
927 // create the ESD output file and tree
928 ffile = TFile::Open("AliESDs.root", "RECREATE");
929 ffile->SetCompressionLevel(2);
930 if (!ffile->IsOpen()) {
931 AliError("opening AliESDs.root failed");
932 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
935 ftree = new TTree("esdTree", "Tree with ESD objects");
936 fesd = new AliESDEvent();
937 fesd->CreateStdContent();
938 fesd->WriteToTree(ftree);
940 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
941 fhltesd = new AliESDEvent();
942 fhltesd->CreateStdContent();
943 fhltesd->WriteToTree(fhlttree);
946 if (fWriteESDfriend) {
947 fesdf = new AliESDfriend();
948 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
949 br->SetFile("AliESDfriends.root");
950 fesd->AddObject(fesdf);
955 if (fRawReader) fRawReader->RewindEvents();
958 gSystem->GetProcInfo(&ProcInfo);
959 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
962 if (fRunQA && fRawReader && fQATasks.Contains(Form("%d", AliQA::kRAWS))) {
963 AliQADataMakerSteer qas ;
964 qas.Run(fQADetectors, fRawReader) ;
965 fSameQACycle = kTRUE ;
967 //Initialize the QA and start of cycle for out-of-cycle QA
969 TString detStr(fQADetectors) ;
970 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
971 if (!IsSelected(fgkDetectorName[iDet], detStr))
973 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
976 AliInfo(Form("Initializing the QA data maker for %s",
977 fgkDetectorName[iDet]));
978 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
979 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
980 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
981 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
983 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
984 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
985 fSameQACycle = kTRUE;
987 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
988 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle);
989 fSameQACycle = kTRUE;
994 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
995 AliInfo(Form("Initializing the global QA data maker"));
996 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
998 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
999 AliTracker::SetResidualsArray(arr);
1001 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1002 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
1005 fSameQACycle = kFALSE;
1006 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1007 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
1008 fSameQACycle = kTRUE;
1010 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1011 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle);
1012 fSameQACycle = kTRUE;
1018 //Initialize the Plane Efficiency framework
1019 if (fRunPlaneEff && !InitPlaneEff()) {
1020 if(fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1023 if (strcmp(gProgName,"alieve") == 0)
1024 fRunAliEVE = InitAliEVE();
1029 //_____________________________________________________________________________
1030 Bool_t AliReconstruction::RunEvent(Int_t iEvent)
1032 // run the reconstruction over a single event
1033 // The event loop is steered in Run method
1035 AliCodeTimerAuto("");
1037 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1038 fRunLoader->SetEventNumber(iEvent);
1039 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1041 //?? fRunLoader->MakeTree("H");
1042 fRunLoader->TreeE()->Fill();
1045 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1046 // copy old ESD to the new one
1048 fesd->ReadFromTree(ftreeOld);
1049 ftreeOld->GetEntry(iEvent);
1053 fhltesd->ReadFromTree(fhlttreeOld);
1054 fhlttreeOld->GetEntry(iEvent);
1060 AliInfo(Form("processing event %d", iEvent));
1062 //Start of cycle for the in-loop QA
1065 fSameQACycle = kFALSE ;
1066 TString detStr(fQADetectors);
1067 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1068 if (!IsSelected(fgkDetectorName[iDet], detStr))
1070 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
1073 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1074 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
1075 fSameQACycle = kTRUE;
1077 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1078 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle) ;
1079 fSameQACycle = kTRUE;
1083 fSameQACycle = kFALSE;
1084 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1085 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1086 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
1087 fSameQACycle = kTRUE;
1089 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1090 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle);
1091 fSameQACycle = kTRUE;
1097 fRunLoader->GetEvent(iEvent);
1099 char aFileName[256];
1100 sprintf(aFileName, "ESD_%d.%d_final.root",
1101 fRunLoader->GetHeader()->GetRun(),
1102 fRunLoader->GetHeader()->GetEventNrInRun());
1103 if (!gSystem->AccessPathName(aFileName)) return kTRUE;
1105 // local single event reconstruction
1106 if (!fRunLocalReconstruction.IsNull()) {
1107 TString detectors=fRunLocalReconstruction;
1108 // run HLT event reconstruction first
1109 // ;-( IsSelected changes the string
1110 if (IsSelected("HLT", detectors) &&
1111 !RunLocalEventReconstruction("HLT")) {
1112 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1114 detectors=fRunLocalReconstruction;
1115 detectors.ReplaceAll("HLT", "");
1116 if (!RunLocalEventReconstruction(detectors)) {
1117 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1121 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1122 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1123 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1124 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1126 // Set magnetic field from the tracker
1127 fesd->SetMagneticField(AliTracker::GetBz());
1128 fhltesd->SetMagneticField(AliTracker::GetBz());
1132 // Fill raw-data error log into the ESD
1133 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1136 if (fRunVertexFinder) {
1137 if (!ReadESD(fesd, "vertex")) {
1138 if (!RunVertexFinder(fesd)) {
1139 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1141 if (fCheckPointLevel > 0) WriteESD(fesd, "vertex");
1146 if (!fRunTracking.IsNull()) {
1147 if (fRunMuonTracking) {
1148 if (!RunMuonTracking(fesd)) {
1149 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1155 if (!fRunTracking.IsNull()) {
1156 if (!ReadESD(fesd, "tracking")) {
1157 if (!RunTracking(fesd)) {
1158 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1160 if (fCheckPointLevel > 0) WriteESD(fesd, "tracking");
1165 if (!fFillESD.IsNull()) {
1166 TString detectors=fFillESD;
1167 // run HLT first and on hltesd
1168 // ;-( IsSelected changes the string
1169 if (IsSelected("HLT", detectors) &&
1170 !FillESD(fhltesd, "HLT")) {
1171 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1174 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1175 if (detectors.Contains("ALL")) {
1177 for (Int_t idet=0; idet<fgkNDetectors; ++idet){
1178 detectors += fgkDetectorName[idet];
1182 detectors.ReplaceAll("HLT", "");
1183 if (!FillESD(fesd, detectors)) {
1184 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1188 // fill Event header information from the RawEventHeader
1189 if (fRawReader){FillRawEventHeaderESD(fesd);}
1192 AliESDpid::MakePID(fesd);
1193 if (fCheckPointLevel > 1) WriteESD(fesd, "PID");
1195 if (fFillTriggerESD) {
1196 if (!ReadESD(fesd, "trigger")) {
1197 if (!FillTriggerESD(fesd)) {
1198 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1200 if (fCheckPointLevel > 1) WriteESD(fesd, "trigger");
1207 // Propagate track to the beam pipe (if not already done by ITS)
1209 const Int_t ntracks = fesd->GetNumberOfTracks();
1210 const Double_t kBz = fesd->GetMagneticField();
1211 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1214 UShort_t *selectedIdx=new UShort_t[ntracks];
1216 for (Int_t itrack=0; itrack<ntracks; itrack++){
1217 const Double_t kMaxStep = 5; //max step over the material
1220 AliESDtrack *track = fesd->GetTrack(itrack);
1221 if (!track) continue;
1223 AliExternalTrackParam *tpcTrack =
1224 (AliExternalTrackParam *)track->GetTPCInnerParam();
1228 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
1231 Int_t n=trkArray.GetEntriesFast();
1232 selectedIdx[n]=track->GetID();
1233 trkArray.AddLast(tpcTrack);
1236 //Tracks refitted by ITS should already be at the SPD vertex
1237 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1240 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1241 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1246 // Improve the reconstructed primary vertex position using the tracks
1248 TObject *obj = fOptions.FindObject("ITS");
1250 TString optITS = obj->GetTitle();
1251 if (optITS.Contains("cosmics") || optITS.Contains("COSMICS"))
1252 fRunVertexFinderTracks=kFALSE;
1254 if (fRunVertexFinderTracks) {
1255 // TPC + ITS primary vertex
1256 ftVertexer->SetITSrefitRequired();
1257 if(fDiamondProfile && fMeanVertexConstraint) {
1258 ftVertexer->SetVtxStart(fDiamondProfile);
1260 ftVertexer->SetConstraintOff();
1262 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1264 if (pvtx->GetStatus()) {
1265 fesd->SetPrimaryVertex(pvtx);
1266 for (Int_t i=0; i<ntracks; i++) {
1267 AliESDtrack *t = fesd->GetTrack(i);
1268 t->RelateToVertex(pvtx, kBz, kVeryBig);
1273 // TPC-only primary vertex
1274 ftVertexer->SetITSrefitNotRequired();
1275 if(fDiamondProfileTPC && fMeanVertexConstraint) {
1276 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1278 ftVertexer->SetConstraintOff();
1280 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1282 if (pvtx->GetStatus()) {
1283 fesd->SetPrimaryVertexTPC(pvtx);
1284 for (Int_t i=0; i<ntracks; i++) {
1285 AliESDtrack *t = fesd->GetTrack(i);
1286 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1292 delete[] selectedIdx;
1294 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1299 AliV0vertexer vtxer;
1300 vtxer.Tracks2V0vertices(fesd);
1302 if (fRunCascadeFinder) {
1304 AliCascadeVertexer cvtxer;
1305 cvtxer.V0sTracks2CascadeVertices(fesd);
1310 if (fCleanESD) CleanESD(fesd);
1314 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1315 if (qadm && fQATasks.Contains(Form("%d", AliQA::kESDS)))
1316 qadm->Exec(AliQA::kESDS, fesd);
1320 if (fWriteESDfriend) {
1321 fesdf->~AliESDfriend();
1322 new (fesdf) AliESDfriend(); // Reset...
1323 fesd->GetESDfriend(fesdf);
1327 // Auto-save the ESD tree in case of prompt reco @P2
1328 if (fRawReader && fRawReader->UseAutoSaveESD())
1329 ftree->AutoSave("SaveSelf");
1335 if (fRunAliEVE) RunAliEVE();
1337 if (fCheckPointLevel > 0) WriteESD(fesd, "final");
1340 if (fWriteESDfriend) {
1341 fesdf->~AliESDfriend();
1342 new (fesdf) AliESDfriend(); // Reset...
1345 ProcInfo_t ProcInfo;
1346 gSystem->GetProcInfo(&ProcInfo);
1347 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1350 // End of cycle for the in-loop
1354 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1355 if (!IsSelected(fgkDetectorName[iDet], fQADetectors))
1357 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1360 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1361 qadm->EndOfCycle(AliQA::kRECPOINTS);
1362 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1363 qadm->EndOfCycle(AliQA::kESDS);
1368 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1370 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1371 qadm->EndOfCycle(AliQA::kRECPOINTS);
1372 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1373 qadm->EndOfCycle(AliQA::kESDS);
1382 //_____________________________________________________________________________
1383 Bool_t AliReconstruction::FinishRun()
1386 // Called after the exit
1387 // from the event loop
1388 AliCodeTimerAuto("");
1390 if (fIsNewRunLoader) { // galice.root didn't exist
1391 fRunLoader->WriteHeader("OVERWRITE");
1392 fRunLoader->CdGAFile();
1393 fRunLoader->Write(0, TObject::kOverwrite);
1396 ftree->GetUserInfo()->Add(fesd);
1397 fhlttree->GetUserInfo()->Add(fhltesd);
1399 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1400 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1402 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1403 cdbMapCopy->SetOwner(1);
1404 cdbMapCopy->SetName("cdbMap");
1405 TIter iter(cdbMap->GetTable());
1408 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1409 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1410 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1411 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1414 TList *cdbListCopy = new TList();
1415 cdbListCopy->SetOwner(1);
1416 cdbListCopy->SetName("cdbList");
1418 TIter iter2(cdbList);
1421 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1422 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1425 ftree->GetUserInfo()->Add(cdbMapCopy);
1426 ftree->GetUserInfo()->Add(cdbListCopy);
1429 if(fESDPar.Contains("ESD.par")){
1430 AliInfo("Attaching ESD.par to Tree");
1431 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1432 ftree->GetUserInfo()->Add(fn);
1438 if (fWriteESDfriend)
1439 ftree->SetBranchStatus("ESDfriend*",0);
1440 // we want to have only one tree version number
1441 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1444 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1445 if (fRunPlaneEff && !FinishPlaneEff()) {
1446 AliWarning("Finish PlaneEff evaluation failed");
1450 CleanUp(ffile, ffileOld);
1453 AliWarning("AOD creation not supported anymore during reconstruction. See ANALYSIS/AliAnalysisTaskESDfilter.cxx instead.");
1456 // Create tags for the events in the ESD tree (the ESD tree is always present)
1457 // In case of empty events the tags will contain dummy values
1458 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1459 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1461 AliWarning("AOD tag creation not supported anymore during reconstruction.");
1464 //Finish QA and end of cycle for out-of-loop QA
1467 AliQADataMakerSteer qas;
1468 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1469 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS, fSameQACycle);
1471 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1472 qas.Run(fRunLocalReconstruction.Data(), AliQA::kESDS, fSameQACycle);
1474 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1476 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1477 qadm->EndOfCycle(AliQA::kRECPOINTS);
1478 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1479 qadm->EndOfCycle(AliQA::kESDS);
1486 // Cleanup of CDB manager: cache and active storages!
1487 AliCDBManager::Instance()->ClearCache();
1493 //_____________________________________________________________________________
1494 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1496 // run the local reconstruction
1497 static Int_t eventNr=0;
1498 AliCodeTimerAuto("")
1500 // AliCDBManager* man = AliCDBManager::Instance();
1501 // Bool_t origCache = man->GetCacheFlag();
1503 // TString detStr = detectors;
1504 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1505 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1506 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1507 // if (!reconstructor) continue;
1508 // if (reconstructor->HasLocalReconstruction()) continue;
1510 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1511 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1513 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1514 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1516 // man->SetCacheFlag(kTRUE);
1517 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1518 // man->GetAll(calibPath); // entries are cached!
1520 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1522 // if (fRawReader) {
1523 // fRawReader->RewindEvents();
1524 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1526 // reconstructor->Reconstruct(fRunLoader);
1529 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1530 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1532 // // unload calibration data
1533 // man->UnloadFromCache(calibPath);
1534 // //man->ClearCache();
1537 // man->SetCacheFlag(origCache);
1539 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1540 // AliError(Form("the following detectors were not found: %s",
1542 // if (fStopOnError) return kFALSE;
1549 //_____________________________________________________________________________
1550 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1552 // run the local reconstruction
1554 static Int_t eventNr=0;
1555 AliCodeTimerAuto("")
1557 TString detStr = detectors;
1558 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1559 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1560 AliReconstructor* reconstructor = GetReconstructor(iDet);
1561 if (!reconstructor) continue;
1562 AliLoader* loader = fLoader[iDet];
1563 // Matthias April 2008: temporary fix to run HLT reconstruction
1564 // although the HLT loader is missing
1565 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1567 reconstructor->Reconstruct(fRawReader, NULL);
1570 reconstructor->Reconstruct(dummy, NULL);
1575 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1578 // conversion of digits
1579 if (fRawReader && reconstructor->HasDigitConversion()) {
1580 AliInfo(Form("converting raw data digits into root objects for %s",
1581 fgkDetectorName[iDet]));
1582 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1583 fgkDetectorName[iDet]));
1584 loader->LoadDigits("update");
1585 loader->CleanDigits();
1586 loader->MakeDigitsContainer();
1587 TTree* digitsTree = loader->TreeD();
1588 reconstructor->ConvertDigits(fRawReader, digitsTree);
1589 loader->WriteDigits("OVERWRITE");
1590 loader->UnloadDigits();
1592 // local reconstruction
1593 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1594 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1595 loader->LoadRecPoints("update");
1596 loader->CleanRecPoints();
1597 loader->MakeRecPointsContainer();
1598 TTree* clustersTree = loader->TreeR();
1599 if (fRawReader && !reconstructor->HasDigitConversion()) {
1600 reconstructor->Reconstruct(fRawReader, clustersTree);
1602 loader->LoadDigits("read");
1603 TTree* digitsTree = loader->TreeD();
1605 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1606 if (fStopOnError) return kFALSE;
1608 reconstructor->Reconstruct(digitsTree, clustersTree);
1610 loader->UnloadDigits();
1613 // In-loop QA for local reconstrucion
1614 TString detQAStr(fQADetectors) ;
1615 if (fRunQA && fInLoopQA && IsSelected(fgkDetectorName[iDet], detQAStr) ) {
1616 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1619 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1621 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1623 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1624 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1626 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1630 loader->WriteRecPoints("OVERWRITE");
1631 loader->UnloadRecPoints();
1632 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1635 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1636 AliError(Form("the following detectors were not found: %s",
1638 if (fStopOnError) return kFALSE;
1644 //_____________________________________________________________________________
1645 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1647 // run the barrel tracking
1649 AliCodeTimerAuto("")
1651 AliESDVertex* vertex = NULL;
1652 Double_t vtxPos[3] = {0, 0, 0};
1653 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1654 TArrayF mcVertex(3);
1655 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1656 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1657 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1661 AliInfo("running the ITS vertex finder");
1663 fLoader[0]->LoadRecPoints();
1664 TTree* cltree = fLoader[0]->TreeR();
1666 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1667 vertex = fVertexer->FindVertexForCurrentEvent(cltree);
1670 AliError("Can't get the ITS cluster tree");
1672 fLoader[0]->UnloadRecPoints();
1675 AliError("Can't get the ITS loader");
1678 AliWarning("Vertex not found");
1679 vertex = new AliESDVertex();
1680 vertex->SetName("default");
1683 vertex->SetName("reconstructed");
1687 AliInfo("getting the primary vertex from MC");
1688 vertex = new AliESDVertex(vtxPos, vtxErr);
1692 vertex->GetXYZ(vtxPos);
1693 vertex->GetSigmaXYZ(vtxErr);
1695 AliWarning("no vertex reconstructed");
1696 vertex = new AliESDVertex(vtxPos, vtxErr);
1698 esd->SetPrimaryVertexSPD(vertex);
1699 // if SPD multiplicity has been determined, it is stored in the ESD
1700 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1701 if(mult)esd->SetMultiplicity(mult);
1703 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1704 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1711 //_____________________________________________________________________________
1712 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1714 // run the HLT barrel tracking
1716 AliCodeTimerAuto("")
1719 AliError("Missing runLoader!");
1723 AliInfo("running HLT tracking");
1725 // Get a pointer to the HLT reconstructor
1726 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1727 if (!reconstructor) return kFALSE;
1730 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1731 TString detName = fgkDetectorName[iDet];
1732 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1733 reconstructor->SetOption(detName.Data());
1734 AliTracker *tracker = reconstructor->CreateTracker();
1736 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1737 if (fStopOnError) return kFALSE;
1741 Double_t vtxErr[3]={0.005,0.005,0.010};
1742 const AliESDVertex *vertex = esd->GetVertex();
1743 vertex->GetXYZ(vtxPos);
1744 tracker->SetVertex(vtxPos,vtxErr);
1746 fLoader[iDet]->LoadRecPoints("read");
1747 TTree* tree = fLoader[iDet]->TreeR();
1749 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1752 tracker->LoadClusters(tree);
1754 if (tracker->Clusters2Tracks(esd) != 0) {
1755 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1759 tracker->UnloadClusters();
1767 //_____________________________________________________________________________
1768 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1770 // run the muon spectrometer tracking
1772 AliCodeTimerAuto("")
1775 AliError("Missing runLoader!");
1778 Int_t iDet = 7; // for MUON
1780 AliInfo("is running...");
1782 // Get a pointer to the MUON reconstructor
1783 AliReconstructor *reconstructor = GetReconstructor(iDet);
1784 if (!reconstructor) return kFALSE;
1787 TString detName = fgkDetectorName[iDet];
1788 AliDebug(1, Form("%s tracking", detName.Data()));
1789 AliTracker *tracker = reconstructor->CreateTracker();
1791 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1796 fLoader[iDet]->LoadRecPoints("read");
1798 tracker->LoadClusters(fLoader[iDet]->TreeR());
1800 Int_t rv = tracker->Clusters2Tracks(esd);
1804 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1808 fLoader[iDet]->UnloadRecPoints();
1810 tracker->UnloadClusters();
1818 //_____________________________________________________________________________
1819 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1821 // run the barrel tracking
1822 static Int_t eventNr=0;
1823 AliCodeTimerAuto("")
1825 AliInfo("running tracking");
1827 //Fill the ESD with the T0 info (will be used by the TOF)
1828 if (fReconstructor[11] && fLoader[11]) {
1829 fLoader[11]->LoadRecPoints("READ");
1830 TTree *treeR = fLoader[11]->TreeR();
1831 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1834 // pass 1: TPC + ITS inwards
1835 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1836 if (!fTracker[iDet]) continue;
1837 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1840 fLoader[iDet]->LoadRecPoints("read");
1841 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1842 TTree* tree = fLoader[iDet]->TreeR();
1844 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1847 fTracker[iDet]->LoadClusters(tree);
1848 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1850 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1851 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1854 if (fCheckPointLevel > 1) {
1855 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1857 // preliminary PID in TPC needed by the ITS tracker
1859 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1860 AliESDpid::MakePID(esd);
1862 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1865 // pass 2: ALL backwards
1867 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1868 if (!fTracker[iDet]) continue;
1869 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1872 if (iDet > 1) { // all except ITS, TPC
1874 fLoader[iDet]->LoadRecPoints("read");
1875 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1876 tree = fLoader[iDet]->TreeR();
1878 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1881 fTracker[iDet]->LoadClusters(tree);
1882 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1886 if (iDet>1) // start filling residuals for the "outer" detectors
1887 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1889 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1890 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1893 if (fCheckPointLevel > 1) {
1894 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1898 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
1899 fTracker[iDet]->UnloadClusters();
1900 fLoader[iDet]->UnloadRecPoints();
1902 // updated PID in TPC needed by the ITS tracker -MI
1904 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1905 AliESDpid::MakePID(esd);
1907 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1909 //stop filling residuals for the "outer" detectors
1910 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1912 // pass 3: TRD + TPC + ITS refit inwards
1914 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1915 if (!fTracker[iDet]) continue;
1916 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1919 if (iDet<2) // start filling residuals for TPC and ITS
1920 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1922 if (fTracker[iDet]->RefitInward(esd) != 0) {
1923 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1926 // run postprocessing
1927 if (fTracker[iDet]->PostProcess(esd) != 0) {
1928 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1931 if (fCheckPointLevel > 1) {
1932 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1934 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1937 // write space-points to the ESD in case alignment data output
1939 if (fWriteAlignmentData)
1940 WriteAlignmentData(esd);
1942 for (Int_t iDet = 3; iDet >= 0; iDet--) {
1943 if (!fTracker[iDet]) continue;
1945 fTracker[iDet]->UnloadClusters();
1946 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1947 fLoader[iDet]->UnloadRecPoints();
1948 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1950 // stop filling residuals for TPC and ITS
1951 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1957 //_____________________________________________________________________________
1958 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1960 // Remove the data which are not needed for the physics analysis.
1963 Int_t nTracks=esd->GetNumberOfTracks();
1964 Int_t nV0s=esd->GetNumberOfV0s();
1966 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1968 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1969 Bool_t rc=esd->Clean(cleanPars);
1971 nTracks=esd->GetNumberOfTracks();
1972 nV0s=esd->GetNumberOfV0s();
1974 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1979 //_____________________________________________________________________________
1980 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1982 // fill the event summary data
1984 AliCodeTimerAuto("")
1985 static Int_t eventNr=0;
1986 TString detStr = detectors;
1988 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1989 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1990 AliReconstructor* reconstructor = GetReconstructor(iDet);
1991 if (!reconstructor) continue;
1992 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1993 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1994 TTree* clustersTree = NULL;
1995 if (fLoader[iDet]) {
1996 fLoader[iDet]->LoadRecPoints("read");
1997 clustersTree = fLoader[iDet]->TreeR();
1998 if (!clustersTree) {
1999 AliError(Form("Can't get the %s clusters tree",
2000 fgkDetectorName[iDet]));
2001 if (fStopOnError) return kFALSE;
2004 if (fRawReader && !reconstructor->HasDigitConversion()) {
2005 reconstructor->FillESD(fRawReader, clustersTree, esd);
2007 TTree* digitsTree = NULL;
2008 if (fLoader[iDet]) {
2009 fLoader[iDet]->LoadDigits("read");
2010 digitsTree = fLoader[iDet]->TreeD();
2012 AliError(Form("Can't get the %s digits tree",
2013 fgkDetectorName[iDet]));
2014 if (fStopOnError) return kFALSE;
2017 reconstructor->FillESD(digitsTree, clustersTree, esd);
2018 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2020 if (fLoader[iDet]) {
2021 fLoader[iDet]->UnloadRecPoints();
2024 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
2028 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2029 AliError(Form("the following detectors were not found: %s",
2031 if (fStopOnError) return kFALSE;
2033 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2038 //_____________________________________________________________________________
2039 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2041 // Reads the trigger decision which is
2042 // stored in Trigger.root file and fills
2043 // the corresponding esd entries
2045 AliCodeTimerAuto("")
2047 AliInfo("Filling trigger information into the ESD");
2049 AliCentralTrigger *aCTP = NULL;
2052 AliCTPRawStream input(fRawReader);
2053 if (!input.Next()) {
2054 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 !");
2055 ULong64_t mask = (((ULong64_t)fRawReader->GetTriggerPattern()[1]) << 32) +
2056 fRawReader->GetTriggerPattern()[0];
2057 esd->SetTriggerMask(mask);
2058 esd->SetTriggerCluster(0);
2061 esd->SetTriggerMask(input.GetClassMask());
2062 esd->SetTriggerCluster(input.GetClusterMask());
2065 aCTP = new AliCentralTrigger();
2066 TString configstr("");
2067 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
2068 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
2074 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
2076 if (!runloader->LoadTrigger()) {
2077 aCTP = runloader->GetTrigger();
2078 esd->SetTriggerMask(aCTP->GetClassMask());
2079 esd->SetTriggerCluster(aCTP->GetClusterMask());
2082 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
2087 AliError("No run loader is available! The trigger information is not stored in the ESD !");
2092 // Now fill the trigger class names into AliESDRun object
2093 AliTriggerConfiguration *config = aCTP->GetConfiguration();
2095 AliError("No trigger configuration has been found! The trigger classes information will not be stored in ESD!");
2096 if (fRawReader) delete aCTP;
2100 const TObjArray& classesArray = config->GetClasses();
2101 Int_t nclasses = classesArray.GetEntriesFast();
2102 for( Int_t j=0; j<nclasses; j++ ) {
2103 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
2104 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
2105 esd->SetTriggerClass(trclass->GetName(),trindex);
2108 if (fRawReader) delete aCTP;
2116 //_____________________________________________________________________________
2117 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2120 // Filling information from RawReader Header
2123 if (!fRawReader) return kFALSE;
2125 AliInfo("Filling information from RawReader Header");
2127 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2128 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2129 esd->SetPeriodNumber(fRawReader->GetPeriod());
2131 esd->SetTimeStamp(fRawReader->GetTimestamp());
2132 esd->SetEventType(fRawReader->GetType());
2138 //_____________________________________________________________________________
2139 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2141 // check whether detName is contained in detectors
2142 // if yes, it is removed from detectors
2144 // check if all detectors are selected
2145 if ((detectors.CompareTo("ALL") == 0) ||
2146 detectors.BeginsWith("ALL ") ||
2147 detectors.EndsWith(" ALL") ||
2148 detectors.Contains(" ALL ")) {
2153 // search for the given detector
2154 Bool_t result = kFALSE;
2155 if ((detectors.CompareTo(detName) == 0) ||
2156 detectors.BeginsWith(detName+" ") ||
2157 detectors.EndsWith(" "+detName) ||
2158 detectors.Contains(" "+detName+" ")) {
2159 detectors.ReplaceAll(detName, "");
2163 // clean up the detectors string
2164 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2165 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2166 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2171 //_____________________________________________________________________________
2172 Bool_t AliReconstruction::InitRunLoader()
2174 // get or create the run loader
2176 if (gAlice) delete gAlice;
2179 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2180 // load all base libraries to get the loader classes
2181 TString libs = gSystem->GetLibraries();
2182 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2183 TString detName = fgkDetectorName[iDet];
2184 if (detName == "HLT") continue;
2185 if (libs.Contains("lib" + detName + "base.so")) continue;
2186 gSystem->Load("lib" + detName + "base.so");
2188 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2190 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2195 fRunLoader->CdGAFile();
2196 fRunLoader->LoadgAlice();
2198 //PH This is a temporary fix to give access to the kinematics
2199 //PH that is needed for the labels of ITS clusters
2200 fRunLoader->LoadHeader();
2201 fRunLoader->LoadKinematics();
2203 } else { // galice.root does not exist
2205 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2209 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2210 AliConfig::GetDefaultEventFolderName(),
2213 AliError(Form("could not create run loader in file %s",
2214 fGAliceFileName.Data()));
2218 fIsNewRunLoader = kTRUE;
2219 fRunLoader->MakeTree("E");
2221 if (fNumberOfEventsPerFile > 0)
2222 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2224 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2230 //_____________________________________________________________________________
2231 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2233 // get the reconstructor object and the loader for a detector
2235 if (fReconstructor[iDet]) return fReconstructor[iDet];
2237 // load the reconstructor object
2238 TPluginManager* pluginManager = gROOT->GetPluginManager();
2239 TString detName = fgkDetectorName[iDet];
2240 TString recName = "Ali" + detName + "Reconstructor";
2242 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2244 AliReconstructor* reconstructor = NULL;
2245 // first check if a plugin is defined for the reconstructor
2246 TPluginHandler* pluginHandler =
2247 pluginManager->FindHandler("AliReconstructor", detName);
2248 // if not, add a plugin for it
2249 if (!pluginHandler) {
2250 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2251 TString libs = gSystem->GetLibraries();
2252 if (libs.Contains("lib" + detName + "base.so") ||
2253 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2254 pluginManager->AddHandler("AliReconstructor", detName,
2255 recName, detName + "rec", recName + "()");
2257 pluginManager->AddHandler("AliReconstructor", detName,
2258 recName, detName, recName + "()");
2260 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2262 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2263 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2265 if (reconstructor) {
2266 TObject* obj = fOptions.FindObject(detName.Data());
2267 if (obj) reconstructor->SetOption(obj->GetTitle());
2268 reconstructor->Init();
2269 fReconstructor[iDet] = reconstructor;
2272 // get or create the loader
2273 if (detName != "HLT") {
2274 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2275 if (!fLoader[iDet]) {
2276 AliConfig::Instance()
2277 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2279 // first check if a plugin is defined for the loader
2281 pluginManager->FindHandler("AliLoader", detName);
2282 // if not, add a plugin for it
2283 if (!pluginHandler) {
2284 TString loaderName = "Ali" + detName + "Loader";
2285 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2286 pluginManager->AddHandler("AliLoader", detName,
2287 loaderName, detName + "base",
2288 loaderName + "(const char*, TFolder*)");
2289 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2291 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2293 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2294 fRunLoader->GetEventFolder());
2296 if (!fLoader[iDet]) { // use default loader
2297 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2299 if (!fLoader[iDet]) {
2300 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2301 if (fStopOnError) return NULL;
2303 fRunLoader->AddLoader(fLoader[iDet]);
2304 fRunLoader->CdGAFile();
2305 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2306 fRunLoader->Write(0, TObject::kOverwrite);
2311 return reconstructor;
2314 //_____________________________________________________________________________
2315 Bool_t AliReconstruction::CreateVertexer()
2317 // create the vertexer
2320 AliReconstructor* itsReconstructor = GetReconstructor(0);
2321 if (itsReconstructor) {
2322 fVertexer = itsReconstructor->CreateVertexer();
2325 AliWarning("couldn't create a vertexer for ITS");
2326 if (fStopOnError) return kFALSE;
2332 //_____________________________________________________________________________
2333 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2335 // create the trackers
2337 TString detStr = detectors;
2338 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2339 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2340 AliReconstructor* reconstructor = GetReconstructor(iDet);
2341 if (!reconstructor) continue;
2342 TString detName = fgkDetectorName[iDet];
2343 if (detName == "HLT") {
2344 fRunHLTTracking = kTRUE;
2347 if (detName == "MUON") {
2348 fRunMuonTracking = kTRUE;
2353 fTracker[iDet] = reconstructor->CreateTracker();
2354 if (!fTracker[iDet] && (iDet < 7)) {
2355 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2356 if (fStopOnError) return kFALSE;
2358 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2364 //_____________________________________________________________________________
2365 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2367 // delete trackers and the run loader and close and delete the file
2369 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2370 delete fReconstructor[iDet];
2371 fReconstructor[iDet] = NULL;
2372 fLoader[iDet] = NULL;
2373 delete fTracker[iDet];
2374 fTracker[iDet] = NULL;
2375 // delete fQADataMaker[iDet];
2376 // fQADataMaker[iDet] = NULL;
2381 if (ftVertexer) delete ftVertexer;
2384 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2385 delete fDiamondProfile;
2386 fDiamondProfile = NULL;
2387 delete fDiamondProfileTPC;
2388 fDiamondProfileTPC = NULL;
2398 if (fParentRawReader) delete fParentRawReader;
2399 fParentRawReader=NULL;
2409 gSystem->Unlink("AliESDs.old.root");
2414 //_____________________________________________________________________________
2416 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2418 // read the ESD event from a file
2420 if (!esd) return kFALSE;
2422 sprintf(fileName, "ESD_%d.%d_%s.root",
2423 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2424 if (gSystem->AccessPathName(fileName)) return kFALSE;
2426 AliInfo(Form("reading ESD from file %s", fileName));
2427 AliDebug(1, Form("reading ESD from file %s", fileName));
2428 TFile* file = TFile::Open(fileName);
2429 if (!file || !file->IsOpen()) {
2430 AliError(Form("opening %s failed", fileName));
2437 esd = (AliESDEvent*) file->Get("ESD");
2446 //_____________________________________________________________________________
2447 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2449 // write the ESD event to a file
2453 sprintf(fileName, "ESD_%d.%d_%s.root",
2454 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2456 AliDebug(1, Form("writing ESD to file %s", fileName));
2457 TFile* file = TFile::Open(fileName, "recreate");
2458 if (!file || !file->IsOpen()) {
2459 AliError(Form("opening %s failed", fileName));
2468 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2470 // Write space-points which are then used in the alignment procedures
2471 // For the moment only ITS, TPC, TRD and TOF
2473 Int_t ntracks = esd->GetNumberOfTracks();
2474 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2476 AliESDtrack *track = esd->GetTrack(itrack);
2479 for (Int_t iDet = 3; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2480 nsp += track->GetNcls(iDet);
2482 if (iDet==0) { // ITS "extra" clusters
2483 track->GetClusters(iDet,idx);
2484 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2489 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2490 track->SetTrackPointArray(sp);
2492 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2493 AliTracker *tracker = fTracker[iDet];
2494 if (!tracker) continue;
2495 Int_t nspdet = track->GetClusters(iDet,idx);
2497 if (iDet==0) // ITS "extra" clusters
2498 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2500 if (nspdet <= 0) continue;
2504 while (isp2 < nspdet) {
2505 Bool_t isvalid=kTRUE;
2507 Int_t index=idx[isp++];
2508 if (index < 0) continue;
2510 TString dets = fgkDetectorName[iDet];
2511 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2512 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2513 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2514 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2515 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2517 isvalid = tracker->GetTrackPoint(index,p);
2520 if (!isvalid) continue;
2521 sp->AddPoint(isptrack,&p); isptrack++;
2528 //_____________________________________________________________________________
2529 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2531 // The method reads the raw-data error log
2532 // accumulated within the rawReader.
2533 // It extracts the raw-data errors related to
2534 // the current event and stores them into
2535 // a TClonesArray inside the esd object.
2537 if (!fRawReader) return;
2539 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2541 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2543 if (iEvent != log->GetEventNumber()) continue;
2545 esd->AddRawDataErrorLog(log);
2550 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2551 // Dump a file content into a char in TNamed
2553 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2554 Int_t kBytes = (Int_t)in.tellg();
2555 printf("Size: %d \n",kBytes);
2558 char* memblock = new char [kBytes];
2559 in.seekg (0, ios::beg);
2560 in.read (memblock, kBytes);
2562 TString fData(memblock,kBytes);
2563 fn = new TNamed(pName,fData);
2564 printf("fData Size: %d \n",fData.Sizeof());
2565 printf("pName Size: %d \n",pName.Sizeof());
2566 printf("fn Size: %d \n",fn->Sizeof());
2570 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2576 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2577 // This is not really needed in AliReconstruction at the moment
2578 // but can serve as a template
2580 TList *fList = fTree->GetUserInfo();
2581 TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2582 printf("fn Size: %d \n",fn->Sizeof());
2584 TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2585 const char* cdata = fn->GetTitle();
2586 printf("fTmp Size %d\n",fTmp.Sizeof());
2588 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2589 printf("calculated size %d\n",size);
2590 ofstream out(pName.Data(),ios::out | ios::binary);
2591 out.write(cdata,size);
2596 //_____________________________________________________________________________
2597 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2599 // get the quality assurance data maker object and the loader for a detector
2601 if (fQADataMaker[iDet])
2602 return fQADataMaker[iDet];
2604 AliQADataMakerRec * qadm = NULL;
2605 if (iDet == fgkNDetectors) { //Global QA
2606 qadm = new AliGlobalQADataMaker();
2607 fQADataMaker[iDet] = qadm;
2611 // load the QA data maker object
2612 TPluginManager* pluginManager = gROOT->GetPluginManager();
2613 TString detName = fgkDetectorName[iDet];
2614 TString qadmName = "Ali" + detName + "QADataMakerRec";
2615 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT"))
2618 // first check if a plugin is defined for the quality assurance data maker
2619 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2620 // if not, add a plugin for it
2621 if (!pluginHandler) {
2622 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
2623 TString libs = gSystem->GetLibraries();
2624 if (libs.Contains("lib" + detName + "base.so") ||
2625 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2626 pluginManager->AddHandler("AliQADataMakerRec", detName,
2627 qadmName, detName + "qadm", qadmName + "()");
2629 pluginManager->AddHandler("AliQADataMakerRec", detName,
2630 qadmName, detName, qadmName + "()");
2632 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2634 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2635 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
2638 fQADataMaker[iDet] = qadm;
2643 //_____________________________________________________________________________
2644 Bool_t AliReconstruction::RunQA(AliESDEvent *& esd)
2646 // run the Quality Assurance data producer
2648 AliCodeTimerAuto("")
2649 TString detStr = fQADetectors ;
2650 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2651 if (!IsSelected(fgkDetectorName[iDet], detStr))
2653 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
2656 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2657 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2659 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
2660 qadm->Exec(AliQA::kESDS, esd) ;
2663 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2665 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2666 AliError(Form("the following detectors were not found: %s",
2676 //_____________________________________________________________________________
2677 void AliReconstruction::CheckQA()
2679 // check the QA of SIM for this run and remove the detectors
2680 // with status Fatal
2682 TString newRunLocalReconstruction ;
2683 TString newRunTracking ;
2684 TString newFillESD ;
2686 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2687 TString detName(AliQA::GetDetName(iDet)) ;
2688 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2689 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2690 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2692 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2693 fRunLocalReconstruction.Contains("ALL") ) {
2694 newRunLocalReconstruction += detName ;
2695 newRunLocalReconstruction += " " ;
2697 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2698 fRunTracking.Contains("ALL") ) {
2699 newRunTracking += detName ;
2700 newRunTracking += " " ;
2702 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2703 fFillESD.Contains("ALL") ) {
2704 newFillESD += detName ;
2709 fRunLocalReconstruction = newRunLocalReconstruction ;
2710 fRunTracking = newRunTracking ;
2711 fFillESD = newFillESD ;
2714 //_____________________________________________________________________________
2715 Int_t AliReconstruction::GetDetIndex(const char* detector)
2717 // return the detector index corresponding to detector
2719 for (index = 0; index < fgkNDetectors ; index++) {
2720 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2725 //_____________________________________________________________________________
2726 Bool_t AliReconstruction::FinishPlaneEff() {
2728 // Here execute all the necessary operationis, at the end of the tracking phase,
2729 // in case that evaluation of PlaneEfficiencies was required for some detector.
2730 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2732 // This Preliminary version works only FOR ITS !!!!!
2733 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2736 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2739 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2740 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2741 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2742 if(fTracker[iDet]) {
2743 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2744 TString name=planeeff->GetName();
2746 TFile* pefile = TFile::Open(name, "RECREATE");
2747 ret=(Bool_t)planeeff->Write();
2749 if(planeeff->GetCreateHistos()) {
2750 TString hname=planeeff->GetName();
2751 hname+="Histo.root";
2752 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
2758 //_____________________________________________________________________________
2759 Bool_t AliReconstruction::InitPlaneEff() {
2761 // Here execute all the necessary operations, before of the tracking phase,
2762 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2763 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2764 // which should be updated/recalculated.
2766 // This Preliminary version will work only FOR ITS !!!!!
2767 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2770 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2772 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2776 //_____________________________________________________________________________
2777 Bool_t AliReconstruction::InitAliEVE()
2779 // This method should be called only in case
2780 // AliReconstruction is run
2781 // within the alieve environment.
2782 // It will initialize AliEVE in a way
2783 // so that it can visualize event processed
2784 // by AliReconstruction.
2785 // The return flag shows whenever the
2786 // AliEVE initialization was successful or not.
2789 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2790 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2791 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2793 gROOT->ProcessLine("if (!gAliEveEvent) {gAliEveEvent = new AliEveEventManager();gAliEveEvent->SetAutoLoad(kTRUE);gAliEveEvent->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(gAliEveEvent);};");
2794 gROOT->ProcessLine("alieve_online_init()");
2799 //_____________________________________________________________________________
2800 void AliReconstruction::RunAliEVE()
2802 // Runs AliEVE visualisation of
2803 // the current event.
2804 // Should be executed only after
2805 // successful initialization of AliEVE.
2807 AliInfo("Running AliEVE...");
2808 gROOT->ProcessLine(Form("gAliEveEvent->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2809 gROOT->ProcessLine("gAliEveEvent->StartStopAutoLoadTimer();");
2813 //_____________________________________________________________________________
2814 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2816 // Allows to run QA for a selected set of detectors
2817 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2818 // all selected detectors run the same selected tasks
2820 if (!detAndAction.Contains(":")) {
2821 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2825 Int_t colon = detAndAction.Index(":") ;
2826 fQADetectors = detAndAction(0, colon) ;
2827 if (fQADetectors.Contains("ALL") )
2828 fQADetectors = fFillESD ;
2829 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2830 if (fQATasks.Contains("ALL") ) {
2831 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2833 fQATasks.ToUpper() ;
2835 if ( fQATasks.Contains("RAW") )
2836 tempo = Form("%d ", AliQA::kRAWS) ;
2837 if ( fQATasks.Contains("RECPOINT") )
2838 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2839 if ( fQATasks.Contains("ESD") )
2840 tempo += Form("%d ", AliQA::kESDS) ;
2842 if (fQATasks.IsNull()) {
2843 AliInfo("No QA requested\n") ;
2848 TString tempo(fQATasks) ;
2849 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2850 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2851 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2852 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;