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),
259 fSameQACycle(kFALSE),
261 fRunPlaneEff(kFALSE),
273 fIsNewRunLoader(kFALSE),
276 // create reconstruction object with default parameters
278 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
279 fReconstructor[iDet] = NULL;
280 fLoader[iDet] = NULL;
281 fTracker[iDet] = NULL;
282 fQADataMaker[iDet] = NULL;
283 fQACycles[iDet] = 999999;
285 fQADataMaker[fgkNDetectors]=NULL; //Global QA
289 //_____________________________________________________________________________
290 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
293 fUniformField(rec.fUniformField),
294 fRunVertexFinder(rec.fRunVertexFinder),
295 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
296 fRunHLTTracking(rec.fRunHLTTracking),
297 fRunMuonTracking(rec.fRunMuonTracking),
298 fRunV0Finder(rec.fRunV0Finder),
299 fRunCascadeFinder(rec.fRunCascadeFinder),
300 fStopOnError(rec.fStopOnError),
301 fWriteAlignmentData(rec.fWriteAlignmentData),
302 fWriteESDfriend(rec.fWriteESDfriend),
303 fWriteAOD(rec.fWriteAOD),
304 fFillTriggerESD(rec.fFillTriggerESD),
306 fCleanESD(rec.fCleanESD),
307 fV0DCAmax(rec.fV0DCAmax),
308 fV0CsPmin(rec.fV0CsPmin),
312 fRunLocalReconstruction(rec.fRunLocalReconstruction),
313 fRunTracking(rec.fRunTracking),
314 fFillESD(rec.fFillESD),
315 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
316 fGAliceFileName(rec.fGAliceFileName),
318 fEquipIdMap(rec.fEquipIdMap),
319 fFirstEvent(rec.fFirstEvent),
320 fLastEvent(rec.fLastEvent),
321 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
324 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
325 fLoadAlignData(rec.fLoadAlignData),
326 fESDPar(rec.fESDPar),
327 fUseHLTData(rec.fUseHLTData),
331 fParentRawReader(NULL),
334 fDiamondProfile(NULL),
335 fDiamondProfileTPC(NULL),
336 fMeanVertexConstraint(rec.fMeanVertexConstraint),
340 fAlignObjArray(rec.fAlignObjArray),
341 fCDBUri(rec.fCDBUri),
343 fInitCDBCalled(rec.fInitCDBCalled),
344 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
346 fRunGlobalQA(rec.fRunGlobalQA),
347 fInLoopQA(rec.fInLoopQA),
348 fSameQACycle(rec.fSameQACycle),
349 fRunPlaneEff(rec.fRunPlaneEff),
361 fIsNewRunLoader(rec.fIsNewRunLoader),
366 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
367 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
369 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
370 fReconstructor[iDet] = NULL;
371 fLoader[iDet] = NULL;
372 fTracker[iDet] = NULL;
373 fQADataMaker[iDet] = NULL;
374 fQACycles[iDet] = rec.fQACycles[iDet];
376 fQADataMaker[fgkNDetectors]=NULL; //Global QA
377 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
378 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
382 //_____________________________________________________________________________
383 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
385 // assignment operator
387 this->~AliReconstruction();
388 new(this) AliReconstruction(rec);
392 //_____________________________________________________________________________
393 AliReconstruction::~AliReconstruction()
399 fSpecCDBUri.Delete();
401 AliCodeTimer::Instance()->Print();
404 //_____________________________________________________________________________
405 void AliReconstruction::InitCDB()
407 // activate a default CDB storage
408 // First check if we have any CDB storage set, because it is used
409 // to retrieve the calibration and alignment constants
411 if (fInitCDBCalled) return;
412 fInitCDBCalled = kTRUE;
414 AliCDBManager* man = AliCDBManager::Instance();
415 if (man->IsDefaultStorageSet())
417 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
418 AliWarning("Default CDB storage has been already set !");
419 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
420 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
421 fCDBUri = man->GetDefaultStorage()->GetURI();
424 if (fCDBUri.Length() > 0)
426 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
427 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
428 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
430 fCDBUri="local://$ALICE_ROOT";
431 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
432 AliWarning("Default CDB storage not yet set !!!!");
433 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
434 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
437 man->SetDefaultStorage(fCDBUri);
440 // Now activate the detector specific CDB storage locations
441 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
442 TObject* obj = fSpecCDBUri[i];
444 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
445 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
446 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
447 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
452 //_____________________________________________________________________________
453 void AliReconstruction::SetDefaultStorage(const char* uri) {
454 // Store the desired default CDB storage location
455 // Activate it later within the Run() method
461 //_____________________________________________________________________________
462 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
463 // Store a detector-specific CDB storage location
464 // Activate it later within the Run() method
466 AliCDBPath aPath(calibType);
467 if(!aPath.IsValid()){
468 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
469 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
470 if(!strcmp(calibType, fgkDetectorName[iDet])) {
471 aPath.SetPath(Form("%s/*", calibType));
472 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
476 if(!aPath.IsValid()){
477 AliError(Form("Not a valid path or detector: %s", calibType));
482 // // check that calibType refers to a "valid" detector name
483 // Bool_t isDetector = kFALSE;
484 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
485 // TString detName = fgkDetectorName[iDet];
486 // if(aPath.GetLevel0() == detName) {
487 // isDetector = kTRUE;
493 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
497 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
498 if (obj) fSpecCDBUri.Remove(obj);
499 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
503 //_____________________________________________________________________________
504 Bool_t AliReconstruction::SetRunNumberFromData()
506 // The method is called in Run() in order
507 // to set a correct run number.
508 // In case of raw data reconstruction the
509 // run number is taken from the raw data header
511 if (fSetRunNumberFromDataCalled) return kTRUE;
512 fSetRunNumberFromDataCalled = kTRUE;
514 AliCDBManager* man = AliCDBManager::Instance();
516 if(man->GetRun() > 0) {
517 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
521 AliError("No run loader is found !");
524 // read run number from gAlice
525 if(fRunLoader->GetAliRun())
526 AliCDBManager::Instance()->SetRun(fRunLoader->GetHeader()->GetRun());
529 if(fRawReader->NextEvent()) {
530 AliCDBManager::Instance()->SetRun(fRawReader->GetRunNumber());
531 fRawReader->RewindEvents();
534 if(man->GetRun() > 0) {
535 AliWarning("No raw events is found ! Using settings in AliCDBManager !");
540 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
546 AliError("Neither gAlice nor RawReader objects are found !");
556 //_____________________________________________________________________________
557 void AliReconstruction::SetCDBLock() {
558 // Set CDB lock: from now on it is forbidden to reset the run number
559 // or the default storage or to activate any further storage!
561 AliCDBManager::Instance()->SetLock(1);
564 //_____________________________________________________________________________
565 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
567 // Read the alignment objects from CDB.
568 // Each detector is supposed to have the
569 // alignment objects in DET/Align/Data CDB path.
570 // All the detector objects are then collected,
571 // sorted by geometry level (starting from ALIC) and
572 // then applied to the TGeo geometry.
573 // Finally an overlaps check is performed.
575 // Load alignment data from CDB and fill fAlignObjArray
576 if(fLoadAlignFromCDB){
578 TString detStr = detectors;
579 TString loadAlObjsListOfDets = "";
581 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
582 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
583 loadAlObjsListOfDets += fgkDetectorName[iDet];
584 loadAlObjsListOfDets += " ";
585 } // end loop over detectors
586 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
587 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
589 // Check if the array with alignment objects was
590 // provided by the user. If yes, apply the objects
591 // to the present TGeo geometry
592 if (fAlignObjArray) {
593 if (gGeoManager && gGeoManager->IsClosed()) {
594 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
595 AliError("The misalignment of one or more volumes failed!"
596 "Compare the list of simulated detectors and the list of detector alignment data!");
601 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
607 delete fAlignObjArray; fAlignObjArray=0;
612 //_____________________________________________________________________________
613 void AliReconstruction::SetGAliceFile(const char* fileName)
615 // set the name of the galice file
617 fGAliceFileName = fileName;
620 //_____________________________________________________________________________
621 void AliReconstruction::SetInput(const char* input)
623 // In case the input string starts with 'mem://', we run in an online mode
624 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
625 // file is assumed. One can give as an input:
626 // mem://: - events taken from DAQ monitoring libs online
628 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
632 //_____________________________________________________________________________
633 void AliReconstruction::SetOption(const char* detector, const char* option)
635 // set options for the reconstruction of a detector
637 TObject* obj = fOptions.FindObject(detector);
638 if (obj) fOptions.Remove(obj);
639 fOptions.Add(new TNamed(detector, option));
642 //_____________________________________________________________________________
643 Bool_t AliReconstruction::Run(const char* input)
646 AliCodeTimerAuto("");
648 if (!InitRun(input)) return kFALSE;
650 //******* The loop over events
652 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
653 (fRawReader && fRawReader->NextEvent())) {
654 if (!RunEvent(iEvent)) return kFALSE;
658 if (!FinishRun()) return kFALSE;
663 //_____________________________________________________________________________
664 Bool_t AliReconstruction::InitRun(const char* input)
666 // Initialize all the stuff before
667 // going into the event loop
668 // If the second argument is given, the first one is ignored and
669 // the reconstruction works in an online mode
670 AliCodeTimerAuto("");
672 // Overwrite the previous setting
673 if (input) fInput = input;
675 // set the input in case of raw data
676 fRawReader = AliRawReader::Create(fInput.Data());
678 AliInfo("Reconstruction will run over digits");
680 if (!fEquipIdMap.IsNull() && fRawReader)
681 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
683 if (!fUseHLTData.IsNull()) {
684 // create the RawReaderHLT which performs redirection of HLT input data for
685 // the specified detectors
686 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
688 fParentRawReader=fRawReader;
689 fRawReader=pRawReader;
691 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
695 AliSysInfo::AddStamp("Start");
696 // get the run loader
697 if (!InitRunLoader()) return kFALSE;
698 AliSysInfo::AddStamp("LoadLoader");
700 // Initialize the CDB storage
703 AliSysInfo::AddStamp("LoadCDB");
705 // Set run number in CDBManager (if it is not already set by the user)
706 if (!SetRunNumberFromData()) if (fStopOnError) return kFALSE;
708 // Set CDB lock: from now on it is forbidden to reset the run number
709 // or the default storage or to activate any further storage!
712 // Import ideal TGeo geometry and apply misalignment
714 TString geom(gSystem->DirName(fGAliceFileName));
715 geom += "/geometry.root";
716 AliGeomManager::LoadGeometry(geom.Data());
718 TString detsToCheck=fRunLocalReconstruction;
719 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data()))
720 AliFatalClass("Current loaded geometry differs in the definition of symbolic names!");
721 if (!gGeoManager) if (fStopOnError) return kFALSE;
724 if (!MisalignGeometry(fLoadAlignData)) if (fStopOnError) return kFALSE;
725 AliSysInfo::AddStamp("LoadGeom");
728 // Get the GRP CDB entry
729 AliCDBEntry* entryGRP = AliCDBManager::Instance()->Get("GRP/GRP/Data");
732 fGRPData = dynamic_cast<TMap*>(entryGRP->GetObject());
735 AliError("No GRP entry found in OCDB!");
740 // Magnetic field map
741 if (!AliTracker::GetFieldMap()) {
742 // Construct the field map out of the information retrieved from GRP.
744 // For the moment, this is a dummy piece of code.
745 // The actual map is expected to be already created in rec.C !
749 Int_t map=AliMagWrapCheb::k5kG;
750 Bool_t dipoleON=kTRUE;
753 TObjString *l3Current=
754 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
756 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
759 TObjString *l3Polarity=
760 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
762 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
767 TObjString *diCurrent=
768 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
770 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
773 TObjString *diPolarity=
774 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
776 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
782 new AliMagWrapCheb("Maps","Maps",2,factor,10.,map,dipoleON);
783 AliTracker::SetFieldMap(field,fUniformField);
786 AliFatal("Please, provide the field map ! Crashing deliberately...");
791 // Get the diamond profile from OCDB
792 AliCDBEntry* entry = AliCDBManager::Instance()
793 ->Get("GRP/Calib/MeanVertex");
796 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
798 AliError("No diamond profile found in OCDB!");
802 entry = AliCDBManager::Instance()
803 ->Get("GRP/Calib/MeanVertexTPC");
806 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
808 AliError("No diamond profile found in OCDB!");
811 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
812 if(fDiamondProfile && fMeanVertexConstraint) ftVertexer->SetVtxStart(fDiamondProfile);
815 if (fRunVertexFinder && !CreateVertexer()) {
821 AliSysInfo::AddStamp("Vertexer");
824 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
830 AliSysInfo::AddStamp("LoadTrackers");
832 // get the possibly already existing ESD file and tree
833 fesd = new AliESDEvent(); fhltesd = new AliESDEvent();
834 if (!gSystem->AccessPathName("AliESDs.root")){
835 gSystem->CopyFile("AliESDs.root", "AliESDs.old.root", kTRUE);
836 ffileOld = TFile::Open("AliESDs.old.root");
837 if (ffileOld && ffileOld->IsOpen()) {
838 ftreeOld = (TTree*) ffileOld->Get("esdTree");
839 if (ftreeOld)fesd->ReadFromTree(ftreeOld);
840 fhlttreeOld = (TTree*) ffileOld->Get("HLTesdTree");
841 if (fhlttreeOld) fhltesd->ReadFromTree(fhlttreeOld);
845 // create the ESD output file and tree
846 ffile = TFile::Open("AliESDs.root", "RECREATE");
847 ffile->SetCompressionLevel(2);
848 if (!ffile->IsOpen()) {
849 AliError("opening AliESDs.root failed");
850 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
853 ftree = new TTree("esdTree", "Tree with ESD objects");
854 fesd = new AliESDEvent();
855 fesd->CreateStdContent();
856 fesd->WriteToTree(ftree);
858 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
859 fhltesd = new AliESDEvent();
860 fhltesd->CreateStdContent();
861 fhltesd->WriteToTree(fhlttree);
864 if (fWriteESDfriend) {
865 fesdf = new AliESDfriend();
866 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
867 br->SetFile("AliESDfriends.root");
868 fesd->AddObject(fesdf);
873 if (fRawReader) fRawReader->RewindEvents();
876 gSystem->GetProcInfo(&ProcInfo);
877 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
881 AliQADataMakerSteer qas ;
882 if (fRunQA && fRawReader) {
883 qas.Run(fRunLocalReconstruction, fRawReader) ;
884 fSameQACycle = kTRUE ;
886 //Initialize the QA and start of cycle for out-of-cycle QA
888 // TString detStr(fFillESD);
889 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
890 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
891 // AliQADataMakerRec *qadm = GetQADataMaker(iDet);
892 // if (!qadm) continue;
893 // AliInfo(Form("Initializing the QA data maker for %s",
894 // fgkDetectorName[iDet]));
895 // qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
896 // qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
898 // qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
899 // qadm->StartOfCycle(AliQA::kESDS,"same");
903 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
904 AliInfo(Form("Initializing the global QA data maker"));
906 qadm->Init(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun());
907 AliTracker::SetResidualsArray(arr);
908 qadm->Init(AliQA::kESDS, AliCDBManager::Instance()->GetRun());
910 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
911 qadm->StartOfCycle(AliQA::kESDS, "same");
915 fSameQACycle = kTRUE;
918 //Initialize the Plane Efficiency framework
919 if (fRunPlaneEff && !InitPlaneEff()) {
920 if(fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
923 if (strcmp(gProgName,"alieve") == 0)
924 fRunAliEVE = InitAliEVE();
929 //_____________________________________________________________________________
930 Bool_t AliReconstruction::RunEvent(Int_t iEvent)
932 // run the reconstruction over a single event
933 // The event loop is steered in Run method
935 AliCodeTimerAuto("");
937 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
938 fRunLoader->SetEventNumber(iEvent);
939 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
941 //?? fRunLoader->MakeTree("H");
942 fRunLoader->TreeE()->Fill();
945 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
946 // copy old ESD to the new one
948 fesd->ReadFromTree(ftreeOld);
949 ftreeOld->GetEntry(iEvent);
953 fhltesd->ReadFromTree(fhlttreeOld);
954 fhlttreeOld->GetEntry(iEvent);
960 AliInfo(Form("processing event %d", iEvent));
962 //Start of cycle for the in-loop QA
965 TString detStr(fFillESD);
966 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
967 if (!IsSelected(fgkDetectorName[iDet], detStr))
969 AliQADataMakerRec *qadm = GetQADataMaker(iDet);
972 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
973 qadm->StartOfCycle(AliQA::kESDS, "same") ;
976 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
977 qadm->StartOfCycle(AliQA::kRECPOINTS, fSameQACycle);
978 qadm->StartOfCycle(AliQA::kESDS, "same");
983 fRunLoader->GetEvent(iEvent);
986 sprintf(aFileName, "ESD_%d.%d_final.root",
987 fRunLoader->GetHeader()->GetRun(),
988 fRunLoader->GetHeader()->GetEventNrInRun());
989 if (!gSystem->AccessPathName(aFileName)) return kTRUE;
991 // local single event reconstruction
992 if (!fRunLocalReconstruction.IsNull()) {
993 TString detectors=fRunLocalReconstruction;
994 // run HLT event reconstruction first
995 // ;-( IsSelected changes the string
996 if (IsSelected("HLT", detectors) &&
997 !RunLocalEventReconstruction("HLT")) {
998 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1000 detectors=fRunLocalReconstruction;
1001 detectors.ReplaceAll("HLT", "");
1002 if (!RunLocalEventReconstruction(detectors)) {
1003 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1007 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1008 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1009 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1010 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1012 // Set magnetic field from the tracker
1013 fesd->SetMagneticField(AliTracker::GetBz());
1014 fhltesd->SetMagneticField(AliTracker::GetBz());
1018 // Fill raw-data error log into the ESD
1019 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1022 if (fRunVertexFinder) {
1023 if (!ReadESD(fesd, "vertex")) {
1024 if (!RunVertexFinder(fesd)) {
1025 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1027 if (fCheckPointLevel > 0) WriteESD(fesd, "vertex");
1032 if (!fRunTracking.IsNull()) {
1033 if (fRunMuonTracking) {
1034 if (!RunMuonTracking(fesd)) {
1035 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1041 if (!fRunTracking.IsNull()) {
1042 if (!ReadESD(fesd, "tracking")) {
1043 if (!RunTracking(fesd)) {
1044 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1046 if (fCheckPointLevel > 0) WriteESD(fesd, "tracking");
1051 if (!fFillESD.IsNull()) {
1052 TString detectors=fFillESD;
1053 // run HLT first and on hltesd
1054 // ;-( IsSelected changes the string
1055 if (IsSelected("HLT", detectors) &&
1056 !FillESD(fhltesd, "HLT")) {
1057 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1060 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1061 if (detectors.Contains("ALL")) {
1063 for (Int_t idet=0; idet<fgkNDetectors; ++idet){
1064 detectors += fgkDetectorName[idet];
1068 detectors.ReplaceAll("HLT", "");
1069 if (!FillESD(fesd, detectors)) {
1070 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1074 // fill Event header information from the RawEventHeader
1075 if (fRawReader){FillRawEventHeaderESD(fesd);}
1078 AliESDpid::MakePID(fesd);
1079 if (fCheckPointLevel > 1) WriteESD(fesd, "PID");
1081 if (fFillTriggerESD) {
1082 if (!ReadESD(fesd, "trigger")) {
1083 if (!FillTriggerESD(fesd)) {
1084 if (fStopOnError) {CleanUp(ffile, ffileOld); return kFALSE;}
1086 if (fCheckPointLevel > 1) WriteESD(fesd, "trigger");
1093 // Propagate track to the beam pipe (if not already done by ITS)
1095 const Int_t ntracks = fesd->GetNumberOfTracks();
1096 const Double_t kBz = fesd->GetMagneticField();
1097 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1100 UShort_t *selectedIdx=new UShort_t[ntracks];
1102 for (Int_t itrack=0; itrack<ntracks; itrack++){
1103 const Double_t kMaxStep = 5; //max step over the material
1106 AliESDtrack *track = fesd->GetTrack(itrack);
1107 if (!track) continue;
1109 AliExternalTrackParam *tpcTrack =
1110 (AliExternalTrackParam *)track->GetTPCInnerParam();
1114 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
1119 Int_t n=trkArray.GetEntriesFast();
1120 selectedIdx[n]=track->GetID();
1121 trkArray.AddLast(tpcTrack);
1124 if (track->GetX() < kRadius) continue;
1127 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1129 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kRadius);
1134 // Improve the reconstructed primary vertex position using the tracks
1136 TObject *obj = fOptions.FindObject("ITS");
1138 TString optITS = obj->GetTitle();
1139 if (optITS.Contains("cosmics") || optITS.Contains("COSMICS"))
1140 fRunVertexFinderTracks=kFALSE;
1142 if (fRunVertexFinderTracks) {
1143 // TPC + ITS primary vertex
1144 ftVertexer->SetITSrefitRequired();
1145 if(fDiamondProfile && fMeanVertexConstraint) {
1146 ftVertexer->SetVtxStart(fDiamondProfile);
1148 ftVertexer->SetConstraintOff();
1150 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1152 if (pvtx->GetStatus()) {
1153 fesd->SetPrimaryVertex(pvtx);
1154 for (Int_t i=0; i<ntracks; i++) {
1155 AliESDtrack *t = fesd->GetTrack(i);
1156 t->RelateToVertex(pvtx, kBz, kRadius);
1161 // TPC-only primary vertex
1162 ftVertexer->SetITSrefitNotRequired();
1163 if(fDiamondProfileTPC && fMeanVertexConstraint) {
1164 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1166 ftVertexer->SetConstraintOff();
1168 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1170 if (pvtx->GetStatus()) {
1171 fesd->SetPrimaryVertexTPC(pvtx);
1172 Int_t nsel=trkArray.GetEntriesFast();
1173 for (Int_t i=0; i<nsel; i++) {
1174 AliExternalTrackParam *t =
1175 (AliExternalTrackParam *)trkArray.UncheckedAt(i);
1176 t->PropagateToDCA(pvtx, kBz, kRadius);
1182 delete[] selectedIdx;
1184 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1189 AliV0vertexer vtxer;
1190 vtxer.Tracks2V0vertices(fesd);
1192 if (fRunCascadeFinder) {
1194 AliCascadeVertexer cvtxer;
1195 cvtxer.V0sTracks2CascadeVertices(fesd);
1200 if (fCleanESD) CleanESD(fesd);
1204 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1205 if (qadm) qadm->Exec(AliQA::kESDS, fesd);
1209 if (fWriteESDfriend) {
1210 fesdf->~AliESDfriend();
1211 new (fesdf) AliESDfriend(); // Reset...
1212 fesd->GetESDfriend(fesdf);
1220 if (fRunAliEVE) RunAliEVE();
1222 if (fCheckPointLevel > 0) WriteESD(fesd, "final");
1225 if (fWriteESDfriend) {
1226 fesdf->~AliESDfriend();
1227 new (fesdf) AliESDfriend(); // Reset...
1230 ProcInfo_t ProcInfo;
1231 gSystem->GetProcInfo(&ProcInfo);
1232 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1235 // End of cycle for the in-loop QA
1238 RunQA(fFillESD.Data(), fesd);
1239 TString detStr(fFillESD);
1240 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1241 if (!IsSelected(fgkDetectorName[iDet], detStr))
1243 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1246 qadm->EndOfCycle(AliQA::kRECPOINTS);
1247 qadm->EndOfCycle(AliQA::kESDS);
1252 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1254 qadm->EndOfCycle(AliQA::kRECPOINTS);
1255 qadm->EndOfCycle(AliQA::kESDS);
1264 //_____________________________________________________________________________
1265 Bool_t AliReconstruction::FinishRun()
1268 // Called after the exit
1269 // from the event loop
1270 AliCodeTimerAuto("");
1272 if (fIsNewRunLoader) { // galice.root didn't exist
1273 fRunLoader->WriteHeader("OVERWRITE");
1274 fRunLoader->CdGAFile();
1275 fRunLoader->Write(0, TObject::kOverwrite);
1278 ftree->GetUserInfo()->Add(fesd);
1279 fhlttree->GetUserInfo()->Add(fhltesd);
1281 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1282 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1284 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1285 cdbMapCopy->SetOwner(1);
1286 cdbMapCopy->SetName("cdbMap");
1287 TIter iter(cdbMap->GetTable());
1290 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1291 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1292 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1293 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1296 TList *cdbListCopy = new TList();
1297 cdbListCopy->SetOwner(1);
1298 cdbListCopy->SetName("cdbList");
1300 TIter iter2(cdbList);
1303 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1304 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1307 ftree->GetUserInfo()->Add(cdbMapCopy);
1308 ftree->GetUserInfo()->Add(cdbListCopy);
1311 if(fESDPar.Contains("ESD.par")){
1312 AliInfo("Attaching ESD.par to Tree");
1313 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1314 ftree->GetUserInfo()->Add(fn);
1320 if (fWriteESDfriend)
1321 ftree->SetBranchStatus("ESDfriend*",0);
1322 // we want to have only one tree version number
1323 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1326 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1327 if (fRunPlaneEff && !FinishPlaneEff()) {
1328 AliWarning("Finish PlaneEff evaluation failed");
1332 CleanUp(ffile, ffileOld);
1335 AliWarning("AOD creation not supported anymore during reconstruction. See ANALYSIS/AliAnalysisTaskESDfilter.cxx instead.");
1338 // Create tags for the events in the ESD tree (the ESD tree is always present)
1339 // In case of empty events the tags will contain dummy values
1340 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1341 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1343 AliWarning("AOD tag creation not supported anymore during reconstruction.");
1346 //Finish QA and end of cycle for out-of-loop QA
1349 AliQADataMakerSteer qas;
1350 qas.Run(fRunLocalReconstruction.Data(), AliQA::kRECPOINTS, fSameQACycle);
1352 qas.Run(fRunLocalReconstruction.Data(), AliQA::kESDS, fSameQACycle);
1354 AliQADataMakerRec *qadm = GetQADataMaker(AliQA::kGLOBAL);
1356 qadm->EndOfCycle(AliQA::kRECPOINTS);
1357 qadm->EndOfCycle(AliQA::kESDS);
1364 // Cleanup of CDB manager: cache and active storages!
1365 AliCDBManager::Instance()->ClearCache();
1371 //_____________________________________________________________________________
1372 Bool_t AliReconstruction::RunLocalReconstruction(const TString& /*detectors*/)
1374 // run the local reconstruction
1375 static Int_t eventNr=0;
1376 AliCodeTimerAuto("")
1378 // AliCDBManager* man = AliCDBManager::Instance();
1379 // Bool_t origCache = man->GetCacheFlag();
1381 // TString detStr = detectors;
1382 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1383 // if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1384 // AliReconstructor* reconstructor = GetReconstructor(iDet);
1385 // if (!reconstructor) continue;
1386 // if (reconstructor->HasLocalReconstruction()) continue;
1388 // AliCodeTimerStart(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1389 // AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1391 // AliCodeTimerStart(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1392 // AliInfo(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1394 // man->SetCacheFlag(kTRUE);
1395 // TString calibPath = Form("%s/Calib/*", fgkDetectorName[iDet]);
1396 // man->GetAll(calibPath); // entries are cached!
1398 // AliCodeTimerStop(Form("Loading calibration data from OCDB for %s", fgkDetectorName[iDet]));
1400 // if (fRawReader) {
1401 // fRawReader->RewindEvents();
1402 // reconstructor->Reconstruct(fRunLoader, fRawReader);
1404 // reconstructor->Reconstruct(fRunLoader);
1407 // AliCodeTimerStop(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1408 // AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr));
1410 // // unload calibration data
1411 // man->UnloadFromCache(calibPath);
1412 // //man->ClearCache();
1415 // man->SetCacheFlag(origCache);
1417 // if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1418 // AliError(Form("the following detectors were not found: %s",
1420 // if (fStopOnError) return kFALSE;
1427 //_____________________________________________________________________________
1428 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1430 // run the local reconstruction
1432 static Int_t eventNr=0;
1433 AliCodeTimerAuto("")
1435 TString detStr = detectors;
1436 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1437 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1438 AliReconstructor* reconstructor = GetReconstructor(iDet);
1439 if (!reconstructor) continue;
1440 AliLoader* loader = fLoader[iDet];
1441 // Matthias April 2008: temporary fix to run HLT reconstruction
1442 // although the HLT loader is missing
1443 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1445 reconstructor->Reconstruct(fRawReader, NULL);
1448 reconstructor->Reconstruct(dummy, NULL);
1453 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1456 // conversion of digits
1457 if (fRawReader && reconstructor->HasDigitConversion()) {
1458 AliInfo(Form("converting raw data digits into root objects for %s",
1459 fgkDetectorName[iDet]));
1460 AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1461 fgkDetectorName[iDet]));
1462 loader->LoadDigits("update");
1463 loader->CleanDigits();
1464 loader->MakeDigitsContainer();
1465 TTree* digitsTree = loader->TreeD();
1466 reconstructor->ConvertDigits(fRawReader, digitsTree);
1467 loader->WriteDigits("OVERWRITE");
1468 loader->UnloadDigits();
1470 // local reconstruction
1471 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1472 AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1473 loader->LoadRecPoints("update");
1474 loader->CleanRecPoints();
1475 loader->MakeRecPointsContainer();
1476 TTree* clustersTree = loader->TreeR();
1477 if (fRawReader && !reconstructor->HasDigitConversion()) {
1478 reconstructor->Reconstruct(fRawReader, clustersTree);
1480 loader->LoadDigits("read");
1481 TTree* digitsTree = loader->TreeD();
1483 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1484 if (fStopOnError) return kFALSE;
1486 reconstructor->Reconstruct(digitsTree, clustersTree);
1488 loader->UnloadDigits();
1491 // In-loop QA for local reconstrucion
1492 if (fRunQA && fInLoopQA) {
1493 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
1496 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1498 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1500 qadm->Exec(AliQA::kRECPOINTS, clustersTree) ;
1503 //(Form("Running QA data maker for %s", fgkDetectorName[iDet]));
1507 loader->WriteRecPoints("OVERWRITE");
1508 loader->UnloadRecPoints();
1509 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1512 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1513 AliError(Form("the following detectors were not found: %s",
1515 if (fStopOnError) return kFALSE;
1521 //_____________________________________________________________________________
1522 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1524 // run the barrel tracking
1526 AliCodeTimerAuto("")
1528 AliESDVertex* vertex = NULL;
1529 Double_t vtxPos[3] = {0, 0, 0};
1530 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1531 TArrayF mcVertex(3);
1532 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1533 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1534 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1538 AliInfo("running the ITS vertex finder");
1540 fLoader[0]->LoadRecPoints();
1541 TTree* cltree = fLoader[0]->TreeR();
1543 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1544 vertex = fVertexer->FindVertexForCurrentEvent(cltree);
1547 AliError("Can't get the ITS cluster tree");
1549 fLoader[0]->UnloadRecPoints();
1552 AliError("Can't get the ITS loader");
1555 AliWarning("Vertex not found");
1556 vertex = new AliESDVertex();
1557 vertex->SetName("default");
1560 vertex->SetName("reconstructed");
1564 AliInfo("getting the primary vertex from MC");
1565 vertex = new AliESDVertex(vtxPos, vtxErr);
1569 vertex->GetXYZ(vtxPos);
1570 vertex->GetSigmaXYZ(vtxErr);
1572 AliWarning("no vertex reconstructed");
1573 vertex = new AliESDVertex(vtxPos, vtxErr);
1575 esd->SetPrimaryVertexSPD(vertex);
1576 // if SPD multiplicity has been determined, it is stored in the ESD
1577 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1578 if(mult)esd->SetMultiplicity(mult);
1580 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1581 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1588 //_____________________________________________________________________________
1589 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1591 // run the HLT barrel tracking
1593 AliCodeTimerAuto("")
1596 AliError("Missing runLoader!");
1600 AliInfo("running HLT tracking");
1602 // Get a pointer to the HLT reconstructor
1603 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1604 if (!reconstructor) return kFALSE;
1607 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1608 TString detName = fgkDetectorName[iDet];
1609 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1610 reconstructor->SetOption(detName.Data());
1611 AliTracker *tracker = reconstructor->CreateTracker();
1613 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1614 if (fStopOnError) return kFALSE;
1618 Double_t vtxErr[3]={0.005,0.005,0.010};
1619 const AliESDVertex *vertex = esd->GetVertex();
1620 vertex->GetXYZ(vtxPos);
1621 tracker->SetVertex(vtxPos,vtxErr);
1623 fLoader[iDet]->LoadRecPoints("read");
1624 TTree* tree = fLoader[iDet]->TreeR();
1626 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1629 tracker->LoadClusters(tree);
1631 if (tracker->Clusters2Tracks(esd) != 0) {
1632 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1636 tracker->UnloadClusters();
1644 //_____________________________________________________________________________
1645 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1647 // run the muon spectrometer tracking
1649 AliCodeTimerAuto("")
1652 AliError("Missing runLoader!");
1655 Int_t iDet = 7; // for MUON
1657 AliInfo("is running...");
1659 // Get a pointer to the MUON reconstructor
1660 AliReconstructor *reconstructor = GetReconstructor(iDet);
1661 if (!reconstructor) return kFALSE;
1664 TString detName = fgkDetectorName[iDet];
1665 AliDebug(1, Form("%s tracking", detName.Data()));
1666 AliTracker *tracker = reconstructor->CreateTracker();
1668 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1673 fLoader[iDet]->LoadRecPoints("read");
1675 tracker->LoadClusters(fLoader[iDet]->TreeR());
1677 Int_t rv = tracker->Clusters2Tracks(esd);
1681 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1685 fLoader[iDet]->UnloadRecPoints();
1687 tracker->UnloadClusters();
1695 //_____________________________________________________________________________
1696 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1698 // run the barrel tracking
1699 static Int_t eventNr=0;
1700 AliCodeTimerAuto("")
1702 AliInfo("running tracking");
1704 //Fill the ESD with the T0 info (will be used by the TOF)
1705 if (fReconstructor[11] && fLoader[11]) {
1706 fLoader[11]->LoadRecPoints("READ");
1707 TTree *treeR = fLoader[11]->TreeR();
1708 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
1711 // pass 1: TPC + ITS inwards
1712 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1713 if (!fTracker[iDet]) continue;
1714 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
1717 fLoader[iDet]->LoadRecPoints("read");
1718 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
1719 TTree* tree = fLoader[iDet]->TreeR();
1721 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1724 fTracker[iDet]->LoadClusters(tree);
1725 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1727 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
1728 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1731 if (fCheckPointLevel > 1) {
1732 WriteESD(esd, Form("%s.tracking", fgkDetectorName[iDet]));
1734 // preliminary PID in TPC needed by the ITS tracker
1736 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1737 AliESDpid::MakePID(esd);
1739 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
1742 // pass 2: ALL backwards
1744 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1745 if (!fTracker[iDet]) continue;
1746 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
1749 if (iDet > 1) { // all except ITS, TPC
1751 fLoader[iDet]->LoadRecPoints("read");
1752 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
1753 tree = fLoader[iDet]->TreeR();
1755 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
1758 fTracker[iDet]->LoadClusters(tree);
1759 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
1763 if (iDet>1) // start filling residuals for the "outer" detectors
1764 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1766 if (fTracker[iDet]->PropagateBack(esd) != 0) {
1767 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
1770 if (fCheckPointLevel > 1) {
1771 WriteESD(esd, Form("%s.back", fgkDetectorName[iDet]));
1775 if (iDet > 2) { // all except ITS, TPC, TRD
1776 fTracker[iDet]->UnloadClusters();
1777 fLoader[iDet]->UnloadRecPoints();
1779 // updated PID in TPC needed by the ITS tracker -MI
1781 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
1782 AliESDpid::MakePID(esd);
1784 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1786 //stop filling residuals for the "outer" detectors
1787 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1789 // write space-points to the ESD in case alignment data output
1791 if (fWriteAlignmentData)
1792 WriteAlignmentData(esd);
1794 // pass 3: TRD + TPC + ITS refit inwards
1796 for (Int_t iDet = 2; iDet >= 0; iDet--) {
1797 if (!fTracker[iDet]) continue;
1798 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
1801 if (iDet<2) // start filling residuals for TPC and ITS
1802 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
1804 if (fTracker[iDet]->RefitInward(esd) != 0) {
1805 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
1808 // run postprocessing
1809 if (fTracker[iDet]->PostProcess(esd) != 0) {
1810 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
1813 if (fCheckPointLevel > 1) {
1814 WriteESD(esd, Form("%s.refit", fgkDetectorName[iDet]));
1816 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
1818 fTracker[iDet]->UnloadClusters();
1819 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
1820 fLoader[iDet]->UnloadRecPoints();
1821 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
1823 // stop filling residuals for TPC and ITS
1824 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
1830 //_____________________________________________________________________________
1831 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
1833 // Remove the data which are not needed for the physics analysis.
1836 Int_t nTracks=esd->GetNumberOfTracks();
1837 Int_t nV0s=esd->GetNumberOfV0s();
1839 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
1841 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
1842 Bool_t rc=esd->Clean(cleanPars);
1844 nTracks=esd->GetNumberOfTracks();
1845 nV0s=esd->GetNumberOfV0s();
1847 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
1852 //_____________________________________________________________________________
1853 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
1855 // fill the event summary data
1857 AliCodeTimerAuto("")
1858 static Int_t eventNr=0;
1859 TString detStr = detectors;
1861 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1862 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1863 AliReconstructor* reconstructor = GetReconstructor(iDet);
1864 if (!reconstructor) continue;
1865 if (!ReadESD(esd, fgkDetectorName[iDet])) {
1866 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
1867 TTree* clustersTree = NULL;
1868 if (fLoader[iDet]) {
1869 fLoader[iDet]->LoadRecPoints("read");
1870 clustersTree = fLoader[iDet]->TreeR();
1871 if (!clustersTree) {
1872 AliError(Form("Can't get the %s clusters tree",
1873 fgkDetectorName[iDet]));
1874 if (fStopOnError) return kFALSE;
1877 if (fRawReader && !reconstructor->HasDigitConversion()) {
1878 reconstructor->FillESD(fRawReader, clustersTree, esd);
1880 TTree* digitsTree = NULL;
1881 if (fLoader[iDet]) {
1882 fLoader[iDet]->LoadDigits("read");
1883 digitsTree = fLoader[iDet]->TreeD();
1885 AliError(Form("Can't get the %s digits tree",
1886 fgkDetectorName[iDet]));
1887 if (fStopOnError) return kFALSE;
1890 reconstructor->FillESD(digitsTree, clustersTree, esd);
1891 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
1893 if (fLoader[iDet]) {
1894 fLoader[iDet]->UnloadRecPoints();
1897 if (fCheckPointLevel > 2) WriteESD(esd, fgkDetectorName[iDet]);
1901 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1902 AliError(Form("the following detectors were not found: %s",
1904 if (fStopOnError) return kFALSE;
1906 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
1911 //_____________________________________________________________________________
1912 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
1914 // Reads the trigger decision which is
1915 // stored in Trigger.root file and fills
1916 // the corresponding esd entries
1918 AliCodeTimerAuto("")
1920 AliInfo("Filling trigger information into the ESD");
1922 AliCentralTrigger *aCTP = NULL;
1925 AliCTPRawStream input(fRawReader);
1926 if (!input.Next()) {
1927 AliError("No valid CTP (trigger) DDL raw data is found ! The trigger information is not stored in the ESD !");
1930 esd->SetTriggerMask(input.GetClassMask());
1931 esd->SetTriggerCluster(input.GetClusterMask());
1933 aCTP = new AliCentralTrigger();
1934 TString configstr("");
1935 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
1936 AliError("No trigger configuration found in OCDB! The trigger classes information will no be stored in ESD!");
1942 AliRunLoader *runloader = AliRunLoader::GetRunLoader();
1944 if (!runloader->LoadTrigger()) {
1945 aCTP = runloader->GetTrigger();
1946 esd->SetTriggerMask(aCTP->GetClassMask());
1947 esd->SetTriggerCluster(aCTP->GetClusterMask());
1950 AliWarning("No trigger can be loaded! The trigger information is not stored in the ESD !");
1955 AliError("No run loader is available! The trigger information is not stored in the ESD !");
1960 // Now fill the trigger class names into AliESDRun object
1961 AliTriggerConfiguration *config = aCTP->GetConfiguration();
1963 AliError("No trigger configuration has been found! The trigger classes information will no be stored in ESD!");
1964 if (fRawReader) delete aCTP;
1968 const TObjArray& classesArray = config->GetClasses();
1969 Int_t nclasses = classesArray.GetEntriesFast();
1970 for( Int_t j=0; j<nclasses; j++ ) {
1971 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At( j );
1972 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
1973 esd->SetTriggerClass(trclass->GetName(),trindex);
1976 if (fRawReader) delete aCTP;
1984 //_____________________________________________________________________________
1985 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
1988 // Filling information from RawReader Header
1991 AliInfo("Filling information from RawReader Header");
1992 esd->SetBunchCrossNumber(0);
1993 esd->SetOrbitNumber(0);
1994 esd->SetPeriodNumber(0);
1995 esd->SetTimeStamp(0);
1996 esd->SetEventType(0);
1997 const AliRawEventHeaderBase * eventHeader = fRawReader->GetEventHeader();
2000 const UInt_t *id = eventHeader->GetP("Id");
2001 esd->SetBunchCrossNumber((id)[1]&0x00000fff);
2002 esd->SetOrbitNumber((((id)[0]<<20)&0xf00000)|(((id)[1]>>12)&0xfffff));
2003 esd->SetPeriodNumber(((id)[0]>>4)&0x0fffffff);
2005 esd->SetTimeStamp((eventHeader->Get("Timestamp")));
2006 esd->SetEventType((eventHeader->Get("Type")));
2013 //_____________________________________________________________________________
2014 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2016 // check whether detName is contained in detectors
2017 // if yes, it is removed from detectors
2019 // check if all detectors are selected
2020 if ((detectors.CompareTo("ALL") == 0) ||
2021 detectors.BeginsWith("ALL ") ||
2022 detectors.EndsWith(" ALL") ||
2023 detectors.Contains(" ALL ")) {
2028 // search for the given detector
2029 Bool_t result = kFALSE;
2030 if ((detectors.CompareTo(detName) == 0) ||
2031 detectors.BeginsWith(detName+" ") ||
2032 detectors.EndsWith(" "+detName) ||
2033 detectors.Contains(" "+detName+" ")) {
2034 detectors.ReplaceAll(detName, "");
2038 // clean up the detectors string
2039 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2040 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2041 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2046 //_____________________________________________________________________________
2047 Bool_t AliReconstruction::InitRunLoader()
2049 // get or create the run loader
2051 if (gAlice) delete gAlice;
2054 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2055 // load all base libraries to get the loader classes
2056 TString libs = gSystem->GetLibraries();
2057 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2058 TString detName = fgkDetectorName[iDet];
2059 if (detName == "HLT") continue;
2060 if (libs.Contains("lib" + detName + "base.so")) continue;
2061 gSystem->Load("lib" + detName + "base.so");
2063 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2065 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2070 fRunLoader->CdGAFile();
2071 fRunLoader->LoadgAlice();
2073 //PH This is a temporary fix to give access to the kinematics
2074 //PH that is needed for the labels of ITS clusters
2075 fRunLoader->LoadHeader();
2076 fRunLoader->LoadKinematics();
2078 } else { // galice.root does not exist
2080 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2084 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2085 AliConfig::GetDefaultEventFolderName(),
2088 AliError(Form("could not create run loader in file %s",
2089 fGAliceFileName.Data()));
2093 fIsNewRunLoader = kTRUE;
2094 fRunLoader->MakeTree("E");
2096 if (fNumberOfEventsPerFile > 0)
2097 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2099 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2105 //_____________________________________________________________________________
2106 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2108 // get the reconstructor object and the loader for a detector
2110 if (fReconstructor[iDet]) return fReconstructor[iDet];
2112 // load the reconstructor object
2113 TPluginManager* pluginManager = gROOT->GetPluginManager();
2114 TString detName = fgkDetectorName[iDet];
2115 TString recName = "Ali" + detName + "Reconstructor";
2117 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2119 AliReconstructor* reconstructor = NULL;
2120 // first check if a plugin is defined for the reconstructor
2121 TPluginHandler* pluginHandler =
2122 pluginManager->FindHandler("AliReconstructor", detName);
2123 // if not, add a plugin for it
2124 if (!pluginHandler) {
2125 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2126 TString libs = gSystem->GetLibraries();
2127 if (libs.Contains("lib" + detName + "base.so") ||
2128 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2129 pluginManager->AddHandler("AliReconstructor", detName,
2130 recName, detName + "rec", recName + "()");
2132 pluginManager->AddHandler("AliReconstructor", detName,
2133 recName, detName, recName + "()");
2135 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2137 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2138 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2140 if (reconstructor) {
2141 TObject* obj = fOptions.FindObject(detName.Data());
2142 if (obj) reconstructor->SetOption(obj->GetTitle());
2143 reconstructor->Init();
2144 fReconstructor[iDet] = reconstructor;
2147 // get or create the loader
2148 if (detName != "HLT") {
2149 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2150 if (!fLoader[iDet]) {
2151 AliConfig::Instance()
2152 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2154 // first check if a plugin is defined for the loader
2156 pluginManager->FindHandler("AliLoader", detName);
2157 // if not, add a plugin for it
2158 if (!pluginHandler) {
2159 TString loaderName = "Ali" + detName + "Loader";
2160 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2161 pluginManager->AddHandler("AliLoader", detName,
2162 loaderName, detName + "base",
2163 loaderName + "(const char*, TFolder*)");
2164 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2166 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2168 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2169 fRunLoader->GetEventFolder());
2171 if (!fLoader[iDet]) { // use default loader
2172 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2174 if (!fLoader[iDet]) {
2175 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2176 if (fStopOnError) return NULL;
2178 fRunLoader->AddLoader(fLoader[iDet]);
2179 fRunLoader->CdGAFile();
2180 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2181 fRunLoader->Write(0, TObject::kOverwrite);
2186 return reconstructor;
2189 //_____________________________________________________________________________
2190 Bool_t AliReconstruction::CreateVertexer()
2192 // create the vertexer
2195 AliReconstructor* itsReconstructor = GetReconstructor(0);
2196 if (itsReconstructor) {
2197 fVertexer = itsReconstructor->CreateVertexer();
2200 AliWarning("couldn't create a vertexer for ITS");
2201 if (fStopOnError) return kFALSE;
2207 //_____________________________________________________________________________
2208 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2210 // create the trackers
2212 TString detStr = detectors;
2213 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2214 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2215 AliReconstructor* reconstructor = GetReconstructor(iDet);
2216 if (!reconstructor) continue;
2217 TString detName = fgkDetectorName[iDet];
2218 if (detName == "HLT") {
2219 fRunHLTTracking = kTRUE;
2222 if (detName == "MUON") {
2223 fRunMuonTracking = kTRUE;
2228 fTracker[iDet] = reconstructor->CreateTracker();
2229 if (!fTracker[iDet] && (iDet < 7)) {
2230 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2231 if (fStopOnError) return kFALSE;
2233 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2239 //_____________________________________________________________________________
2240 void AliReconstruction::CleanUp(TFile* file, TFile* fileOld)
2242 // delete trackers and the run loader and close and delete the file
2244 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2245 delete fReconstructor[iDet];
2246 fReconstructor[iDet] = NULL;
2247 fLoader[iDet] = NULL;
2248 delete fTracker[iDet];
2249 fTracker[iDet] = NULL;
2250 // delete fQADataMaker[iDet];
2251 // fQADataMaker[iDet] = NULL;
2256 if (ftVertexer) delete ftVertexer;
2259 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2260 delete fDiamondProfile;
2261 fDiamondProfile = NULL;
2262 delete fDiamondProfileTPC;
2263 fDiamondProfileTPC = NULL;
2273 if (fParentRawReader) delete fParentRawReader;
2274 fParentRawReader=NULL;
2284 gSystem->Unlink("AliESDs.old.root");
2289 //_____________________________________________________________________________
2291 Bool_t AliReconstruction::ReadESD(AliESDEvent*& esd, const char* recStep) const
2293 // read the ESD event from a file
2295 if (!esd) return kFALSE;
2297 sprintf(fileName, "ESD_%d.%d_%s.root",
2298 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2299 if (gSystem->AccessPathName(fileName)) return kFALSE;
2301 AliInfo(Form("reading ESD from file %s", fileName));
2302 AliDebug(1, Form("reading ESD from file %s", fileName));
2303 TFile* file = TFile::Open(fileName);
2304 if (!file || !file->IsOpen()) {
2305 AliError(Form("opening %s failed", fileName));
2312 esd = (AliESDEvent*) file->Get("ESD");
2321 //_____________________________________________________________________________
2322 void AliReconstruction::WriteESD(AliESDEvent* esd, const char* recStep) const
2324 // write the ESD event to a file
2328 sprintf(fileName, "ESD_%d.%d_%s.root",
2329 esd->GetRunNumber(), esd->GetEventNumberInFile(), recStep);
2331 AliDebug(1, Form("writing ESD to file %s", fileName));
2332 TFile* file = TFile::Open(fileName, "recreate");
2333 if (!file || !file->IsOpen()) {
2334 AliError(Form("opening %s failed", fileName));
2343 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2345 // Write space-points which are then used in the alignment procedures
2346 // For the moment only ITS, TRD and TPC
2348 // Load TOF clusters
2350 fLoader[3]->LoadRecPoints("read");
2351 TTree* tree = fLoader[3]->TreeR();
2353 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[3]));
2356 fTracker[3]->LoadClusters(tree);
2358 Int_t ntracks = esd->GetNumberOfTracks();
2359 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2361 AliESDtrack *track = esd->GetTrack(itrack);
2364 for (Int_t iDet = 3; iDet >= 0; iDet--)
2365 nsp += track->GetNcls(iDet);
2367 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2368 track->SetTrackPointArray(sp);
2370 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2371 AliTracker *tracker = fTracker[iDet];
2372 if (!tracker) continue;
2373 Int_t nspdet = track->GetNcls(iDet);
2374 if (nspdet <= 0) continue;
2375 track->GetClusters(iDet,idx);
2379 while (isp2 < nspdet) {
2381 TString dets = fgkDetectorName[iDet];
2382 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2383 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2384 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2385 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2386 isvalid = tracker->GetTrackPointTrackingError(idx[isp2],p,track);
2388 isvalid = tracker->GetTrackPoint(idx[isp2],p);
2391 const Int_t kNTPCmax = 159;
2392 if (iDet==1 && isp2>kNTPCmax) break; // to be fixed
2393 if (!isvalid) continue;
2394 sp->AddPoint(isptrack,&p); isptrack++; isp++;
2400 fTracker[3]->UnloadClusters();
2401 fLoader[3]->UnloadRecPoints();
2405 //_____________________________________________________________________________
2406 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2408 // The method reads the raw-data error log
2409 // accumulated within the rawReader.
2410 // It extracts the raw-data errors related to
2411 // the current event and stores them into
2412 // a TClonesArray inside the esd object.
2414 if (!fRawReader) return;
2416 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2418 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2420 if (iEvent != log->GetEventNumber()) continue;
2422 esd->AddRawDataErrorLog(log);
2427 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2428 // Dump a file content into a char in TNamed
2430 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2431 Int_t kBytes = (Int_t)in.tellg();
2432 printf("Size: %d \n",kBytes);
2435 char* memblock = new char [kBytes];
2436 in.seekg (0, ios::beg);
2437 in.read (memblock, kBytes);
2439 TString fData(memblock,kBytes);
2440 fn = new TNamed(pName,fData);
2441 printf("fData Size: %d \n",fData.Sizeof());
2442 printf("pName Size: %d \n",pName.Sizeof());
2443 printf("fn Size: %d \n",fn->Sizeof());
2447 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2453 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2454 // This is not really needed in AliReconstruction at the moment
2455 // but can serve as a template
2457 TList *fList = fTree->GetUserInfo();
2458 TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2459 printf("fn Size: %d \n",fn->Sizeof());
2461 TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2462 const char* cdata = fn->GetTitle();
2463 printf("fTmp Size %d\n",fTmp.Sizeof());
2465 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2466 printf("calculated size %d\n",size);
2467 ofstream out(pName.Data(),ios::out | ios::binary);
2468 out.write(cdata,size);
2473 //_____________________________________________________________________________
2474 AliQADataMakerRec * AliReconstruction::GetQADataMaker(Int_t iDet)
2476 // get the quality assurance data maker object and the loader for a detector
2478 if (fQADataMaker[iDet])
2479 return fQADataMaker[iDet];
2481 AliQADataMakerRec * qadm = NULL;
2482 if (iDet == fgkNDetectors) { //Global QA
2483 qadm = new AliGlobalQADataMaker();
2484 fQADataMaker[iDet] = qadm;
2488 // load the QA data maker object
2489 TPluginManager* pluginManager = gROOT->GetPluginManager();
2490 TString detName = fgkDetectorName[iDet];
2491 TString qadmName = "Ali" + detName + "QADataMakerRec";
2492 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT"))
2495 // first check if a plugin is defined for the quality assurance data maker
2496 TPluginHandler* pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2497 // if not, add a plugin for it
2498 if (!pluginHandler) {
2499 AliDebug(1, Form("defining plugin for %s", qadmName.Data()));
2500 TString libs = gSystem->GetLibraries();
2501 if (libs.Contains("lib" + detName + "base.so") ||
2502 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2503 pluginManager->AddHandler("AliQADataMakerRec", detName,
2504 qadmName, detName + "qadm", qadmName + "()");
2506 pluginManager->AddHandler("AliQADataMakerRec", detName,
2507 qadmName, detName, qadmName + "()");
2509 pluginHandler = pluginManager->FindHandler("AliQADataMakerRec", detName);
2511 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2512 qadm = (AliQADataMakerRec *) pluginHandler->ExecPlugin(0);
2515 fQADataMaker[iDet] = qadm;
2520 //_____________________________________________________________________________
2521 Bool_t AliReconstruction::RunQA(const char* detectors, AliESDEvent *& esd)
2523 // run the Quality Assurance data producer
2525 AliCodeTimerAuto("")
2526 TString detStr = detectors;
2527 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2528 if (!IsSelected(fgkDetectorName[iDet], detStr))
2530 AliQADataMakerRec * qadm = GetQADataMaker(iDet);
2533 AliCodeTimerStart(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2534 AliInfo(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2536 qadm->Exec(AliQA::kESDS, esd) ;
2539 AliCodeTimerStop(Form("running quality assurance data maker for %s", fgkDetectorName[iDet]));
2541 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2542 AliError(Form("the following detectors were not found: %s",
2552 //_____________________________________________________________________________
2553 void AliReconstruction::CheckQA()
2555 // check the QA of SIM for this run and remove the detectors
2556 // with status Fatal
2558 TString newRunLocalReconstruction ;
2559 TString newRunTracking ;
2560 TString newFillESD ;
2562 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2563 TString detName(AliQA::GetDetName(iDet)) ;
2564 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2565 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2566 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2568 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2569 fRunLocalReconstruction.Contains("ALL") ) {
2570 newRunLocalReconstruction += detName ;
2571 newRunLocalReconstruction += " " ;
2573 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2574 fRunTracking.Contains("ALL") ) {
2575 newRunTracking += detName ;
2576 newRunTracking += " " ;
2578 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2579 fFillESD.Contains("ALL") ) {
2580 newFillESD += detName ;
2585 fRunLocalReconstruction = newRunLocalReconstruction ;
2586 fRunTracking = newRunTracking ;
2587 fFillESD = newFillESD ;
2590 //_____________________________________________________________________________
2591 Int_t AliReconstruction::GetDetIndex(const char* detector)
2593 // return the detector index corresponding to detector
2595 for (index = 0; index < fgkNDetectors ; index++) {
2596 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2601 //_____________________________________________________________________________
2602 Bool_t AliReconstruction::FinishPlaneEff() {
2604 // Here execute all the necessary operationis, at the end of the tracking phase,
2605 // in case that evaluation of PlaneEfficiencies was required for some detector.
2606 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2608 // This Preliminary version works only FOR ITS !!!!!
2609 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2612 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2615 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2616 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2617 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2618 if(fTracker[iDet]) {
2619 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2620 ret=planeeff->WriteIntoCDB();
2621 if(planeeff->GetCreateHistos()) {
2622 TString name="PlaneEffHisto";
2623 name+=fgkDetectorName[iDet];
2625 ret*=planeeff->WriteHistosToFile(name,"RECREATE");
2631 //_____________________________________________________________________________
2632 Bool_t AliReconstruction::InitPlaneEff() {
2634 // Here execute all the necessary operations, before of the tracking phase,
2635 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2636 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2637 // which should be updated/recalculated.
2639 // This Preliminary version will work only FOR ITS !!!!!
2640 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2643 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2645 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2649 //_____________________________________________________________________________
2650 Bool_t AliReconstruction::InitAliEVE()
2652 // This method should be called only in case
2653 // AliReconstruction is run
2654 // within the alieve environment.
2655 // It will initialize AliEVE in a way
2656 // so that it can visualize event processed
2657 // by AliReconstruction.
2658 // The return flag shows whenever the
2659 // AliEVE initialization was successful or not.
2662 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2663 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2664 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2666 gROOT->ProcessLine("if (!gAliEveEvent) {gAliEveEvent = new AliEveEventManager();gAliEveEvent->SetAutoLoad(kTRUE);gAliEveEvent->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(gAliEveEvent);};");
2667 gROOT->ProcessLine("alieve_online_init()");
2672 //_____________________________________________________________________________
2673 void AliReconstruction::RunAliEVE()
2675 // Runs AliEVE visualisation of
2676 // the current event.
2677 // Should be executed only after
2678 // successful initialization of AliEVE.
2680 AliInfo("Running AliEVE...");
2681 gROOT->ProcessLine(Form("gAliEveEvent->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2682 gROOT->ProcessLine("gAliEveEvent->StartStopAutoLoadTimer();");