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 fRunVertexFinder(kTRUE),
206 fRunVertexFinderTracks(kTRUE),
207 fRunHLTTracking(kFALSE),
208 fRunMuonTracking(kFALSE),
210 fRunCascadeFinder(kTRUE),
211 fStopOnError(kFALSE),
212 fWriteAlignmentData(kFALSE),
213 fWriteESDfriend(kFALSE),
215 fFillTriggerESD(kTRUE),
223 fRunLocalReconstruction("ALL"),
226 fUseTrackingErrorsForAlignment(""),
227 fGAliceFileName(gAliceFilename),
232 fNumberOfEventsPerFile(1),
235 fLoadAlignFromCDB(kTRUE),
236 fLoadAlignData("ALL"),
242 fParentRawReader(NULL),
245 fDiamondProfile(NULL),
246 fDiamondProfileTPC(NULL),
247 fMeanVertexConstraint(kTRUE),
251 fAlignObjArray(NULL),
254 fInitCDBCalled(kFALSE),
255 fSetRunNumberFromDataCalled(kFALSE),
261 fSameQACycle(kFALSE),
263 fRunPlaneEff(kFALSE),
275 fIsNewRunLoader(kFALSE),
278 // create reconstruction object with default parameters
280 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
281 fReconstructor[iDet] = NULL;
282 fLoader[iDet] = NULL;
283 fTracker[iDet] = NULL;
284 fQADataMaker[iDet] = NULL;
285 fQACycles[iDet] = 999999;
287 fQADataMaker[fgkNDetectors]=NULL; //Global QA
291 //_____________________________________________________________________________
292 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
295 fUniformField(rec.fUniformField),
296 fRunVertexFinder(rec.fRunVertexFinder),
297 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
298 fRunHLTTracking(rec.fRunHLTTracking),
299 fRunMuonTracking(rec.fRunMuonTracking),
300 fRunV0Finder(rec.fRunV0Finder),
301 fRunCascadeFinder(rec.fRunCascadeFinder),
302 fStopOnError(rec.fStopOnError),
303 fWriteAlignmentData(rec.fWriteAlignmentData),
304 fWriteESDfriend(rec.fWriteESDfriend),
305 fWriteAOD(rec.fWriteAOD),
306 fFillTriggerESD(rec.fFillTriggerESD),
308 fCleanESD(rec.fCleanESD),
309 fV0DCAmax(rec.fV0DCAmax),
310 fV0CsPmin(rec.fV0CsPmin),
314 fRunLocalReconstruction(rec.fRunLocalReconstruction),
315 fRunTracking(rec.fRunTracking),
316 fFillESD(rec.fFillESD),
317 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
318 fGAliceFileName(rec.fGAliceFileName),
320 fEquipIdMap(rec.fEquipIdMap),
321 fFirstEvent(rec.fFirstEvent),
322 fLastEvent(rec.fLastEvent),
323 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
326 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
327 fLoadAlignData(rec.fLoadAlignData),
328 fESDPar(rec.fESDPar),
329 fUseHLTData(rec.fUseHLTData),
333 fParentRawReader(NULL),
336 fDiamondProfile(NULL),
337 fDiamondProfileTPC(NULL),
338 fMeanVertexConstraint(rec.fMeanVertexConstraint),
342 fAlignObjArray(rec.fAlignObjArray),
343 fCDBUri(rec.fCDBUri),
345 fInitCDBCalled(rec.fInitCDBCalled),
346 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
347 fQADetectors(rec.fQADetectors),
348 fQATasks(rec.fQATasks),
350 fRunGlobalQA(rec.fRunGlobalQA),
351 fInLoopQA(rec.fInLoopQA),
352 fSameQACycle(rec.fSameQACycle),
353 fRunPlaneEff(rec.fRunPlaneEff),
365 fIsNewRunLoader(rec.fIsNewRunLoader),
370 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
371 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
373 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
374 fReconstructor[iDet] = NULL;
375 fLoader[iDet] = NULL;
376 fTracker[iDet] = NULL;
377 fQADataMaker[iDet] = NULL;
378 fQACycles[iDet] = rec.fQACycles[iDet];
380 fQADataMaker[fgkNDetectors]=NULL; //Global QA
381 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
382 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
386 //_____________________________________________________________________________
387 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
389 // assignment operator
391 this->~AliReconstruction();
392 new(this) AliReconstruction(rec);
396 //_____________________________________________________________________________
397 AliReconstruction::~AliReconstruction()
403 fSpecCDBUri.Delete();
405 AliCodeTimer::Instance()->Print();
408 //_____________________________________________________________________________
409 void AliReconstruction::InitCDB()
411 // activate a default CDB storage
412 // First check if we have any CDB storage set, because it is used
413 // to retrieve the calibration and alignment constants
415 if (fInitCDBCalled) return;
416 fInitCDBCalled = kTRUE;
418 AliCDBManager* man = AliCDBManager::Instance();
419 if (man->IsDefaultStorageSet())
421 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
422 AliWarning("Default CDB storage has been already set !");
423 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
424 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
425 fCDBUri = man->GetDefaultStorage()->GetURI();
428 if (fCDBUri.Length() > 0)
430 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
431 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
432 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
434 fCDBUri="local://$ALICE_ROOT";
435 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
436 AliWarning("Default CDB storage not yet set !!!!");
437 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
438 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
441 man->SetDefaultStorage(fCDBUri);
444 // Now activate the detector specific CDB storage locations
445 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
446 TObject* obj = fSpecCDBUri[i];
448 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
449 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
450 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
451 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
456 //_____________________________________________________________________________
457 void AliReconstruction::SetDefaultStorage(const char* uri) {
458 // Store the desired default CDB storage location
459 // Activate it later within the Run() method
465 //_____________________________________________________________________________
466 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
467 // Store a detector-specific CDB storage location
468 // Activate it later within the Run() method
470 AliCDBPath aPath(calibType);
471 if(!aPath.IsValid()){
472 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
473 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
474 if(!strcmp(calibType, fgkDetectorName[iDet])) {
475 aPath.SetPath(Form("%s/*", calibType));
476 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
480 if(!aPath.IsValid()){
481 AliError(Form("Not a valid path or detector: %s", calibType));
486 // // check that calibType refers to a "valid" detector name
487 // Bool_t isDetector = kFALSE;
488 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
489 // TString detName = fgkDetectorName[iDet];
490 // if(aPath.GetLevel0() == detName) {
491 // isDetector = kTRUE;
497 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
501 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
502 if (obj) fSpecCDBUri.Remove(obj);
503 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
507 //_____________________________________________________________________________
508 Bool_t AliReconstruction::SetRunNumberFromData()
510 // The method is called in Run() in order
511 // to set a correct run number.
512 // In case of raw data reconstruction the
513 // run number is taken from the raw data header
515 if (fSetRunNumberFromDataCalled) return kTRUE;
516 fSetRunNumberFromDataCalled = kTRUE;
518 AliCDBManager* man = AliCDBManager::Instance();
520 if(man->GetRun() > 0) {
521 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
525 AliError("No run loader is found !");
528 // read run number from gAlice
529 if(fRunLoader->GetAliRun())
530 AliCDBManager::Instance()->SetRun(fRunLoader->GetHeader()->GetRun());
533 if(fRawReader->NextEvent()) {
534 AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
535 fRawReader->RewindEvents();
538 if(man->GetRun() > 0) {
539 AliWarning("No raw events is found ! Using settings in AliCDBManager !");
544 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
550 AliError("Neither gAlice nor RawReader objects are found !");
560 //_____________________________________________________________________________
561 void AliReconstruction::SetCDBLock() {
562 // Set CDB lock: from now on it is forbidden to reset the run number
563 // or the default storage or to activate any further storage!
565 AliCDBManager::Instance()->SetLock(1);
568 //_____________________________________________________________________________
569 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
571 // Read the alignment objects from CDB.
572 // Each detector is supposed to have the
573 // alignment objects in DET/Align/Data CDB path.
574 // All the detector objects are then collected,
575 // sorted by geometry level (starting from ALIC) and
576 // then applied to the TGeo geometry.
577 // Finally an overlaps check is performed.
579 // Load alignment data from CDB and fill fAlignObjArray
580 if(fLoadAlignFromCDB){
582 TString detStr = detectors;
583 TString loadAlObjsListOfDets = "";
585 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
586 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
587 loadAlObjsListOfDets += fgkDetectorName[iDet];
588 loadAlObjsListOfDets += " ";
589 } // end loop over detectors
590 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
591 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
593 // Check if the array with alignment objects was
594 // provided by the user. If yes, apply the objects
595 // to the present TGeo geometry
596 if (fAlignObjArray) {
597 if (gGeoManager && gGeoManager->IsClosed()) {
598 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
599 AliError("The misalignment of one or more volumes failed!"
600 "Compare the list of simulated detectors and the list of detector alignment data!");
605 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
611 delete fAlignObjArray; fAlignObjArray=0;
616 //_____________________________________________________________________________
617 void AliReconstruction::SetGAliceFile(const char* fileName)
619 // set the name of the galice file
621 fGAliceFileName = fileName;
624 //_____________________________________________________________________________
625 void AliReconstruction::SetInput(const char* input)
627 // In case the input string starts with 'mem://', we run in an online mode
628 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
629 // file is assumed. One can give as an input:
630 // mem://: - events taken from DAQ monitoring libs online
632 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
636 //_____________________________________________________________________________
637 void AliReconstruction::SetOption(const char* detector, const char* option)
639 // set options for the reconstruction of a detector
641 TObject* obj = fOptions.FindObject(detector);
642 if (obj) fOptions.Remove(obj);
643 fOptions.Add(new TNamed(detector, option));
646 //_____________________________________________________________________________
647 Bool_t AliReconstruction::Run(const char* input)
650 AliCodeTimerAuto("");
652 if (!InitRun(input)) return kFALSE;
653 //******* The loop over events
655 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
656 (fRawReader && fRawReader->NextEvent())) {
657 if (!RunEvent(iEvent)) return kFALSE;
661 if (!FinishRun()) return kFALSE;
666 //_____________________________________________________________________________
667 Bool_t AliReconstruction::InitRun(const char* input)
669 // Initialize all the stuff before
670 // going into the event loop
671 // If the second argument is given, the first one is ignored and
672 // the reconstruction works in an online mode
673 AliCodeTimerAuto("");
675 // Overwrite the previous setting
676 if (input) fInput = input;
678 // set the input in case of raw data
679 fRawReader = AliRawReader::Create(fInput.Data());
681 AliInfo("Reconstruction will run over digits");
683 if (!fEquipIdMap.IsNull() && fRawReader)
684 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
686 if (!fUseHLTData.IsNull()) {
687 // create the RawReaderHLT which performs redirection of HLT input data for
688 // the specified detectors
689 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
691 fParentRawReader=fRawReader;
692 fRawReader=pRawReader;
694 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
698 AliSysInfo::AddStamp("Start");
699 // get the run loader
700 if (!InitRunLoader()) return kFALSE;
701 AliSysInfo::AddStamp("LoadLoader");
703 // Initialize the CDB storage
706 AliSysInfo::AddStamp("LoadCDB");
708 // Set run number in CDBManager (if it is not already set by the user)
709 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
711 // Set CDB lock: from now on it is forbidden to reset the run number
712 // or the default storage or to activate any further storage!
715 // Import ideal TGeo geometry and apply misalignment
717 TString geom(gSystem->DirName(fGAliceFileName));
718 geom += "/geometry.root";
719 AliGeomManager::LoadGeometry(geom.Data());
721 TString detsToCheck=fRunLocalReconstruction;
722 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data()))
723 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
724 if (!gGeoManager) if (fStopOnError) return kFALSE;
727 if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
728 AliSysInfo::AddStamp("LoadGeom");
731 // Get the GRP CDB entry
732 AliCDBEntry* entryGRP = AliCDBManager::Instance()->Get("GRP/GRP/Data");
735 fGRPData = dynamic_cast<TMap*>(entryGRP->GetObject());
738 AliError("No GRP entry found in OCDB!");
743 // Magnetic field map
744 if (!AliTracker::GetFieldMap()) {
745 // Construct the field map out of the information retrieved from GRP.
747 // For the moment, this is a dummy piece of code.
748 // The actual map is expected to be already created in rec.C !
752 Int_t map=AliMagWrapCheb::k5kG;
753 Bool_t dipoleON=kTRUE;
756 TObjString *l3Current=
757 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
759 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
762 TObjString *l3Polarity=
763 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
765 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
770 TObjString *diCurrent=
771 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
773 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
776 TObjString *diPolarity=
777 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
779 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
785 new AliMagWrapCheb("Maps","Maps",2,factor,10.,map,dipoleON);
786 AliTracker::SetFieldMap(field,fUniformField);
789 AliFatal("Please, provide the field map ! Crashing deliberately...");
794 // Get the diamond profile from OCDB
795 AliCDBEntry* entry = AliCDBManager::Instance()
796 ->Get("GRP/Calib/MeanVertex");
799 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
801 AliError("No diamond profile found in OCDB!");
805 entry = AliCDBManager::Instance()
806 ->Get("GRP/Calib/MeanVertexTPC");
809 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
811 AliError("No diamond profile found in OCDB!");
814 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
815 if(fDiamondProfile && fMeanVertexConstraint) ftVertexer->SetVtxStart(fDiamondProfile);
818 if (fRunVertexFinder && !CreateVertexer()) {
824 AliSysInfo::AddStamp("Vertexer");
827 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
833 AliSysInfo::AddStamp("LoadTrackers");
835 // get the possibly already existing ESD file and tree
836 fesd = new AliESDEvent(); fhltesd = new AliESDEvent();
837 if (!gSystem->AccessPathName("AliESDs.root")){
838 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
839 ffileOld = TFile::Open("AliESDs.old.root");
840 if (ffileOld && ffileOld->IsOpen()) {
841 ftreeOld = (TTree*) ffileOld->Get("esdTree");
842 if (ftreeOld)fesd->ReadFromTree(ftreeOld);
843 fhlttreeOld = (TTree*) ffileOld->Get("HLTesdTree");
844 if (fhlttreeOld) fhltesd->ReadFromTree(fhlttreeOld);
848 // create the ESD output file and tree
849 ffile = TFile::Open("AliESDs.root", "RECREATE");
850 ffile->SetCompressionLevel(2);
851 if (!ffile->IsOpen()) {
852 AliError("opening AliESDs.root failed");
853 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
856 ftree = new TTree("esdTree", "Tree with ESD objects");
857 fesd = new AliESDEvent();
858 fesd->CreateStdContent();
859 fesd->WriteToTree(ftree);
861 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
862 fhltesd = new AliESDEvent();
863 fhltesd->CreateStdContent();
864 fhltesd->WriteToTree(fhlttree);
867 if (fWriteESDfriend) {
868 fesdf = new AliESDfriend();
869 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
870 br->SetFile("AliESDfriends.root");
871 fesd->AddObject(fesdf);
876 if (fRawReader) fRawReader->RewindEvents();
879 gSystem->GetProcInfo(&ProcInfo);
880 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
883 if (fRunQA && fRawReader && fQATasks.Contains(Form("%d", AliQA::kRAWS))) {
884 AliQADataMakerSteer qas ;
885 qas.Run(fRunLocalReconstruction, fRawReader) ;
886 fSameQACycle = kTRUE ;
888 //Initialize the QA and start of cycle for out-of-cycle QA
890 TString detStr(fQADetectors) ;
891 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
892 if (!IsSelected(fgkDetectorName[iDet], detStr))
894 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
897 AliInfo(Form("Initializing the QA data maker for %s",
898 fgkDetectorName[iDet]));
899 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
900 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
901 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
902 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
904 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
905 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
906 fSameQACycle = kTRUE;
908 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
909 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle);
910 fSameQACycle = kTRUE;
915 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
916 AliInfo(Form("Initializing the global QA data maker"));
917 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
919 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
920 AliTracker::SetResidualsArray(arr);
922 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
923 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
926 fSameQACycle = kFALSE;
927 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
928 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
929 fSameQACycle = kTRUE;
931 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
932 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle);
933 fSameQACycle = kTRUE;
939 //Initialize the Plane Efficiency framework
940 if (fRunPlaneEff && !InitPlaneEff()) {
941 if(fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
944 if (strcmp(gProgName,"alieve") == 0)
945 fRunAliEVE = InitAliEVE();
950 //_____________________________________________________________________________
951 Bool_t AliReconstruction::RunEvent(Int_t iEvent)
953 // run the reconstruction over a single event
954 // The event loop is steered in Run method
956 AliCodeTimerAuto("");
958 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
959 fRunLoader->SetEventNumber(iEvent);
960 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
962 //?? fRunLoader->MakeTree("H");
963 fRunLoader->TreeE()->Fill();
966 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
967 // copy old ESD to the new one
969 fesd->ReadFromTree(ftreeOld);
970 ftreeOld->GetEntry(iEvent);
974 fhltesd->ReadFromTree(fhlttreeOld);
975 fhlttreeOld->GetEntry(iEvent);
981 AliInfo(Form("processing event %d", iEvent));
983 //Start of cycle for the in-loop QA
986 fSameQACycle = kFALSE ;
987 TString detStr(fQADetectors);
988 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
989 if (!IsSelected(fgkDetectorName[iDet], detStr))
991 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
994 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
995 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
996 fSameQACycle = kTRUE;
998 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
999 qadm->StartOfCycle(AliQA::kESDS, fSameQACycle) ;
1000 fSameQACycle = kTRUE;
1004 fSameQACycle = kFALSE;
1005 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
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 fRunLoader->GetEvent(iEvent);
1020 char aFileName[256];
1021 sprintf(aFileName, "ESD_%d.%d_final.root",
1022 fRunLoader->GetHeader()->GetRun(),
1023 fRunLoader->GetHeader()->GetEventNrInRun());
1024 if (!gSystem->AccessPathName(aFileName)) return kTRUE;
1026 // local single event reconstruction
1027 if (!fRunLocalReconstruction.IsNull()) {
1028 TString detectors=fRunLocalReconstruction;
1029 // run HLT event reconstruction first
1030 // ;-( IsSelected changes the string
1031 if (IsSelected("HLT", detectors) &&
1032 !RunLocalEventReconstruction("HLT")) {
1033 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1035 detectors=fRunLocalReconstruction;
1036 detectors.ReplaceAll("HLT", "");
1037 if (!RunLocalEventReconstruction(detectors)) {
1038 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1042 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1043 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1044 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1045 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1047 // Set magnetic field from the tracker
1048 fesd->SetMagneticField(AliTracker::GetBz());
1049 fhltesd->SetMagneticField(AliTracker::GetBz());
1053 // Fill raw-data error log into the ESD
1054 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1057 if (fRunVertexFinder) {
1058 if (!ReadESD(fesd, "vertex")) {
1059 if (!RunVertexFinder(fesd)) {
1060 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1062 if (fCheckPointLevel > 0) WriteESD(fesd, "vertex");
1067 if (!fRunTracking.IsNull()) {
1068 if (fRunMuonTracking) {
1069 if (!RunMuonTracking(fesd)) {
1070 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1076 if (!fRunTracking.IsNull()) {
1077 if (!ReadESD(fesd, "tracking")) {
1078 if (!RunTracking(fesd)) {
1079 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1081 if (fCheckPointLevel > 0) WriteESD(fesd, "tracking");
1086 if (!fFillESD.IsNull()) {
1087 TString detectors=fFillESD;
1088 // run HLT first and on hltesd
1089 // ;-( IsSelected changes the string
1090 if (IsSelected("HLT", detectors) &&
1091 !FillESD(fhltesd, "HLT")) {
1092 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1095 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1096 if (detectors.Contains("ALL")) {
1098 for (Int_t idet=0; idet<fgkNDetectors; ++idet){
1099 detectors += fgkDetectorName[idet];
1103 detectors.ReplaceAll("HLT", "");
1104 if (!FillESD(fesd, detectors)) {
1105 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1109 // fill Event header information from the RawEventHeader
1110 if (fRawReader){FillRawEventHeaderESD(fesd);}
1113 AliESDpid::MakePID(fesd);
1114 if (fCheckPointLevel > 1) WriteESD(fesd, "PID");
1116 if (fFillTriggerESD) {
1117 if (!ReadESD(fesd, "trigger")) {
1118 if (!FillTriggerESD(fesd)) {
1119 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1121 if (fCheckPointLevel > 1) WriteESD(fesd, "trigger");
1128 // Propagate track to the beam pipe (if not already done by ITS)
1130 const Int_t ntracks = fesd->GetNumberOfTracks();
1131 const Double_t kBz = fesd->GetMagneticField();
1132 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1135 UShort_t *selectedIdx=new UShort_t[ntracks];
1137 for (Int_t itrack=0; itrack<ntracks; itrack++){
1138 const Double_t kMaxStep = 5; //max step over the material
1141 AliESDtrack *track = fesd->GetTrack(itrack);
1142 if (!track) continue;
1144 AliExternalTrackParam *tpcTrack =
1145 (AliExternalTrackParam *)track->GetTPCInnerParam();
1149 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
1154 Int_t n=trkArray.GetEntriesFast();
1155 selectedIdx[n]=track->GetID();
1156 trkArray.AddLast(tpcTrack);
1159 if (track->GetX() < kRadius) continue;
1162 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1164 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kRadius);
1169 // Improve the reconstructed primary vertex position using the tracks
1171 TObject *obj = fOptions.FindObject("ITS");
1173 TString optITS = obj->GetTitle();
1174 if (optITS.Contains("cosmics") || optITS.Contains("COSMICS"))
1175 fRunVertexFinderTracks=kFALSE;
1177 if (fRunVertexFinderTracks) {
1178 // TPC + ITS primary vertex
1179 ftVertexer->SetITSrefitRequired();
1180 if(fDiamondProfile && fMeanVertexConstraint) {
1181 ftVertexer->SetVtxStart(fDiamondProfile);
1183 ftVertexer->SetConstraintOff();
1185 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1187 if (pvtx->GetStatus()) {
1188 fesd->SetPrimaryVertex(pvtx);
1189 for (Int_t i=0; i<ntracks; i++) {
1190 AliESDtrack *t = fesd->GetTrack(i);
1191 t->RelateToVertex(pvtx, kBz, kRadius);
1196 // TPC-only primary vertex
1197 ftVertexer->SetITSrefitNotRequired();
1198 if(fDiamondProfileTPC && fMeanVertexConstraint) {
1199 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1201 ftVertexer->SetConstraintOff();
1203 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1205 if (pvtx->GetStatus()) {
1206 fesd->SetPrimaryVertexTPC(pvtx);
1207 Int_t nsel=trkArray.GetEntriesFast();
1208 for (Int_t i=0; i<nsel; i++) {
1209 AliExternalTrackParam *t =
1210 (AliExternalTrackParam *)trkArray.UncheckedAt(i);
1211 t->PropagateToDCA(pvtx, kBz, kRadius);
1217 delete[] selectedIdx;
1219 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1224 AliV0vertexer vtxer;
1225 vtxer.Tracks2V0vertices(fesd);
1227 if (fRunCascadeFinder) {
1229 AliCascadeVertexer cvtxer;
1230 cvtxer.V0sTracks2CascadeVertices(fesd);
1235 if (fCleanESD) CleanESD(fesd);
1239 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1240 if (qadm && fQATasks.Contains(Form("%d", AliQA::kESDS)))
1241 qadm->Exec(AliQA::kESDS, fesd);
1245 if (fWriteESDfriend) {
1246 fesdf->~AliESDfriend();
1247 new (fesdf) AliESDfriend(); // Reset...
1248 fesd->GetESDfriend(fesdf);
1256 if (fRunAliEVE) RunAliEVE();
1258 if (fCheckPointLevel > 0) WriteESD(fesd, "final");
1261 if (fWriteESDfriend) {
1262 fesdf->~AliESDfriend();
1263 new (fesdf) AliESDfriend(); // Reset...
1266 ProcInfo_t ProcInfo;
1267 gSystem->GetProcInfo(&ProcInfo);
1268 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1271 // End of cycle for the in-loop
1275 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1276 if (!IsSelected(fgkDetectorName[iDet], fQADetectors))
1278 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1281 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1282 qadm->EndOfCycle(AliQA::kRECPOINTS);
1283 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1284 qadm->EndOfCycle(AliQA::kESDS);
1289 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1291 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1292 qadm->EndOfCycle(AliQA::kRECPOINTS);
1293 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1294 qadm->EndOfCycle(AliQA::kESDS);
1303 //_____________________________________________________________________________
1304 Bool_t AliReconstruction::FinishRun()
1307 // Called after the exit
1308 // from the event loop
1309 AliCodeTimerAuto("");
1311 if (fIsNewRunLoader) { // galice.root didn't exist
1312 fRunLoader->WriteHeader("OVERWRITE");
1313 fRunLoader->CdGAFile();
1314 fRunLoader->Write(0, TObject::kOverwrite);
1317 ftree->GetUserInfo()->Add(fesd);
1318 fhlttree->GetUserInfo()->Add(fhltesd);
1320 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1321 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1323 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1324 cdbMapCopy->SetOwner(1);
1325 cdbMapCopy->SetName("cdbMap");
1326 TIter iter(cdbMap->GetTable());
1329 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1330 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1331 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1332 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1335 TList *cdbListCopy = new TList();
1336 cdbListCopy->SetOwner(1);
1337 cdbListCopy->SetName("cdbList");
1339 TIter iter2(cdbList);
1342 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1343 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1346 ftree->GetUserInfo()->Add(cdbMapCopy);
1347 ftree->GetUserInfo()->Add(cdbListCopy);
1350 if(fESDPar.Contains("ESD.par")){
1351 AliInfo("Attaching ESD.par to Tree");
1352 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1353 ftree->GetUserInfo()->Add(fn);
1359 if (fWriteESDfriend)
1360 ftree->SetBranchStatus("ESDfriend*",0);
1361 // we want to have only one tree version number
1362 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1365 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1366 if (fRunPlaneEff && !FinishPlaneEff()) {
1367 AliWarning("Finish PlaneEff evaluation failed");
1371 CleanUp(ffile, ffileOld);
1374 AliWarning("AOD creation not supported anymore during reconstruction. See ANALYSIS/AliAnalysisTaskESDfilter.cxx instead.");
1377 // Create tags for the events in the ESD tree (the ESD tree is always present)
1378 // In case of empty events the tags will contain dummy values
1379 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1380 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1382 AliWarning("AOD tag creation not supported anymore during reconstruction.");
1385 //Finish QA and end of cycle for out-of-loop QA
1388 AliQADataMakerSteer qas;
1389 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1390 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS, fSameQACycle);
1392 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1393 qas.Run(fRunLocalReconstruction.Data(), AliQA::kESDS, fSameQACycle);
1395 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1397 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1398 qadm->EndOfCycle(AliQA::kRECPOINTS);
1399 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1400 qadm->EndOfCycle(AliQA::kESDS);
1407 // Cleanup of CDB manager: cache and active storages!
1408 AliCDBManager::Instance()->ClearCache();
1414 //_____________________________________________________________________________
1415 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1417 // run the local reconstruction
1418 static Int_t eventNr=0;
1419 AliCodeTimerAuto("")
1421 // AliCDBManager* man = AliCDBManager::Instance();
1422 // Bool_t origCache = man->GetCacheFlag();
1424 // TString detStr = detectors;
1425 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1426 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1427 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1428 // if (!reconstructor) continue;
1429 // if (reconstructor->HasLocalReconstruction()) continue;
1431 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1432 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1434 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1435 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1437 // man->SetCacheFlag(kTRUE);
1438 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1439 // man->GetAll(calibPath); // entries are cached!
1441 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1443 // if (fRawReader) {
1444 // fRawReader->RewindEvents();
1445 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1447 // reconstructor->Reconstruct(fRunLoader);
1450 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1451 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1453 // // unload calibration data
1454 // man->UnloadFromCache(calibPath);
1455 // //man->ClearCache();
1458 // man->SetCacheFlag(origCache);
1460 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1461 // AliError(Form("the following detectors were not found: %s",
1463 // if (fStopOnError) return kFALSE;
1470 //_____________________________________________________________________________
1471 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1473 // run the local reconstruction
1475 static Int_t eventNr=0;
1476 AliCodeTimerAuto("")
1478 TString detStr = detectors;
1479 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1480 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1481 AliReconstructor* reconstructor = GetReconstructor(iDet);
1482 if (!reconstructor) continue;
1483 AliLoader* loader = fLoader[iDet];
1484 // Matthias April 2008: temporary fix to run HLT reconstruction
1485 // although the HLT loader is missing
1486 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1488 reconstructor->Reconstruct(fRawReader, NULL);
1491 reconstructor->Reconstruct(dummy, NULL);
1496 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1499 // conversion of digits
1500 if (fRawReader && reconstructor->HasDigitConversion()) {
1501 AliInfo(Form("converting raw data digits into root objects for %s",
1502 fgkDetectorName[iDet]));
1503 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1504 fgkDetectorName[iDet]));
1505 loader->LoadDigits("update");
1506 loader->CleanDigits();
1507 loader->MakeDigitsContainer();
1508 TTree* digitsTree = loader->TreeD();
1509 reconstructor->ConvertDigits(fRawReader, digitsTree);
1510 loader->WriteDigits("OVERWRITE");
1511 loader->UnloadDigits();
1513 // local reconstruction
1514 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1515 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1516 loader->LoadRecPoints("update");
1517 loader->CleanRecPoints();
1518 loader->MakeRecPointsContainer();
1519 TTree* clustersTree = loader->TreeR();
1520 if (fRawReader && !reconstructor->HasDigitConversion()) {
1521 reconstructor->Reconstruct(fRawReader, clustersTree);
1523 loader->LoadDigits("read");
1524 TTree* digitsTree = loader->TreeD();
1526 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1527 if (fStopOnError) return kFALSE;
1529 reconstructor->Reconstruct(digitsTree, clustersTree);
1531 loader->UnloadDigits();
1534 // In-loop QA for local reconstrucion
1535 if (fRunQA && fInLoopQA) {
1536 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1539 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1541 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1543 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1544 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1546 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1550 loader->WriteRecPoints("OVERWRITE");
1551 loader->UnloadRecPoints();
1552 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1555 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1556 AliError(Form("the following detectors were not found: %s",
1558 if (fStopOnError) return kFALSE;
1564 //_____________________________________________________________________________
1565 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1567 // run the barrel tracking
1569 AliCodeTimerAuto("")
1571 AliESDVertex* vertex = NULL;
1572 Double_t vtxPos[3] = {0, 0, 0};
1573 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1574 TArrayF mcVertex(3);
1575 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1576 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1577 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1581 AliInfo("running the ITS vertex finder");
1583 fLoader[0]->LoadRecPoints();
1584 TTree* cltree = fLoader[0]->TreeR();
1586 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1587 vertex = fVertexer->FindVertexForCurrentEvent(cltree);
1590 AliError("Can't get the ITS cluster tree");
1592 fLoader[0]->UnloadRecPoints();
1595 AliError("Can't get the ITS loader");
1598 AliWarning("Vertex not found");
1599 vertex = new AliESDVertex();
1600 vertex->SetName("default");
1603 vertex->SetName("reconstructed");
1607 AliInfo("getting the primary vertex from MC");
1608 vertex = new AliESDVertex(vtxPos, vtxErr);
1612 vertex->GetXYZ(vtxPos);
1613 vertex->GetSigmaXYZ(vtxErr);
1615 AliWarning("no vertex reconstructed");
1616 vertex = new AliESDVertex(vtxPos, vtxErr);
1618 esd->SetPrimaryVertexSPD(vertex);
1619 // if SPD multiplicity has been determined, it is stored in the ESD
1620 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1621 if(mult)esd->SetMultiplicity(mult);
1623 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1624 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1631 //_____________________________________________________________________________
1632 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1634 // run the HLT barrel tracking
1636 AliCodeTimerAuto("")
1639 AliError("Missing runLoader!");
1643 AliInfo("running HLT tracking");
1645 // Get a pointer to the HLT reconstructor
1646 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1647 if (!reconstructor) return kFALSE;
1650 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1651 TString detName = fgkDetectorName[iDet];
1652 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1653 reconstructor->SetOption(detName.Data());
1654 AliTracker *tracker = reconstructor->CreateTracker();
1656 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1657 if (fStopOnError) return kFALSE;
1661 Double_t vtxErr[3]={0.005,0.005,0.010};
1662 const AliESDVertex *vertex = esd->GetVertex();
1663 vertex->GetXYZ(vtxPos);
1664 tracker->SetVertex(vtxPos,vtxErr);
1666 fLoader[iDet]->LoadRecPoints("read");
1667 TTree* tree = fLoader[iDet]->TreeR();
1669 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1672 tracker->LoadClusters(tree);
1674 if (tracker->Clusters2Tracks(esd) != 0) {
1675 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1679 tracker->UnloadClusters();
1687 //_____________________________________________________________________________
1688 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1690 // run the muon spectrometer tracking
1692 AliCodeTimerAuto("")
1695 AliError("Missing runLoader!");
1698 Int_t iDet = 7; // for MUON
1700 AliInfo("is running...");
1702 // Get a pointer to the MUON reconstructor
1703 AliReconstructor *reconstructor = GetReconstructor(iDet);
1704 if (!reconstructor) return kFALSE;
1707 TString detName = fgkDetectorName[iDet];
1708 AliDebug(1, Form("%s tracking", detName.Data()));
1709 AliTracker *tracker = reconstructor->CreateTracker();
1711 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1716 fLoader[iDet]->LoadRecPoints("read");
1718 tracker->LoadClusters(fLoader[iDet]->TreeR());
1720 Int_t rv = tracker->Clusters2Tracks(esd);
1724 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1728 fLoader[iDet]->UnloadRecPoints();
1730 tracker->UnloadClusters();
1738 //_____________________________________________________________________________
1739 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1741 // run the barrel tracking
1742 static Int_t eventNr=0;
1743 AliCodeTimerAuto("")
1745 AliInfo("running tracking");
1747 //Fill the ESD with the T0 info (will be used by the TOF)
1748 if (fReconstructor[11] && fLoader[11]) {
1749 fLoader[11]->LoadRecPoints("READ");
1750 TTree *treeR = fLoader[11]->TreeR();
1751 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1754 // pass 1: TPC + ITS inwards
1755 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1756 if (!fTracker[iDet]) continue;
1757 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1760 fLoader[iDet]->LoadRecPoints("read");
1761 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1762 TTree* tree = fLoader[iDet]->TreeR();
1764 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1767 fTracker[iDet]->LoadClusters(tree);
1768 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1770 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1771 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1774 if (fCheckPointLevel > 1) {
1775 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1777 // preliminary PID in TPC needed by the ITS tracker
1779 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1780 AliESDpid::MakePID(esd);
1782 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1785 // pass 2: ALL backwards
1787 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1788 if (!fTracker[iDet]) continue;
1789 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1792 if (iDet > 1) { // all except ITS, TPC
1794 fLoader[iDet]->LoadRecPoints("read");
1795 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1796 tree = fLoader[iDet]->TreeR();
1798 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1801 fTracker[iDet]->LoadClusters(tree);
1802 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1806 if (iDet>1) // start filling residuals for the "outer" detectors
1807 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1809 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1810 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1813 if (fCheckPointLevel > 1) {
1814 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1818 if (iDet > 2) { // all except ITS, TPC, TRD
1819 fTracker[iDet]->UnloadClusters();
1820 fLoader[iDet]->UnloadRecPoints();
1822 // updated PID in TPC needed by the ITS tracker -MI
1824 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1825 AliESDpid::MakePID(esd);
1827 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1829 //stop filling residuals for the "outer" detectors
1830 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1832 // write space-points to the ESD in case alignment data output
1834 if (fWriteAlignmentData)
1835 WriteAlignmentData(esd);
1837 // pass 3: TRD + TPC + ITS refit inwards
1839 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1840 if (!fTracker[iDet]) continue;
1841 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1844 if (iDet<2) // start filling residuals for TPC and ITS
1845 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1847 if (fTracker[iDet]->RefitInward(esd) != 0) {
1848 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1851 // run postprocessing
1852 if (fTracker[iDet]->PostProcess(esd) != 0) {
1853 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1856 if (fCheckPointLevel > 1) {
1857 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1859 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1861 fTracker[iDet]->UnloadClusters();
1862 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1863 fLoader[iDet]->UnloadRecPoints();
1864 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1866 // stop filling residuals for TPC and ITS
1867 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1873 //_____________________________________________________________________________
1874 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1876 // Remove the data which are not needed for the physics analysis.
1879 Int_t nTracks=esd->GetNumberOfTracks();
1880 Int_t nV0s=esd->GetNumberOfV0s();
1882 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1884 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1885 Bool_t rc=esd->Clean(cleanPars);
1887 nTracks=esd->GetNumberOfTracks();
1888 nV0s=esd->GetNumberOfV0s();
1890 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1895 //_____________________________________________________________________________
1896 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1898 // fill the event summary data
1900 AliCodeTimerAuto("")
1901 static Int_t eventNr=0;
1902 TString detStr = detectors;
1904 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1905 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1906 AliReconstructor* reconstructor = GetReconstructor(iDet);
1907 if (!reconstructor) continue;
1908 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1909 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1910 TTree* clustersTree = NULL;
1911 if (fLoader[iDet]) {
1912 fLoader[iDet]->LoadRecPoints("read");
1913 clustersTree = fLoader[iDet]->TreeR();
1914 if (!clustersTree) {
1915 AliError(Form("Can't get the %s clusters tree",
1916 fgkDetectorName[iDet]));
1917 if (fStopOnError) return kFALSE;
1920 if (fRawReader && !reconstructor->HasDigitConversion()) {
1921 reconstructor->FillESD(fRawReader, clustersTree, esd);
1923 TTree* digitsTree = NULL;
1924 if (fLoader[iDet]) {
1925 fLoader[iDet]->LoadDigits("read");
1926 digitsTree = fLoader[iDet]->TreeD();
1928 AliError(Form("Can't get the %s digits tree",
1929 fgkDetectorName[iDet]));
1930 if (fStopOnError) return kFALSE;
1933 reconstructor->FillESD(digitsTree, clustersTree, esd);
1934 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1936 if (fLoader[iDet]) {
1937 fLoader[iDet]->UnloadRecPoints();
1940 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
1944 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1945 AliError(Form("the following detectors were not found: %s",
1947 if (fStopOnError) return kFALSE;
1949 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
1954 //_____________________________________________________________________________
1955 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
1957 // Reads the trigger decision which is
1958 // stored in Trigger.root file and fills
1959 // the corresponding esd entries
1961 AliCodeTimerAuto("")
1963 AliInfo("Filling trigger information into the ESD");
1965 AliCentralTrigger *aCTP = NULL;
1968 AliCTPRawStream input(fRawReader);
1969 if (!input.Next()) {
1970 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 !");
1971 ULong_t mask = (((ULong_t)fRawReader->GetTriggerPattern()[1]) << 32) +
1972 fRawReader->GetTriggerPattern()[0];
1973 esd->SetTriggerMask(mask);
1974 esd->SetTriggerCluster(0);
1977 esd->SetTriggerMask(input.GetClassMask());
1978 esd->SetTriggerCluster(input.GetClusterMask());
1981 aCTP = new AliCentralTrigger();
1982 TString configstr("");
1983 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
1984 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
1990 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1992 if (!runloader->LoadTrigger()) {
1993 aCTP = runloader->GetTrigger();
1994 esd->SetTriggerMask(aCTP->GetClassMask());
1995 esd->SetTriggerCluster(aCTP->GetClusterMask());
1998 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
2003 AliError("No run loader is available! The trigger information is not stored in the ESD !");
2008 // Now fill the trigger class names into AliESDRun object
2009 AliTriggerConfiguration *config = aCTP->GetConfiguration();
2011 AliError("No trigger configuration has been found! The trigger classes information will not be stored in ESD!");
2012 if (fRawReader) delete aCTP;
2016 const TObjArray& classesArray = config->GetClasses();
2017 Int_t nclasses = classesArray.GetEntriesFast();
2018 for( Int_t j=0; j<nclasses; j++ ) {
2019 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
2020 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
2021 esd->SetTriggerClass(trclass->GetName(),trindex);
2024 if (fRawReader) delete aCTP;
2032 //_____________________________________________________________________________
2033 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2036 // Filling information from RawReader Header
2039 AliInfo("Filling information from RawReader Header");
2040 esd->SetBunchCrossNumber(0);
2041 esd->SetOrbitNumber(0);
2042 esd->SetPeriodNumber(0);
2043 esd->SetTimeStamp(0);
2044 esd->SetEventType(0);
2045 const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
2048 const UInt_t *id = eventHeader->GetP("Id");
2049 esd->SetBunchCrossNumber((id)[1]&0x00000fff);
2050 esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
2051 esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
2053 esd->SetTimeStamp((eventHeader->Get("Timestamp")));
2054 esd->SetEventType((eventHeader->Get("Type")));
2061 //_____________________________________________________________________________
2062 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2064 // check whether detName is contained in detectors
2065 // if yes, it is removed from detectors
2067 // check if all detectors are selected
2068 if ((detectors.CompareTo("ALL") == 0) ||
2069 detectors.BeginsWith("ALL ") ||
2070 detectors.EndsWith(" ALL") ||
2071 detectors.Contains(" ALL ")) {
2076 // search for the given detector
2077 Bool_t result = kFALSE;
2078 if ((detectors.CompareTo(detName) == 0) ||
2079 detectors.BeginsWith(detName+" ") ||
2080 detectors.EndsWith(" "+detName) ||
2081 detectors.Contains(" "+detName+" ")) {
2082 detectors.ReplaceAll(detName, "");
2086 // clean up the detectors string
2087 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2088 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2089 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2094 //_____________________________________________________________________________
2095 Bool_t AliReconstruction::InitRunLoader()
2097 // get or create the run loader
2099 if (gAlice) delete gAlice;
2102 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2103 // load all base libraries to get the loader classes
2104 TString libs = gSystem->GetLibraries();
2105 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2106 TString detName = fgkDetectorName[iDet];
2107 if (detName == "HLT") continue;
2108 if (libs.Contains("lib" + detName + "base.so")) continue;
2109 gSystem->Load("lib" + detName + "base.so");
2111 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2113 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2118 fRunLoader->CdGAFile();
2119 fRunLoader->LoadgAlice();
2121 //PH This is a temporary fix to give access to the kinematics
2122 //PH that is needed for the labels of ITS clusters
2123 fRunLoader->LoadHeader();
2124 fRunLoader->LoadKinematics();
2126 } else { // galice.root does not exist
2128 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2132 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2133 AliConfig::GetDefaultEventFolderName(),
2136 AliError(Form("could not create run loader in file %s",
2137 fGAliceFileName.Data()));
2141 fIsNewRunLoader = kTRUE;
2142 fRunLoader->MakeTree("E");
2144 if (fNumberOfEventsPerFile > 0)
2145 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2147 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2153 //_____________________________________________________________________________
2154 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2156 // get the reconstructor object and the loader for a detector
2158 if (fReconstructor[iDet]) return fReconstructor[iDet];
2160 // load the reconstructor object
2161 TPluginManager* pluginManager = gROOT->GetPluginManager();
2162 TString detName = fgkDetectorName[iDet];
2163 TString recName = "Ali" + detName + "Reconstructor";
2165 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2167 AliReconstructor* reconstructor = NULL;
2168 // first check if a plugin is defined for the reconstructor
2169 TPluginHandler* pluginHandler =
2170 pluginManager->FindHandler("AliReconstructor", detName);
2171 // if not, add a plugin for it
2172 if (!pluginHandler) {
2173 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2174 TString libs = gSystem->GetLibraries();
2175 if (libs.Contains("lib" + detName + "base.so") ||
2176 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2177 pluginManager->AddHandler("AliReconstructor", detName,
2178 recName, detName + "rec", recName + "()");
2180 pluginManager->AddHandler("AliReconstructor", detName,
2181 recName, detName, recName + "()");
2183 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2185 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2186 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2188 if (reconstructor) {
2189 TObject* obj = fOptions.FindObject(detName.Data());
2190 if (obj) reconstructor->SetOption(obj->GetTitle());
2191 reconstructor->Init();
2192 fReconstructor[iDet] = reconstructor;
2195 // get or create the loader
2196 if (detName != "HLT") {
2197 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2198 if (!fLoader[iDet]) {
2199 AliConfig::Instance()
2200 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2202 // first check if a plugin is defined for the loader
2204 pluginManager->FindHandler("AliLoader", detName);
2205 // if not, add a plugin for it
2206 if (!pluginHandler) {
2207 TString loaderName = "Ali" + detName + "Loader";
2208 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2209 pluginManager->AddHandler("AliLoader", detName,
2210 loaderName, detName + "base",
2211 loaderName + "(const char*, TFolder*)");
2212 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2214 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2216 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2217 fRunLoader->GetEventFolder());
2219 if (!fLoader[iDet]) { // use default loader
2220 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2222 if (!fLoader[iDet]) {
2223 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2224 if (fStopOnError) return NULL;
2226 fRunLoader->AddLoader(fLoader[iDet]);
2227 fRunLoader->CdGAFile();
2228 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2229 fRunLoader->Write(0, TObject::kOverwrite);
2234 return reconstructor;
2237 //_____________________________________________________________________________
2238 Bool_t AliReconstruction::CreateVertexer()
2240 // create the vertexer
2243 AliReconstructor* itsReconstructor = GetReconstructor(0);
2244 if (itsReconstructor) {
2245 fVertexer = itsReconstructor->CreateVertexer();
2248 AliWarning("couldn't create a vertexer for ITS");
2249 if (fStopOnError) return kFALSE;
2255 //_____________________________________________________________________________
2256 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2258 // create the trackers
2260 TString detStr = detectors;
2261 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2262 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2263 AliReconstructor* reconstructor = GetReconstructor(iDet);
2264 if (!reconstructor) continue;
2265 TString detName = fgkDetectorName[iDet];
2266 if (detName == "HLT") {
2267 fRunHLTTracking = kTRUE;
2270 if (detName == "MUON") {
2271 fRunMuonTracking = kTRUE;
2276 fTracker[iDet] = reconstructor->CreateTracker();
2277 if (!fTracker[iDet] && (iDet < 7)) {
2278 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2279 if (fStopOnError) return kFALSE;
2281 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2287 //_____________________________________________________________________________
2288 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2290 // delete trackers and the run loader and close and delete the file
2292 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2293 delete fReconstructor[iDet];
2294 fReconstructor[iDet] = NULL;
2295 fLoader[iDet] = NULL;
2296 delete fTracker[iDet];
2297 fTracker[iDet] = NULL;
2298 // delete fQADataMaker[iDet];
2299 // fQADataMaker[iDet] = NULL;
2304 if (ftVertexer) delete ftVertexer;
2307 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2308 delete fDiamondProfile;
2309 fDiamondProfile = NULL;
2310 delete fDiamondProfileTPC;
2311 fDiamondProfileTPC = NULL;
2321 if (fParentRawReader) delete fParentRawReader;
2322 fParentRawReader=NULL;
2332 gSystem->Unlink("AliESDs.old.root");
2337 //_____________________________________________________________________________
2339 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2341 // read the ESD event from a file
2343 if (!esd) return kFALSE;
2345 sprintf(fileName, "ESD_%d.%d_%s.root",
2346 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2347 if (gSystem->AccessPathName(fileName)) return kFALSE;
2349 AliInfo(Form("reading ESD from file %s", fileName));
2350 AliDebug(1, Form("reading ESD from file %s", fileName));
2351 TFile* file = TFile::Open(fileName);
2352 if (!file || !file->IsOpen()) {
2353 AliError(Form("opening %s failed", fileName));
2360 esd = (AliESDEvent*) file->Get("ESD");
2369 //_____________________________________________________________________________
2370 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2372 // write the ESD event to a file
2376 sprintf(fileName, "ESD_%d.%d_%s.root",
2377 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2379 AliDebug(1, Form("writing ESD to file %s", fileName));
2380 TFile* file = TFile::Open(fileName, "recreate");
2381 if (!file || !file->IsOpen()) {
2382 AliError(Form("opening %s failed", fileName));
2391 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2393 // Write space-points which are then used in the alignment procedures
2394 // For the moment only ITS, TRD and TPC
2396 // Load TOF clusters
2398 fLoader[3]->LoadRecPoints("read");
2399 TTree* tree = fLoader[3]->TreeR();
2401 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2404 fTracker[3]->LoadClusters(tree);
2406 Int_t ntracks = esd->GetNumberOfTracks();
2407 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2409 AliESDtrack *track = esd->GetTrack(itrack);
2412 for (Int_t iDet = 3; iDet >= 0; iDet--)
2413 nsp += track->GetNcls(iDet);
2415 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2416 track->SetTrackPointArray(sp);
2418 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2419 AliTracker *tracker = fTracker[iDet];
2420 if (!tracker) continue;
2421 Int_t nspdet = track->GetNcls(iDet);
2422 if (nspdet <= 0) continue;
2423 track->GetClusters(iDet,idx);
2427 while (isp2 < nspdet) {
2429 TString dets = fgkDetectorName[iDet];
2430 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2431 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2432 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2433 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2434 isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2436 isvalid = tracker->GetTrackPoint(idx[isp2],p);
2439 const Int_t kNTPCmax = 159;
2440 if (iDet==1 && isp2>kNTPCmax) break; // to be fixed
2441 if (!isvalid) continue;
2442 sp->AddPoint(isptrack,&p); isptrack++; isp++;
2448 fTracker[3]->UnloadClusters();
2449 fLoader[3]->UnloadRecPoints();
2453 //_____________________________________________________________________________
2454 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2456 // The method reads the raw-data error log
2457 // accumulated within the rawReader.
2458 // It extracts the raw-data errors related to
2459 // the current event and stores them into
2460 // a TClonesArray inside the esd object.
2462 if (!fRawReader) return;
2464 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2466 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2468 if (iEvent != log->GetEventNumber()) continue;
2470 esd->AddRawDataErrorLog(log);
2475 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2476 // Dump a file content into a char in TNamed
2478 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2479 Int_t kBytes = (Int_t)in.tellg();
2480 printf("Size: %d \n",kBytes);
2483 char* memblock = new char [kBytes];
2484 in.seekg (0, ios::beg);
2485 in.read (memblock, kBytes);
2487 TString fData(memblock,kBytes);
2488 fn = new TNamed(pName,fData);
2489 printf("fData Size: %d \n",fData.Sizeof());
2490 printf("pName Size: %d \n",pName.Sizeof());
2491 printf("fn Size: %d \n",fn->Sizeof());
2495 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2501 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2502 // This is not really needed in AliReconstruction at the moment
2503 // but can serve as a template
2505 TList *fList = fTree->GetUserInfo();
2506 TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2507 printf("fn Size: %d \n",fn->Sizeof());
2509 TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2510 const char* cdata = fn->GetTitle();
2511 printf("fTmp Size %d\n",fTmp.Sizeof());
2513 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2514 printf("calculated size %d\n",size);
2515 ofstream out(pName.Data(),ios::out | ios::binary);
2516 out.write(cdata,size);
2521 //_____________________________________________________________________________
2522 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2524 // get the quality assurance data maker object and the loader for a detector
2526 if (fQADataMaker[iDet])
2527 return fQADataMaker[iDet];
2529 AliQADataMakerRec * qadm = NULL;
2530 if (iDet == fgkNDetectors) { //Global QA
2531 qadm = new AliGlobalQADataMaker();
2532 fQADataMaker[iDet] = qadm;
2536 // load the QA data maker object
2537 TPluginManager* pluginManager = gROOT->GetPluginManager();
2538 TString detName = fgkDetectorName[iDet];
2539 TString qadmName = "Ali" + detName + "QADataMakerRec";
2540 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT"))
2543 // first check if a plugin is defined for the quality assurance data maker
2544 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2545 // if not, add a plugin for it
2546 if (!pluginHandler) {
2547 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
2548 TString libs = gSystem->GetLibraries();
2549 if (libs.Contains("lib" + detName + "base.so") ||
2550 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2551 pluginManager->AddHandler("AliQADataMakerRec", detName,
2552 qadmName, detName + "qadm", qadmName + "()");
2554 pluginManager->AddHandler("AliQADataMakerRec", detName,
2555 qadmName, detName, qadmName + "()");
2557 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2559 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2560 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
2563 fQADataMaker[iDet] = qadm;
2568 //_____________________________________________________________________________
2569 Bool_t AliReconstruction::RunQA(AliESDEvent *& esd)
2571 // run the Quality Assurance data producer
2573 AliCodeTimerAuto("")
2574 TString detStr = fQADetectors ;
2575 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2576 if (!IsSelected(fgkDetectorName[iDet], detStr))
2578 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
2581 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2582 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2584 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
2585 qadm->Exec(AliQA::kESDS, esd) ;
2588 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2590 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2591 AliError(Form("the following detectors were not found: %s",
2601 //_____________________________________________________________________________
2602 void AliReconstruction::CheckQA()
2604 // check the QA of SIM for this run and remove the detectors
2605 // with status Fatal
2607 TString newRunLocalReconstruction ;
2608 TString newRunTracking ;
2609 TString newFillESD ;
2611 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2612 TString detName(AliQA::GetDetName(iDet)) ;
2613 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2614 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2615 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2617 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2618 fRunLocalReconstruction.Contains("ALL") ) {
2619 newRunLocalReconstruction += detName ;
2620 newRunLocalReconstruction += " " ;
2622 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2623 fRunTracking.Contains("ALL") ) {
2624 newRunTracking += detName ;
2625 newRunTracking += " " ;
2627 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2628 fFillESD.Contains("ALL") ) {
2629 newFillESD += detName ;
2634 fRunLocalReconstruction = newRunLocalReconstruction ;
2635 fRunTracking = newRunTracking ;
2636 fFillESD = newFillESD ;
2639 //_____________________________________________________________________________
2640 Int_t AliReconstruction::GetDetIndex(const char* detector)
2642 // return the detector index corresponding to detector
2644 for (index = 0; index < fgkNDetectors ; index++) {
2645 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2650 //_____________________________________________________________________________
2651 Bool_t AliReconstruction::FinishPlaneEff() {
2653 // Here execute all the necessary operationis, at the end of the tracking phase,
2654 // in case that evaluation of PlaneEfficiencies was required for some detector.
2655 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2657 // This Preliminary version works only FOR ITS !!!!!
2658 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2661 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2664 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2665 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2666 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2667 if(fTracker[iDet]) {
2668 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2669 ret=planeeff->WriteIntoCDB();
2670 if(planeeff->GetCreateHistos()) {
2671 TString name="PlaneEffHisto";
2672 name+=fgkDetectorName[iDet];
2674 ret*=planeeff->WriteHistosToFile(name,"RECREATE");
2680 //_____________________________________________________________________________
2681 Bool_t AliReconstruction::InitPlaneEff() {
2683 // Here execute all the necessary operations, before of the tracking phase,
2684 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2685 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2686 // which should be updated/recalculated.
2688 // This Preliminary version will work only FOR ITS !!!!!
2689 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2692 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2694 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2698 //_____________________________________________________________________________
2699 Bool_t AliReconstruction::InitAliEVE()
2701 // This method should be called only in case
2702 // AliReconstruction is run
2703 // within the alieve environment.
2704 // It will initialize AliEVE in a way
2705 // so that it can visualize event processed
2706 // by AliReconstruction.
2707 // The return flag shows whenever the
2708 // AliEVE initialization was successful or not.
2711 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2712 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2713 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2715 gROOT->ProcessLine("if (!gAliEveEvent) {gAliEveEvent = new AliEveEventManager();gAliEveEvent->SetAutoLoad(kTRUE);gAliEveEvent->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(gAliEveEvent);};");
2716 gROOT->ProcessLine("alieve_online_init()");
2721 //_____________________________________________________________________________
2722 void AliReconstruction::RunAliEVE()
2724 // Runs AliEVE visualisation of
2725 // the current event.
2726 // Should be executed only after
2727 // successful initialization of AliEVE.
2729 AliInfo("Running AliEVE...");
2730 gROOT->ProcessLine(Form("gAliEveEvent->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2731 gROOT->ProcessLine("gAliEveEvent->StartStopAutoLoadTimer();");
2735 //_____________________________________________________________________________
2736 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2738 // Allows to run QA for a selected set of detectors
2739 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2740 // all selected detectors run the same selected tasks
2742 if (!detAndAction.Contains(":")) {
2743 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2747 Int_t colon = detAndAction.Index(":") ;
2748 fQADetectors = detAndAction(0, colon) ;
2749 if (fQADetectors.Contains("ALL") )
2750 fQADetectors = fFillESD ;
2751 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2752 if (fQATasks.Contains("ALL") ) {
2753 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2755 fQATasks.ToUpper() ;
2757 if ( fQATasks.Contains("RAW") )
2758 tempo = Form("%d ", AliQA::kRAWS) ;
2759 if ( fQATasks.Contains("RECPOINT") )
2760 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2761 if ( fQATasks.Contains("ESD") )
2762 tempo += Form("%d ", AliQA::kESDS) ;
2764 if (fQATasks.IsNull()) {
2765 AliInfo("No QA requested\n") ;
2770 TString tempo(fQATasks) ;
2771 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2772 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2773 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2774 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;