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);
1331 if (fRunAliEVE) RunAliEVE();
1333 if (fCheckPointLevel > 0) WriteESD(fesd, "final");
1336 if (fWriteESDfriend) {
1337 fesdf->~AliESDfriend();
1338 new (fesdf) AliESDfriend(); // Reset...
1341 ProcInfo_t ProcInfo;
1342 gSystem->GetProcInfo(&ProcInfo);
1343 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1346 // End of cycle for the in-loop
1350 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1351 if (!IsSelected(fgkDetectorName[iDet], fQADetectors))
1353 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1356 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1357 qadm->EndOfCycle(AliQA::kRECPOINTS);
1358 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1359 qadm->EndOfCycle(AliQA::kESDS);
1364 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1366 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1367 qadm->EndOfCycle(AliQA::kRECPOINTS);
1368 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1369 qadm->EndOfCycle(AliQA::kESDS);
1378 //_____________________________________________________________________________
1379 Bool_t AliReconstruction::FinishRun()
1382 // Called after the exit
1383 // from the event loop
1384 AliCodeTimerAuto("");
1386 if (fIsNewRunLoader) { // galice.root didn't exist
1387 fRunLoader->WriteHeader("OVERWRITE");
1388 fRunLoader->CdGAFile();
1389 fRunLoader->Write(0, TObject::kOverwrite);
1392 ftree->GetUserInfo()->Add(fesd);
1393 fhlttree->GetUserInfo()->Add(fhltesd);
1395 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1396 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1398 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1399 cdbMapCopy->SetOwner(1);
1400 cdbMapCopy->SetName("cdbMap");
1401 TIter iter(cdbMap->GetTable());
1404 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1405 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1406 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1407 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1410 TList *cdbListCopy = new TList();
1411 cdbListCopy->SetOwner(1);
1412 cdbListCopy->SetName("cdbList");
1414 TIter iter2(cdbList);
1417 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1418 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1421 ftree->GetUserInfo()->Add(cdbMapCopy);
1422 ftree->GetUserInfo()->Add(cdbListCopy);
1425 if(fESDPar.Contains("ESD.par")){
1426 AliInfo("Attaching ESD.par to Tree");
1427 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1428 ftree->GetUserInfo()->Add(fn);
1434 if (fWriteESDfriend)
1435 ftree->SetBranchStatus("ESDfriend*",0);
1436 // we want to have only one tree version number
1437 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1440 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1441 if (fRunPlaneEff && !FinishPlaneEff()) {
1442 AliWarning("Finish PlaneEff evaluation failed");
1446 CleanUp(ffile, ffileOld);
1449 AliWarning("AOD creation not supported anymore during reconstruction. See ANALYSIS/AliAnalysisTaskESDfilter.cxx instead.");
1452 // Create tags for the events in the ESD tree (the ESD tree is always present)
1453 // In case of empty events the tags will contain dummy values
1454 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1455 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1457 AliWarning("AOD tag creation not supported anymore during reconstruction.");
1460 //Finish QA and end of cycle for out-of-loop QA
1463 AliQADataMakerSteer qas;
1464 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1465 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS, fSameQACycle);
1467 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1468 qas.Run(fRunLocalReconstruction.Data(), AliQA::kESDS, fSameQACycle);
1470 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1472 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1473 qadm->EndOfCycle(AliQA::kRECPOINTS);
1474 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1475 qadm->EndOfCycle(AliQA::kESDS);
1482 // Cleanup of CDB manager: cache and active storages!
1483 AliCDBManager::Instance()->ClearCache();
1489 //_____________________________________________________________________________
1490 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1492 // run the local reconstruction
1493 static Int_t eventNr=0;
1494 AliCodeTimerAuto("")
1496 // AliCDBManager* man = AliCDBManager::Instance();
1497 // Bool_t origCache = man->GetCacheFlag();
1499 // TString detStr = detectors;
1500 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1501 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1502 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1503 // if (!reconstructor) continue;
1504 // if (reconstructor->HasLocalReconstruction()) continue;
1506 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1507 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1509 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1510 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1512 // man->SetCacheFlag(kTRUE);
1513 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1514 // man->GetAll(calibPath); // entries are cached!
1516 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1518 // if (fRawReader) {
1519 // fRawReader->RewindEvents();
1520 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1522 // reconstructor->Reconstruct(fRunLoader);
1525 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1526 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1528 // // unload calibration data
1529 // man->UnloadFromCache(calibPath);
1530 // //man->ClearCache();
1533 // man->SetCacheFlag(origCache);
1535 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1536 // AliError(Form("the following detectors were not found: %s",
1538 // if (fStopOnError) return kFALSE;
1545 //_____________________________________________________________________________
1546 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1548 // run the local reconstruction
1550 static Int_t eventNr=0;
1551 AliCodeTimerAuto("")
1553 TString detStr = detectors;
1554 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1555 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1556 AliReconstructor* reconstructor = GetReconstructor(iDet);
1557 if (!reconstructor) continue;
1558 AliLoader* loader = fLoader[iDet];
1559 // Matthias April 2008: temporary fix to run HLT reconstruction
1560 // although the HLT loader is missing
1561 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1563 reconstructor->Reconstruct(fRawReader, NULL);
1566 reconstructor->Reconstruct(dummy, NULL);
1571 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1574 // conversion of digits
1575 if (fRawReader && reconstructor->HasDigitConversion()) {
1576 AliInfo(Form("converting raw data digits into root objects for %s",
1577 fgkDetectorName[iDet]));
1578 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1579 fgkDetectorName[iDet]));
1580 loader->LoadDigits("update");
1581 loader->CleanDigits();
1582 loader->MakeDigitsContainer();
1583 TTree* digitsTree = loader->TreeD();
1584 reconstructor->ConvertDigits(fRawReader, digitsTree);
1585 loader->WriteDigits("OVERWRITE");
1586 loader->UnloadDigits();
1588 // local reconstruction
1589 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1590 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1591 loader->LoadRecPoints("update");
1592 loader->CleanRecPoints();
1593 loader->MakeRecPointsContainer();
1594 TTree* clustersTree = loader->TreeR();
1595 if (fRawReader && !reconstructor->HasDigitConversion()) {
1596 reconstructor->Reconstruct(fRawReader, clustersTree);
1598 loader->LoadDigits("read");
1599 TTree* digitsTree = loader->TreeD();
1601 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1602 if (fStopOnError) return kFALSE;
1604 reconstructor->Reconstruct(digitsTree, clustersTree);
1606 loader->UnloadDigits();
1609 // In-loop QA for local reconstrucion
1610 TString detQAStr(fQADetectors) ;
1611 if (fRunQA && fInLoopQA && IsSelected(fgkDetectorName[iDet], detQAStr) ) {
1612 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1615 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1617 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1619 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1620 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1622 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1626 loader->WriteRecPoints("OVERWRITE");
1627 loader->UnloadRecPoints();
1628 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1631 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1632 AliError(Form("the following detectors were not found: %s",
1634 if (fStopOnError) return kFALSE;
1640 //_____________________________________________________________________________
1641 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1643 // run the barrel tracking
1645 AliCodeTimerAuto("")
1647 AliESDVertex* vertex = NULL;
1648 Double_t vtxPos[3] = {0, 0, 0};
1649 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1650 TArrayF mcVertex(3);
1651 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1652 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1653 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1657 AliInfo("running the ITS vertex finder");
1659 fLoader[0]->LoadRecPoints();
1660 TTree* cltree = fLoader[0]->TreeR();
1662 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1663 vertex = fVertexer->FindVertexForCurrentEvent(cltree);
1666 AliError("Can't get the ITS cluster tree");
1668 fLoader[0]->UnloadRecPoints();
1671 AliError("Can't get the ITS loader");
1674 AliWarning("Vertex not found");
1675 vertex = new AliESDVertex();
1676 vertex->SetName("default");
1679 vertex->SetName("reconstructed");
1683 AliInfo("getting the primary vertex from MC");
1684 vertex = new AliESDVertex(vtxPos, vtxErr);
1688 vertex->GetXYZ(vtxPos);
1689 vertex->GetSigmaXYZ(vtxErr);
1691 AliWarning("no vertex reconstructed");
1692 vertex = new AliESDVertex(vtxPos, vtxErr);
1694 esd->SetPrimaryVertexSPD(vertex);
1695 // if SPD multiplicity has been determined, it is stored in the ESD
1696 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1697 if(mult)esd->SetMultiplicity(mult);
1699 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1700 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1707 //_____________________________________________________________________________
1708 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1710 // run the HLT barrel tracking
1712 AliCodeTimerAuto("")
1715 AliError("Missing runLoader!");
1719 AliInfo("running HLT tracking");
1721 // Get a pointer to the HLT reconstructor
1722 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1723 if (!reconstructor) return kFALSE;
1726 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1727 TString detName = fgkDetectorName[iDet];
1728 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1729 reconstructor->SetOption(detName.Data());
1730 AliTracker *tracker = reconstructor->CreateTracker();
1732 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1733 if (fStopOnError) return kFALSE;
1737 Double_t vtxErr[3]={0.005,0.005,0.010};
1738 const AliESDVertex *vertex = esd->GetVertex();
1739 vertex->GetXYZ(vtxPos);
1740 tracker->SetVertex(vtxPos,vtxErr);
1742 fLoader[iDet]->LoadRecPoints("read");
1743 TTree* tree = fLoader[iDet]->TreeR();
1745 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1748 tracker->LoadClusters(tree);
1750 if (tracker->Clusters2Tracks(esd) != 0) {
1751 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1755 tracker->UnloadClusters();
1763 //_____________________________________________________________________________
1764 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1766 // run the muon spectrometer tracking
1768 AliCodeTimerAuto("")
1771 AliError("Missing runLoader!");
1774 Int_t iDet = 7; // for MUON
1776 AliInfo("is running...");
1778 // Get a pointer to the MUON reconstructor
1779 AliReconstructor *reconstructor = GetReconstructor(iDet);
1780 if (!reconstructor) return kFALSE;
1783 TString detName = fgkDetectorName[iDet];
1784 AliDebug(1, Form("%s tracking", detName.Data()));
1785 AliTracker *tracker = reconstructor->CreateTracker();
1787 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1792 fLoader[iDet]->LoadRecPoints("read");
1794 tracker->LoadClusters(fLoader[iDet]->TreeR());
1796 Int_t rv = tracker->Clusters2Tracks(esd);
1800 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1804 fLoader[iDet]->UnloadRecPoints();
1806 tracker->UnloadClusters();
1814 //_____________________________________________________________________________
1815 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1817 // run the barrel tracking
1818 static Int_t eventNr=0;
1819 AliCodeTimerAuto("")
1821 AliInfo("running tracking");
1823 //Fill the ESD with the T0 info (will be used by the TOF)
1824 if (fReconstructor[11] && fLoader[11]) {
1825 fLoader[11]->LoadRecPoints("READ");
1826 TTree *treeR = fLoader[11]->TreeR();
1827 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1830 // pass 1: TPC + ITS inwards
1831 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1832 if (!fTracker[iDet]) continue;
1833 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1836 fLoader[iDet]->LoadRecPoints("read");
1837 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1838 TTree* tree = fLoader[iDet]->TreeR();
1840 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1843 fTracker[iDet]->LoadClusters(tree);
1844 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1846 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1847 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1850 if (fCheckPointLevel > 1) {
1851 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1853 // preliminary PID in TPC needed by the ITS tracker
1855 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1856 AliESDpid::MakePID(esd);
1858 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1861 // pass 2: ALL backwards
1863 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1864 if (!fTracker[iDet]) continue;
1865 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1868 if (iDet > 1) { // all except ITS, TPC
1870 fLoader[iDet]->LoadRecPoints("read");
1871 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1872 tree = fLoader[iDet]->TreeR();
1874 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1877 fTracker[iDet]->LoadClusters(tree);
1878 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1882 if (iDet>1) // start filling residuals for the "outer" detectors
1883 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1885 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1886 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1889 if (fCheckPointLevel > 1) {
1890 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1894 if (iDet > 2) { // all except ITS, TPC, TRD
1895 fTracker[iDet]->UnloadClusters();
1896 fLoader[iDet]->UnloadRecPoints();
1898 // updated PID in TPC needed by the ITS tracker -MI
1900 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1901 AliESDpid::MakePID(esd);
1903 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1905 //stop filling residuals for the "outer" detectors
1906 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1908 // write space-points to the ESD in case alignment data output
1910 if (fWriteAlignmentData)
1911 WriteAlignmentData(esd);
1913 // pass 3: TRD + TPC + ITS refit inwards
1915 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1916 if (!fTracker[iDet]) continue;
1917 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1920 if (iDet<2) // start filling residuals for TPC and ITS
1921 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1923 if (fTracker[iDet]->RefitInward(esd) != 0) {
1924 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1927 // run postprocessing
1928 if (fTracker[iDet]->PostProcess(esd) != 0) {
1929 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1932 if (fCheckPointLevel > 1) {
1933 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1935 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1937 fTracker[iDet]->UnloadClusters();
1938 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1939 fLoader[iDet]->UnloadRecPoints();
1940 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1942 // stop filling residuals for TPC and ITS
1943 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1949 //_____________________________________________________________________________
1950 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1952 // Remove the data which are not needed for the physics analysis.
1955 Int_t nTracks=esd->GetNumberOfTracks();
1956 Int_t nV0s=esd->GetNumberOfV0s();
1958 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1960 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1961 Bool_t rc=esd->Clean(cleanPars);
1963 nTracks=esd->GetNumberOfTracks();
1964 nV0s=esd->GetNumberOfV0s();
1966 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1971 //_____________________________________________________________________________
1972 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1974 // fill the event summary data
1976 AliCodeTimerAuto("")
1977 static Int_t eventNr=0;
1978 TString detStr = detectors;
1980 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1981 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1982 AliReconstructor* reconstructor = GetReconstructor(iDet);
1983 if (!reconstructor) continue;
1984 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1985 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1986 TTree* clustersTree = NULL;
1987 if (fLoader[iDet]) {
1988 fLoader[iDet]->LoadRecPoints("read");
1989 clustersTree = fLoader[iDet]->TreeR();
1990 if (!clustersTree) {
1991 AliError(Form("Can't get the %s clusters tree",
1992 fgkDetectorName[iDet]));
1993 if (fStopOnError) return kFALSE;
1996 if (fRawReader && !reconstructor->HasDigitConversion()) {
1997 reconstructor->FillESD(fRawReader, clustersTree, esd);
1999 TTree* digitsTree = NULL;
2000 if (fLoader[iDet]) {
2001 fLoader[iDet]->LoadDigits("read");
2002 digitsTree = fLoader[iDet]->TreeD();
2004 AliError(Form("Can't get the %s digits tree",
2005 fgkDetectorName[iDet]));
2006 if (fStopOnError) return kFALSE;
2009 reconstructor->FillESD(digitsTree, clustersTree, esd);
2010 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2012 if (fLoader[iDet]) {
2013 fLoader[iDet]->UnloadRecPoints();
2016 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
2020 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2021 AliError(Form("the following detectors were not found: %s",
2023 if (fStopOnError) return kFALSE;
2025 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2030 //_____________________________________________________________________________
2031 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2033 // Reads the trigger decision which is
2034 // stored in Trigger.root file and fills
2035 // the corresponding esd entries
2037 AliCodeTimerAuto("")
2039 AliInfo("Filling trigger information into the ESD");
2041 AliCentralTrigger *aCTP = NULL;
2044 AliCTPRawStream input(fRawReader);
2045 if (!input.Next()) {
2046 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 !");
2047 ULong64_t mask = (((ULong64_t)fRawReader->GetTriggerPattern()[1]) << 32) +
2048 fRawReader->GetTriggerPattern()[0];
2049 esd->SetTriggerMask(mask);
2050 esd->SetTriggerCluster(0);
2053 esd->SetTriggerMask(input.GetClassMask());
2054 esd->SetTriggerCluster(input.GetClusterMask());
2057 aCTP = new AliCentralTrigger();
2058 TString configstr("");
2059 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
2060 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
2066 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
2068 if (!runloader->LoadTrigger()) {
2069 aCTP = runloader->GetTrigger();
2070 esd->SetTriggerMask(aCTP->GetClassMask());
2071 esd->SetTriggerCluster(aCTP->GetClusterMask());
2074 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
2079 AliError("No run loader is available! The trigger information is not stored in the ESD !");
2084 // Now fill the trigger class names into AliESDRun object
2085 AliTriggerConfiguration *config = aCTP->GetConfiguration();
2087 AliError("No trigger configuration has been found! The trigger classes information will not be stored in ESD!");
2088 if (fRawReader) delete aCTP;
2092 const TObjArray& classesArray = config->GetClasses();
2093 Int_t nclasses = classesArray.GetEntriesFast();
2094 for( Int_t j=0; j<nclasses; j++ ) {
2095 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
2096 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
2097 esd->SetTriggerClass(trclass->GetName(),trindex);
2100 if (fRawReader) delete aCTP;
2108 //_____________________________________________________________________________
2109 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2112 // Filling information from RawReader Header
2115 AliInfo("Filling information from RawReader Header");
2116 esd->SetBunchCrossNumber(0);
2117 esd->SetOrbitNumber(0);
2118 esd->SetPeriodNumber(0);
2119 esd->SetTimeStamp(0);
2120 esd->SetEventType(0);
2121 const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
2124 const UInt_t *id = eventHeader->GetP("Id");
2125 esd->SetBunchCrossNumber((id)[1]&0x00000fff);
2126 esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
2127 esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
2129 esd->SetTimeStamp((eventHeader->Get("Timestamp")));
2130 esd->SetEventType((eventHeader->Get("Type")));
2137 //_____________________________________________________________________________
2138 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2140 // check whether detName is contained in detectors
2141 // if yes, it is removed from detectors
2143 // check if all detectors are selected
2144 if ((detectors.CompareTo("ALL") == 0) ||
2145 detectors.BeginsWith("ALL ") ||
2146 detectors.EndsWith(" ALL") ||
2147 detectors.Contains(" ALL ")) {
2152 // search for the given detector
2153 Bool_t result = kFALSE;
2154 if ((detectors.CompareTo(detName) == 0) ||
2155 detectors.BeginsWith(detName+" ") ||
2156 detectors.EndsWith(" "+detName) ||
2157 detectors.Contains(" "+detName+" ")) {
2158 detectors.ReplaceAll(detName, "");
2162 // clean up the detectors string
2163 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2164 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2165 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2170 //_____________________________________________________________________________
2171 Bool_t AliReconstruction::InitRunLoader()
2173 // get or create the run loader
2175 if (gAlice) delete gAlice;
2178 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2179 // load all base libraries to get the loader classes
2180 TString libs = gSystem->GetLibraries();
2181 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2182 TString detName = fgkDetectorName[iDet];
2183 if (detName == "HLT") continue;
2184 if (libs.Contains("lib" + detName + "base.so")) continue;
2185 gSystem->Load("lib" + detName + "base.so");
2187 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2189 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2194 fRunLoader->CdGAFile();
2195 fRunLoader->LoadgAlice();
2197 //PH This is a temporary fix to give access to the kinematics
2198 //PH that is needed for the labels of ITS clusters
2199 fRunLoader->LoadHeader();
2200 fRunLoader->LoadKinematics();
2202 } else { // galice.root does not exist
2204 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2208 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2209 AliConfig::GetDefaultEventFolderName(),
2212 AliError(Form("could not create run loader in file %s",
2213 fGAliceFileName.Data()));
2217 fIsNewRunLoader = kTRUE;
2218 fRunLoader->MakeTree("E");
2220 if (fNumberOfEventsPerFile > 0)
2221 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2223 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2229 //_____________________________________________________________________________
2230 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2232 // get the reconstructor object and the loader for a detector
2234 if (fReconstructor[iDet]) return fReconstructor[iDet];
2236 // load the reconstructor object
2237 TPluginManager* pluginManager = gROOT->GetPluginManager();
2238 TString detName = fgkDetectorName[iDet];
2239 TString recName = "Ali" + detName + "Reconstructor";
2241 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2243 AliReconstructor* reconstructor = NULL;
2244 // first check if a plugin is defined for the reconstructor
2245 TPluginHandler* pluginHandler =
2246 pluginManager->FindHandler("AliReconstructor", detName);
2247 // if not, add a plugin for it
2248 if (!pluginHandler) {
2249 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2250 TString libs = gSystem->GetLibraries();
2251 if (libs.Contains("lib" + detName + "base.so") ||
2252 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2253 pluginManager->AddHandler("AliReconstructor", detName,
2254 recName, detName + "rec", recName + "()");
2256 pluginManager->AddHandler("AliReconstructor", detName,
2257 recName, detName, recName + "()");
2259 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2261 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2262 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2264 if (reconstructor) {
2265 TObject* obj = fOptions.FindObject(detName.Data());
2266 if (obj) reconstructor->SetOption(obj->GetTitle());
2267 reconstructor->Init();
2268 fReconstructor[iDet] = reconstructor;
2271 // get or create the loader
2272 if (detName != "HLT") {
2273 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2274 if (!fLoader[iDet]) {
2275 AliConfig::Instance()
2276 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2278 // first check if a plugin is defined for the loader
2280 pluginManager->FindHandler("AliLoader", detName);
2281 // if not, add a plugin for it
2282 if (!pluginHandler) {
2283 TString loaderName = "Ali" + detName + "Loader";
2284 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2285 pluginManager->AddHandler("AliLoader", detName,
2286 loaderName, detName + "base",
2287 loaderName + "(const char*, TFolder*)");
2288 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2290 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2292 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2293 fRunLoader->GetEventFolder());
2295 if (!fLoader[iDet]) { // use default loader
2296 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2298 if (!fLoader[iDet]) {
2299 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2300 if (fStopOnError) return NULL;
2302 fRunLoader->AddLoader(fLoader[iDet]);
2303 fRunLoader->CdGAFile();
2304 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2305 fRunLoader->Write(0, TObject::kOverwrite);
2310 return reconstructor;
2313 //_____________________________________________________________________________
2314 Bool_t AliReconstruction::CreateVertexer()
2316 // create the vertexer
2319 AliReconstructor* itsReconstructor = GetReconstructor(0);
2320 if (itsReconstructor) {
2321 fVertexer = itsReconstructor->CreateVertexer();
2324 AliWarning("couldn't create a vertexer for ITS");
2325 if (fStopOnError) return kFALSE;
2331 //_____________________________________________________________________________
2332 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2334 // create the trackers
2336 TString detStr = detectors;
2337 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2338 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2339 AliReconstructor* reconstructor = GetReconstructor(iDet);
2340 if (!reconstructor) continue;
2341 TString detName = fgkDetectorName[iDet];
2342 if (detName == "HLT") {
2343 fRunHLTTracking = kTRUE;
2346 if (detName == "MUON") {
2347 fRunMuonTracking = kTRUE;
2352 fTracker[iDet] = reconstructor->CreateTracker();
2353 if (!fTracker[iDet] && (iDet < 7)) {
2354 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2355 if (fStopOnError) return kFALSE;
2357 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2363 //_____________________________________________________________________________
2364 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2366 // delete trackers and the run loader and close and delete the file
2368 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2369 delete fReconstructor[iDet];
2370 fReconstructor[iDet] = NULL;
2371 fLoader[iDet] = NULL;
2372 delete fTracker[iDet];
2373 fTracker[iDet] = NULL;
2374 // delete fQADataMaker[iDet];
2375 // fQADataMaker[iDet] = NULL;
2380 if (ftVertexer) delete ftVertexer;
2383 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2384 delete fDiamondProfile;
2385 fDiamondProfile = NULL;
2386 delete fDiamondProfileTPC;
2387 fDiamondProfileTPC = NULL;
2397 if (fParentRawReader) delete fParentRawReader;
2398 fParentRawReader=NULL;
2408 gSystem->Unlink("AliESDs.old.root");
2413 //_____________________________________________________________________________
2415 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2417 // read the ESD event from a file
2419 if (!esd) return kFALSE;
2421 sprintf(fileName, "ESD_%d.%d_%s.root",
2422 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2423 if (gSystem->AccessPathName(fileName)) return kFALSE;
2425 AliInfo(Form("reading ESD from file %s", fileName));
2426 AliDebug(1, Form("reading ESD from file %s", fileName));
2427 TFile* file = TFile::Open(fileName);
2428 if (!file || !file->IsOpen()) {
2429 AliError(Form("opening %s failed", fileName));
2436 esd = (AliESDEvent*) file->Get("ESD");
2445 //_____________________________________________________________________________
2446 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2448 // write the ESD event to a file
2452 sprintf(fileName, "ESD_%d.%d_%s.root",
2453 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2455 AliDebug(1, Form("writing ESD to file %s", fileName));
2456 TFile* file = TFile::Open(fileName, "recreate");
2457 if (!file || !file->IsOpen()) {
2458 AliError(Form("opening %s failed", fileName));
2467 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2469 // Write space-points which are then used in the alignment procedures
2470 // For the moment only ITS, TRD and TPC
2472 // Load TOF clusters
2474 fLoader[3]->LoadRecPoints("read");
2475 TTree* tree = fLoader[3]->TreeR();
2477 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2480 fTracker[3]->LoadClusters(tree);
2482 Int_t ntracks = esd->GetNumberOfTracks();
2483 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2485 AliESDtrack *track = esd->GetTrack(itrack);
2488 for (Int_t iDet = 3; iDet >= 0; iDet--)
2489 nsp += track->GetNcls(iDet);
2491 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2492 track->SetTrackPointArray(sp);
2494 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2495 AliTracker *tracker = fTracker[iDet];
2496 if (!tracker) continue;
2497 Int_t nspdet = track->GetNcls(iDet);
2498 if (nspdet <= 0) continue;
2499 track->GetClusters(iDet,idx);
2503 while (isp2 < nspdet) {
2504 Bool_t isvalid=kTRUE;
2506 Int_t index=idx[isp++];
2507 if (index < 0) continue;
2509 TString dets = fgkDetectorName[iDet];
2510 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2511 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2512 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2513 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2514 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2516 isvalid = tracker->GetTrackPoint(index,p);
2519 if (!isvalid) continue;
2520 sp->AddPoint(isptrack,&p); isptrack++;
2526 fTracker[3]->UnloadClusters();
2527 fLoader[3]->UnloadRecPoints();
2531 //_____________________________________________________________________________
2532 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2534 // The method reads the raw-data error log
2535 // accumulated within the rawReader.
2536 // It extracts the raw-data errors related to
2537 // the current event and stores them into
2538 // a TClonesArray inside the esd object.
2540 if (!fRawReader) return;
2542 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2544 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2546 if (iEvent != log->GetEventNumber()) continue;
2548 esd->AddRawDataErrorLog(log);
2553 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2554 // Dump a file content into a char in TNamed
2556 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2557 Int_t kBytes = (Int_t)in.tellg();
2558 printf("Size: %d \n",kBytes);
2561 char* memblock = new char [kBytes];
2562 in.seekg (0, ios::beg);
2563 in.read (memblock, kBytes);
2565 TString fData(memblock,kBytes);
2566 fn = new TNamed(pName,fData);
2567 printf("fData Size: %d \n",fData.Sizeof());
2568 printf("pName Size: %d \n",pName.Sizeof());
2569 printf("fn Size: %d \n",fn->Sizeof());
2573 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2579 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2580 // This is not really needed in AliReconstruction at the moment
2581 // but can serve as a template
2583 TList *fList = fTree->GetUserInfo();
2584 TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2585 printf("fn Size: %d \n",fn->Sizeof());
2587 TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2588 const char* cdata = fn->GetTitle();
2589 printf("fTmp Size %d\n",fTmp.Sizeof());
2591 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2592 printf("calculated size %d\n",size);
2593 ofstream out(pName.Data(),ios::out | ios::binary);
2594 out.write(cdata,size);
2599 //_____________________________________________________________________________
2600 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2602 // get the quality assurance data maker object and the loader for a detector
2604 if (fQADataMaker[iDet])
2605 return fQADataMaker[iDet];
2607 AliQADataMakerRec * qadm = NULL;
2608 if (iDet == fgkNDetectors) { //Global QA
2609 qadm = new AliGlobalQADataMaker();
2610 fQADataMaker[iDet] = qadm;
2614 // load the QA data maker object
2615 TPluginManager* pluginManager = gROOT->GetPluginManager();
2616 TString detName = fgkDetectorName[iDet];
2617 TString qadmName = "Ali" + detName + "QADataMakerRec";
2618 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT"))
2621 // first check if a plugin is defined for the quality assurance data maker
2622 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2623 // if not, add a plugin for it
2624 if (!pluginHandler) {
2625 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
2626 TString libs = gSystem->GetLibraries();
2627 if (libs.Contains("lib" + detName + "base.so") ||
2628 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2629 pluginManager->AddHandler("AliQADataMakerRec", detName,
2630 qadmName, detName + "qadm", qadmName + "()");
2632 pluginManager->AddHandler("AliQADataMakerRec", detName,
2633 qadmName, detName, qadmName + "()");
2635 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2637 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2638 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
2641 fQADataMaker[iDet] = qadm;
2646 //_____________________________________________________________________________
2647 Bool_t AliReconstruction::RunQA(AliESDEvent *& esd)
2649 // run the Quality Assurance data producer
2651 AliCodeTimerAuto("")
2652 TString detStr = fQADetectors ;
2653 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2654 if (!IsSelected(fgkDetectorName[iDet], detStr))
2656 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
2659 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2660 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2662 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
2663 qadm->Exec(AliQA::kESDS, esd) ;
2666 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2668 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2669 AliError(Form("the following detectors were not found: %s",
2679 //_____________________________________________________________________________
2680 void AliReconstruction::CheckQA()
2682 // check the QA of SIM for this run and remove the detectors
2683 // with status Fatal
2685 TString newRunLocalReconstruction ;
2686 TString newRunTracking ;
2687 TString newFillESD ;
2689 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2690 TString detName(AliQA::GetDetName(iDet)) ;
2691 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2692 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2693 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2695 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2696 fRunLocalReconstruction.Contains("ALL") ) {
2697 newRunLocalReconstruction += detName ;
2698 newRunLocalReconstruction += " " ;
2700 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2701 fRunTracking.Contains("ALL") ) {
2702 newRunTracking += detName ;
2703 newRunTracking += " " ;
2705 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2706 fFillESD.Contains("ALL") ) {
2707 newFillESD += detName ;
2712 fRunLocalReconstruction = newRunLocalReconstruction ;
2713 fRunTracking = newRunTracking ;
2714 fFillESD = newFillESD ;
2717 //_____________________________________________________________________________
2718 Int_t AliReconstruction::GetDetIndex(const char* detector)
2720 // return the detector index corresponding to detector
2722 for (index = 0; index < fgkNDetectors ; index++) {
2723 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2728 //_____________________________________________________________________________
2729 Bool_t AliReconstruction::FinishPlaneEff() {
2731 // Here execute all the necessary operationis, at the end of the tracking phase,
2732 // in case that evaluation of PlaneEfficiencies was required for some detector.
2733 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2735 // This Preliminary version works only FOR ITS !!!!!
2736 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2739 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2742 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2743 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2744 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2745 if(fTracker[iDet]) {
2746 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2747 TString name=planeeff->GetName();
2749 TFile* pefile = TFile::Open(name, "RECREATE");
2750 ret=(Bool_t)planeeff->Write();
2752 if(planeeff->GetCreateHistos()) {
2753 TString hname=planeeff->GetName();
2754 hname+="Histo.root";
2755 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
2761 //_____________________________________________________________________________
2762 Bool_t AliReconstruction::InitPlaneEff() {
2764 // Here execute all the necessary operations, before of the tracking phase,
2765 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2766 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2767 // which should be updated/recalculated.
2769 // This Preliminary version will work only FOR ITS !!!!!
2770 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2773 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2775 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2779 //_____________________________________________________________________________
2780 Bool_t AliReconstruction::InitAliEVE()
2782 // This method should be called only in case
2783 // AliReconstruction is run
2784 // within the alieve environment.
2785 // It will initialize AliEVE in a way
2786 // so that it can visualize event processed
2787 // by AliReconstruction.
2788 // The return flag shows whenever the
2789 // AliEVE initialization was successful or not.
2792 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2793 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2794 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2796 gROOT->ProcessLine("if (!gAliEveEvent) {gAliEveEvent = new AliEveEventManager();gAliEveEvent->SetAutoLoad(kTRUE);gAliEveEvent->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(gAliEveEvent);};");
2797 gROOT->ProcessLine("alieve_online_init()");
2802 //_____________________________________________________________________________
2803 void AliReconstruction::RunAliEVE()
2805 // Runs AliEVE visualisation of
2806 // the current event.
2807 // Should be executed only after
2808 // successful initialization of AliEVE.
2810 AliInfo("Running AliEVE...");
2811 gROOT->ProcessLine(Form("gAliEveEvent->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2812 gROOT->ProcessLine("gAliEveEvent->StartStopAutoLoadTimer();");
2816 //_____________________________________________________________________________
2817 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2819 // Allows to run QA for a selected set of detectors
2820 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2821 // all selected detectors run the same selected tasks
2823 if (!detAndAction.Contains(":")) {
2824 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2828 Int_t colon = detAndAction.Index(":") ;
2829 fQADetectors = detAndAction(0, colon) ;
2830 if (fQADetectors.Contains("ALL") )
2831 fQADetectors = fFillESD ;
2832 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2833 if (fQATasks.Contains("ALL") ) {
2834 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2836 fQATasks.ToUpper() ;
2838 if ( fQATasks.Contains("RAW") )
2839 tempo = Form("%d ", AliQA::kRAWS) ;
2840 if ( fQATasks.Contains("RECPOINT") )
2841 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2842 if ( fQATasks.Contains("ESD") )
2843 tempo += Form("%d ", AliQA::kESDS) ;
2845 if (fQATasks.IsNull()) {
2846 AliInfo("No QA requested\n") ;
2851 TString tempo(fQATasks) ;
2852 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2853 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2854 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2855 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;