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("..."); //
109 ///////////////////////////////////////////////////////////////////////////////
116 #include <TPluginManager.h>
117 #include <TGeoManager.h>
118 #include <TLorentzVector.h>
121 #include <TObjArray.h>
125 #include <TProofOutputFile.h>
126 #include <TParameter.h>
128 #include "AliReconstruction.h"
129 #include "AliCodeTimer.h"
130 #include "AliReconstructor.h"
132 #include "AliRunLoader.h"
134 #include "AliRawReaderFile.h"
135 #include "AliRawReaderDate.h"
136 #include "AliRawReaderRoot.h"
137 #include "AliRawEventHeaderBase.h"
138 #include "AliRawEvent.h"
139 #include "AliESDEvent.h"
140 #include "AliESDMuonTrack.h"
141 #include "AliESDfriend.h"
142 #include "AliESDVertex.h"
143 #include "AliESDcascade.h"
144 #include "AliESDkink.h"
145 #include "AliESDtrack.h"
146 #include "AliESDCaloCluster.h"
147 #include "AliESDCaloCells.h"
148 #include "AliMultiplicity.h"
149 #include "AliTracker.h"
150 #include "AliVertexer.h"
151 #include "AliVertexerTracks.h"
152 #include "AliV0vertexer.h"
153 #include "AliCascadeVertexer.h"
154 #include "AliHeader.h"
155 #include "AliGenEventHeader.h"
157 #include "AliESDpid.h"
158 #include "AliESDtrack.h"
159 #include "AliESDPmdTrack.h"
161 #include "AliESDTagCreator.h"
163 #include "AliGeomManager.h"
164 #include "AliTrackPointArray.h"
165 #include "AliCDBManager.h"
166 #include "AliCDBStorage.h"
167 #include "AliCDBEntry.h"
168 #include "AliAlignObj.h"
170 #include "AliCentralTrigger.h"
171 #include "AliTriggerConfiguration.h"
172 #include "AliTriggerClass.h"
173 #include "AliTriggerCluster.h"
174 #include "AliCTPRawStream.h"
176 #include "AliQADataMakerRec.h"
177 #include "AliGlobalQADataMaker.h"
179 #include "AliQADataMakerSteer.h"
181 #include "AliPlaneEff.h"
183 #include "AliSysInfo.h" // memory snapshots
184 #include "AliRawHLTManager.h"
186 #include "AliMagWrapCheb.h"
188 #include "AliDetectorRecoParam.h"
189 #include "AliGRPRecoParam.h"
190 #include "AliRunInfo.h"
191 #include "AliEventInfo.h"
195 #include "AliGRPObject.h"
197 ClassImp(AliReconstruction)
199 //_____________________________________________________________________________
200 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
202 //_____________________________________________________________________________
203 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
205 fUniformField(kFALSE),
206 fForcedFieldMap(NULL),
207 fRunVertexFinder(kTRUE),
208 fRunVertexFinderTracks(kTRUE),
209 fRunHLTTracking(kFALSE),
210 fRunMuonTracking(kFALSE),
212 fRunCascadeFinder(kTRUE),
213 fStopOnError(kFALSE),
214 fWriteAlignmentData(kFALSE),
215 fWriteESDfriend(kFALSE),
216 fFillTriggerESD(kTRUE),
224 fRunLocalReconstruction("ALL"),
228 fUseTrackingErrorsForAlignment(""),
229 fGAliceFileName(gAliceFilename),
234 fNumberOfEventsPerFile(1),
236 fLoadAlignFromCDB(kTRUE),
237 fLoadAlignData("ALL"),
245 fParentRawReader(NULL),
250 fDiamondProfileSPD(NULL),
251 fDiamondProfile(NULL),
252 fDiamondProfileTPC(NULL),
256 fAlignObjArray(NULL),
259 fInitCDBCalled(kFALSE),
260 fSetRunNumberFromDataCalled(kFALSE),
266 fSameQACycle(kFALSE),
268 fRunPlaneEff(kFALSE),
277 fIsNewRunLoader(kFALSE),
281 // create reconstruction object with default parameters
284 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
285 fReconstructor[iDet] = NULL;
286 fLoader[iDet] = NULL;
287 fTracker[iDet] = NULL;
289 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
290 fQACycles[iDet] = 999999 ;
291 fQAWriteExpert[iDet] = kFALSE ;
297 //_____________________________________________________________________________
298 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
300 fUniformField(rec.fUniformField),
301 fForcedFieldMap(NULL),
302 fRunVertexFinder(rec.fRunVertexFinder),
303 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
304 fRunHLTTracking(rec.fRunHLTTracking),
305 fRunMuonTracking(rec.fRunMuonTracking),
306 fRunV0Finder(rec.fRunV0Finder),
307 fRunCascadeFinder(rec.fRunCascadeFinder),
308 fStopOnError(rec.fStopOnError),
309 fWriteAlignmentData(rec.fWriteAlignmentData),
310 fWriteESDfriend(rec.fWriteESDfriend),
311 fFillTriggerESD(rec.fFillTriggerESD),
313 fCleanESD(rec.fCleanESD),
314 fV0DCAmax(rec.fV0DCAmax),
315 fV0CsPmin(rec.fV0CsPmin),
319 fRunLocalReconstruction(rec.fRunLocalReconstruction),
320 fRunTracking(rec.fRunTracking),
321 fFillESD(rec.fFillESD),
322 fLoadCDB(rec.fLoadCDB),
323 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
324 fGAliceFileName(rec.fGAliceFileName),
325 fRawInput(rec.fRawInput),
326 fEquipIdMap(rec.fEquipIdMap),
327 fFirstEvent(rec.fFirstEvent),
328 fLastEvent(rec.fLastEvent),
329 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
331 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
332 fLoadAlignData(rec.fLoadAlignData),
333 fESDPar(rec.fESDPar),
334 fUseHLTData(rec.fUseHLTData),
340 fParentRawReader(NULL),
342 fRecoParam(rec.fRecoParam),
345 fDiamondProfileSPD(rec.fDiamondProfileSPD),
346 fDiamondProfile(rec.fDiamondProfile),
347 fDiamondProfileTPC(rec.fDiamondProfileTPC),
351 fAlignObjArray(rec.fAlignObjArray),
352 fCDBUri(rec.fCDBUri),
354 fInitCDBCalled(rec.fInitCDBCalled),
355 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
356 fQADetectors(rec.fQADetectors),
358 fQATasks(rec.fQATasks),
360 fRunGlobalQA(rec.fRunGlobalQA),
361 fSameQACycle(rec.fSameQACycle),
362 fRunPlaneEff(rec.fRunPlaneEff),
371 fIsNewRunLoader(rec.fIsNewRunLoader),
377 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
378 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
380 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
381 fReconstructor[iDet] = NULL;
382 fLoader[iDet] = NULL;
383 fTracker[iDet] = NULL;
386 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
387 fQACycles[iDet] = rec.fQACycles[iDet];
388 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
391 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
392 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
396 //_____________________________________________________________________________
397 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
399 // assignment operator
400 // Used in PROOF mode
401 // Be very careful while modifing it!
402 // Simple rules to follow:
403 // for persistent data members - use their assignment operators
404 // for non-persistent ones - do nothing or take the default values from constructor
405 // TSelector members should not be touched
406 if(&rec == this) return *this;
408 fUniformField = rec.fUniformField;
409 fForcedFieldMap = NULL;
410 fRunVertexFinder = rec.fRunVertexFinder;
411 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
412 fRunHLTTracking = rec.fRunHLTTracking;
413 fRunMuonTracking = rec.fRunMuonTracking;
414 fRunV0Finder = rec.fRunV0Finder;
415 fRunCascadeFinder = rec.fRunCascadeFinder;
416 fStopOnError = rec.fStopOnError;
417 fWriteAlignmentData = rec.fWriteAlignmentData;
418 fWriteESDfriend = rec.fWriteESDfriend;
419 fFillTriggerESD = rec.fFillTriggerESD;
421 fCleanESD = rec.fCleanESD;
422 fV0DCAmax = rec.fV0DCAmax;
423 fV0CsPmin = rec.fV0CsPmin;
427 fRunLocalReconstruction = rec.fRunLocalReconstruction;
428 fRunTracking = rec.fRunTracking;
429 fFillESD = rec.fFillESD;
430 fLoadCDB = rec.fLoadCDB;
431 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
432 fGAliceFileName = rec.fGAliceFileName;
433 fRawInput = rec.fRawInput;
434 fEquipIdMap = rec.fEquipIdMap;
435 fFirstEvent = rec.fFirstEvent;
436 fLastEvent = rec.fLastEvent;
437 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
439 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
440 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
443 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
444 fLoadAlignData = rec.fLoadAlignData;
445 fESDPar = rec.fESDPar;
446 fUseHLTData = rec.fUseHLTData;
448 delete fRunInfo; fRunInfo = NULL;
449 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
451 fEventInfo = rec.fEventInfo;
455 fParentRawReader = NULL;
457 fRecoParam = rec.fRecoParam;
459 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
460 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
461 delete fLoader[iDet]; fLoader[iDet] = NULL;
462 delete fTracker[iDet]; fTracker[iDet] = NULL;
465 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
466 fQACycles[iDet] = rec.fQACycles[iDet];
467 fQAWriteExpert[iDet] = rec.fQAWriteExpert[iDet] ;
471 delete fDiamondProfileSPD; fDiamondProfileSPD = NULL;
472 if (rec.fDiamondProfileSPD) fDiamondProfileSPD = new AliESDVertex(*rec.fDiamondProfileSPD);
473 delete fDiamondProfile; fDiamondProfile = NULL;
474 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
475 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
476 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
478 delete fGRPData; fGRPData = NULL;
479 // if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
480 if (rec.fGRPData) fGRPData = (AliGRPObject*)((rec.fGRPData)->Clone());
482 delete fAlignObjArray; fAlignObjArray = NULL;
485 fSpecCDBUri.Delete();
486 fInitCDBCalled = rec.fInitCDBCalled;
487 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
488 fQADetectors = rec.fQADetectors;
490 fQATasks = rec.fQATasks;
492 fRunGlobalQA = rec.fRunGlobalQA;
493 fSameQACycle = rec.fSameQACycle;
494 fRunPlaneEff = rec.fRunPlaneEff;
503 fIsNewRunLoader = rec.fIsNewRunLoader;
510 //_____________________________________________________________________________
511 AliReconstruction::~AliReconstruction()
517 delete fForcedFieldMap;
519 if (fAlignObjArray) {
520 fAlignObjArray->Delete();
521 delete fAlignObjArray;
523 fSpecCDBUri.Delete();
525 AliCodeTimer::Instance()->Print();
528 //_____________________________________________________________________________
529 void AliReconstruction::InitCDB()
531 // activate a default CDB storage
532 // First check if we have any CDB storage set, because it is used
533 // to retrieve the calibration and alignment constants
534 AliCodeTimerAuto("");
536 if (fInitCDBCalled) return;
537 fInitCDBCalled = kTRUE;
539 AliCDBManager* man = AliCDBManager::Instance();
540 if (man->IsDefaultStorageSet())
542 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
543 AliWarning("Default CDB storage has been already set !");
544 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
545 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
546 fCDBUri = man->GetDefaultStorage()->GetURI();
549 if (fCDBUri.Length() > 0)
551 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
552 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
553 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
555 fCDBUri="local://$ALICE_ROOT";
556 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
557 AliWarning("Default CDB storage not yet set !!!!");
558 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
559 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
562 man->SetDefaultStorage(fCDBUri);
565 // Now activate the detector specific CDB storage locations
566 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
567 TObject* obj = fSpecCDBUri[i];
569 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
570 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
571 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
572 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
574 AliSysInfo::AddStamp("InitCDB");
577 //_____________________________________________________________________________
578 void AliReconstruction::SetDefaultStorage(const char* uri) {
579 // Store the desired default CDB storage location
580 // Activate it later within the Run() method
586 //_____________________________________________________________________________
587 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
588 // Store a detector-specific CDB storage location
589 // Activate it later within the Run() method
591 AliCDBPath aPath(calibType);
592 if(!aPath.IsValid()){
593 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
594 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
595 if(!strcmp(calibType, fgkDetectorName[iDet])) {
596 aPath.SetPath(Form("%s/*", calibType));
597 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
601 if(!aPath.IsValid()){
602 AliError(Form("Not a valid path or detector: %s", calibType));
607 // // check that calibType refers to a "valid" detector name
608 // Bool_t isDetector = kFALSE;
609 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
610 // TString detName = fgkDetectorName[iDet];
611 // if(aPath.GetLevel0() == detName) {
612 // isDetector = kTRUE;
618 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
622 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
623 if (obj) fSpecCDBUri.Remove(obj);
624 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
628 //_____________________________________________________________________________
629 Bool_t AliReconstruction::SetRunNumberFromData()
631 // The method is called in Run() in order
632 // to set a correct run number.
633 // In case of raw data reconstruction the
634 // run number is taken from the raw data header
636 if (fSetRunNumberFromDataCalled) return kTRUE;
637 fSetRunNumberFromDataCalled = kTRUE;
639 AliCDBManager* man = AliCDBManager::Instance();
642 if(fRawReader->NextEvent()) {
643 if(man->GetRun() > 0) {
644 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
646 man->SetRun(fRawReader->GetRunNumber());
647 fRawReader->RewindEvents();
650 if(man->GetRun() > 0) {
651 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
654 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
660 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
662 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
667 // read run number from gAlice
668 if(rl->GetHeader()) {
669 man->SetRun(rl->GetHeader()->GetRun());
674 AliError("Neither run-loader header nor RawReader objects are found !");
686 //_____________________________________________________________________________
687 void AliReconstruction::SetCDBLock() {
688 // Set CDB lock: from now on it is forbidden to reset the run number
689 // or the default storage or to activate any further storage!
691 AliCDBManager::Instance()->SetLock(1);
694 //_____________________________________________________________________________
695 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
697 // Read the alignment objects from CDB.
698 // Each detector is supposed to have the
699 // alignment objects in DET/Align/Data CDB path.
700 // All the detector objects are then collected,
701 // sorted by geometry level (starting from ALIC) and
702 // then applied to the TGeo geometry.
703 // Finally an overlaps check is performed.
705 // Load alignment data from CDB and fill fAlignObjArray
706 if(fLoadAlignFromCDB){
708 TString detStr = detectors;
709 TString loadAlObjsListOfDets = "";
711 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
712 if(!IsSelected(fgkDetectorName[iDet], detStr)) continue;
713 if(!strcmp(fgkDetectorName[iDet],"HLT")) continue;
714 if(AliGeomManager::IsModuleInGeom(fgkDetectorName[iDet]))
716 loadAlObjsListOfDets += fgkDetectorName[iDet];
717 loadAlObjsListOfDets += " ";
719 } // end loop over detectors
720 if(AliGeomManager::IsModuleInGeom("FRAME"))
721 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
722 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
723 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
725 // Check if the array with alignment objects was
726 // provided by the user. If yes, apply the objects
727 // to the present TGeo geometry
728 if (fAlignObjArray) {
729 if (gGeoManager && gGeoManager->IsClosed()) {
730 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
731 AliError("The misalignment of one or more volumes failed!"
732 "Compare the list of simulated detectors and the list of detector alignment data!");
737 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
743 if (fAlignObjArray) {
744 fAlignObjArray->Delete();
745 delete fAlignObjArray; fAlignObjArray=NULL;
751 //_____________________________________________________________________________
752 void AliReconstruction::SetGAliceFile(const char* fileName)
754 // set the name of the galice file
756 fGAliceFileName = fileName;
759 //_____________________________________________________________________________
760 void AliReconstruction::SetInput(const char* input)
762 // In case the input string starts with 'mem://', we run in an online mode
763 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
764 // file is assumed. One can give as an input:
765 // mem://: - events taken from DAQ monitoring libs online
767 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
768 if (input) fRawInput = input;
771 //_____________________________________________________________________________
772 void AliReconstruction::SetOption(const char* detector, const char* option)
774 // set options for the reconstruction of a detector
776 TObject* obj = fOptions.FindObject(detector);
777 if (obj) fOptions.Remove(obj);
778 fOptions.Add(new TNamed(detector, option));
781 //_____________________________________________________________________________
782 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
784 // Set custom reconstruction parameters for a given detector
785 // Single set of parameters for all the events
787 // First check if the reco-params are global
788 if(!strcmp(detector, "GRP")) {
790 fRecoParam.AddDetRecoParam(fgkNDetectors,par);
794 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
795 if(!strcmp(detector, fgkDetectorName[iDet])) {
797 fRecoParam.AddDetRecoParam(iDet,par);
804 //_____________________________________________________________________________
805 Bool_t AliReconstruction::SetFieldMap(Float_t l3Current, Float_t diCurrent, Float_t factor, const char *path) {
806 //------------------------------------------------
807 // The magnetic field map, defined externally...
808 // L3 current 30000 A -> 0.5 T
809 // L3 current 12000 A -> 0.2 T
810 // dipole current 6000 A
811 // The polarities must be the same
812 //------------------------------------------------
813 const Float_t l3NominalCurrent1=30000.; // (A)
814 const Float_t l3NominalCurrent2=12000.; // (A)
815 const Float_t diNominalCurrent =6000. ; // (A)
817 const Float_t tolerance=0.03; // relative current tolerance
818 const Float_t zero=77.; // "zero" current (A)
821 Bool_t dipoleON=kFALSE;
823 TString s=(factor < 0) ? "L3: -" : "L3: +";
825 l3Current = TMath::Abs(l3Current);
826 if (TMath::Abs(l3Current-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
827 map=AliMagWrapCheb::k5kG;
830 if (TMath::Abs(l3Current-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
831 map=AliMagWrapCheb::k2kG;
834 if (l3Current < zero) {
835 map=AliMagWrapCheb::k2kG;
837 factor=0.; // in fact, this is a global factor...
838 fUniformField=kTRUE; // track with the uniform (zero) B field
840 AliError(Form("Wrong L3 current (%f A)!",l3Current));
844 diCurrent = TMath::Abs(diCurrent);
845 if (TMath::Abs(diCurrent-diNominalCurrent)/diNominalCurrent < tolerance) {
846 // 3% current tolerance...
850 if (diCurrent < zero) { // some small current..
854 AliError(Form("Wrong dipole current (%f A)!",diCurrent));
858 delete fForcedFieldMap;
860 new AliMagWrapCheb("B field map ",s,2,factor,10.,map,dipoleON,path);
862 fForcedFieldMap->Print();
864 AliTracker::SetFieldMap(fForcedFieldMap,fUniformField);
870 Bool_t AliReconstruction::InitGRP() {
871 //------------------------------------
872 // Initialization of the GRP entry
873 //------------------------------------
874 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
878 TMap* m = dynamic_cast<TMap*>(entry->GetObject()); // old GRP entry
881 AliInfo("Found a TMap in GRP/GRP/Data, converting it into an AliGRPObject");
883 fGRPData = new AliGRPObject();
884 fGRPData->ReadValuesFromMap(m);
888 AliInfo("Found an AliGRPObject in GRP/GRP/Data, reading it");
889 fGRPData = dynamic_cast<AliGRPObject*>(entry->GetObject()); // new GRP entry
893 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
897 AliError("No GRP entry found in OCDB!");
901 TString lhcState = fGRPData->GetLHCState();
902 if (lhcState==AliGRPObject::GetInvalidString()) {
903 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
904 lhcState = "UNKNOWN";
907 TString beamType = fGRPData->GetBeamType();
908 if (beamType==AliGRPObject::GetInvalidString()) {
909 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
910 beamType = "UNKNOWN";
913 Float_t beamEnergy = fGRPData->GetBeamEnergy();
914 if (beamEnergy==AliGRPObject::GetInvalidFloat()) {
915 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
919 TString runType = fGRPData->GetRunType();
920 if (runType==AliGRPObject::GetInvalidString()) {
921 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
925 Int_t activeDetectors = fGRPData->GetDetectorMask();
926 if (activeDetectors==AliGRPObject::GetInvalidUInt()) {
927 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
928 activeDetectors = 1074790399;
931 fRunInfo = new AliRunInfo(lhcState, beamType, beamEnergy, runType, activeDetectors);
936 // Process the list of active detectors
937 if (activeDetectors) {
938 UInt_t detMask = activeDetectors;
939 fLoadCDB.Form("%s %s %s %s",
940 fRunLocalReconstruction.Data(),
943 fQADetectors.Data());
944 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
945 fRunTracking = MatchDetectorList(fRunTracking,detMask);
946 fFillESD = MatchDetectorList(fFillESD,detMask);
947 fQADetectors = MatchDetectorList(fQADetectors,detMask);
948 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
949 if (!((detMask >> AliDAQ::DetectorID("ITSSPD")) & 0x1)) {
950 // switch off the vertexer
951 AliInfo("SPD is not in the list of active detectors. Vertexer switched off.");
952 fRunVertexFinder = kFALSE;
956 AliInfo("===================================================================================");
957 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
958 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
959 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
960 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
961 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
962 AliInfo("===================================================================================");
964 //*** Dealing with the magnetic field map
965 if (AliTracker::GetFieldMap()) {
966 AliInfo("Running with the externally set B field !");
968 // Construct the field map out of the information retrieved from GRP.
973 Float_t l3Current = fGRPData->GetL3Current((AliGRPObject::Stats)0);
974 if (l3Current == AliGRPObject::GetInvalidFloat()) {
975 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
979 Char_t l3Polarity = fGRPData->GetL3Polarity();
980 if (l3Polarity == AliGRPObject::GetInvalidChar()) {
981 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
986 Float_t diCurrent = fGRPData->GetDipoleCurrent((AliGRPObject::Stats)0);
987 if (diCurrent == AliGRPObject::GetInvalidFloat()) {
988 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
992 Char_t diPolarity = fGRPData->GetDipolePolarity();
993 if (diPolarity == AliGRPObject::GetInvalidChar()) {
994 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
999 TObjString *l3Current=
1000 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
1002 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
1005 TObjString *l3Polarity=
1006 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
1008 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
1013 TObjString *diCurrent=
1014 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
1016 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
1019 TObjString *diPolarity=
1020 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
1022 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
1028 Float_t l3Cur=TMath::Abs(l3Current);
1029 Float_t diCur=TMath::Abs(diCurrent);
1030 Float_t l3Pol=l3Polarity;
1031 // Float_t l3Cur=TMath::Abs(atof(l3Current->GetName()));
1032 //Float_t diCur=TMath::Abs(atof(diCurrent->GetName()));
1033 //Float_t l3Pol=atof(l3Polarity->GetName());
1035 if (l3Pol != 0.) factor=-1.;
1038 if (!SetFieldMap(l3Cur, diCur, factor)) {
1039 AliFatal("Failed to creat a B field map ! Exiting...");
1041 AliInfo("Running with the B field constructed out of GRP !");
1044 AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
1049 //*** Get the diamond profiles from OCDB
1050 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
1052 fDiamondProfileSPD = dynamic_cast<AliESDVertex*> (entry->GetObject());
1054 AliError("No SPD diamond profile found in OCDB!");
1057 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
1059 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
1061 AliError("No diamond profile found in OCDB!");
1064 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
1066 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
1068 AliError("No TPC diamond profile found in OCDB!");
1074 //_____________________________________________________________________________
1075 Bool_t AliReconstruction::LoadCDB()
1077 AliCodeTimerAuto("");
1079 AliCDBManager::Instance()->Get("GRP/CTP/Config");
1081 TString detStr = fLoadCDB;
1082 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1083 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1084 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1089 //_____________________________________________________________________________
1090 Bool_t AliReconstruction::Run(const char* input)
1093 AliCodeTimerAuto("");
1096 if (GetAbort() != TSelector::kContinue) return kFALSE;
1098 TChain *chain = NULL;
1099 if (fRawReader && (chain = fRawReader->GetChain())) {
1102 gProof->AddInput(this);
1104 outputFile.SetProtocol("root",kTRUE);
1105 outputFile.SetHost(gSystem->HostName());
1106 outputFile.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1107 AliInfo(Form("Output file with ESDs is %s",outputFile.GetUrl()));
1108 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.GetUrl()));
1110 chain->Process("AliReconstruction");
1113 chain->Process(this);
1118 if (GetAbort() != TSelector::kContinue) return kFALSE;
1120 if (GetAbort() != TSelector::kContinue) return kFALSE;
1121 //******* The loop over events
1122 AliInfo("Starting looping over events");
1124 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1125 (fRawReader && fRawReader->NextEvent())) {
1126 if (!ProcessEvent(iEvent)) {
1127 Abort("ProcessEvent",TSelector::kAbortFile);
1133 if (GetAbort() != TSelector::kContinue) return kFALSE;
1135 if (GetAbort() != TSelector::kContinue) return kFALSE;
1141 //_____________________________________________________________________________
1142 void AliReconstruction::InitRawReader(const char* input)
1144 AliCodeTimerAuto("");
1146 // Init raw-reader and
1147 // set the input in case of raw data
1148 if (input) fRawInput = input;
1149 fRawReader = AliRawReader::Create(fRawInput.Data());
1151 AliInfo("Reconstruction will run over digits");
1153 if (!fEquipIdMap.IsNull() && fRawReader)
1154 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1156 if (!fUseHLTData.IsNull()) {
1157 // create the RawReaderHLT which performs redirection of HLT input data for
1158 // the specified detectors
1159 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1161 fParentRawReader=fRawReader;
1162 fRawReader=pRawReader;
1164 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1167 AliSysInfo::AddStamp("CreateRawReader");
1170 //_____________________________________________________________________________
1171 void AliReconstruction::InitRun(const char* input)
1173 // Initialization of raw-reader,
1174 // run number, CDB etc.
1175 AliCodeTimerAuto("");
1176 AliSysInfo::AddStamp("Start");
1178 // Initialize raw-reader if any
1179 InitRawReader(input);
1181 // Initialize the CDB storage
1184 // Set run number in CDBManager (if it is not already set by the user)
1185 if (!SetRunNumberFromData()) {
1186 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1190 // Set CDB lock: from now on it is forbidden to reset the run number
1191 // or the default storage or to activate any further storage!
1196 //_____________________________________________________________________________
1197 void AliReconstruction::Begin(TTree *)
1199 // Initialize AlReconstruction before
1200 // going into the event loop
1201 // Should follow the TSelector convention
1202 // i.e. initialize only the object on the client side
1203 AliCodeTimerAuto("");
1205 AliReconstruction *reco = NULL;
1207 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1210 AliSysInfo::AddStamp("ReadInputInBegin");
1213 // Import ideal TGeo geometry and apply misalignment
1215 TString geom(gSystem->DirName(fGAliceFileName));
1216 geom += "/geometry.root";
1217 AliGeomManager::LoadGeometry(geom.Data());
1219 Abort("LoadGeometry", TSelector::kAbortProcess);
1222 AliSysInfo::AddStamp("LoadGeom");
1223 TString detsToCheck=fRunLocalReconstruction;
1224 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1225 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1228 AliSysInfo::AddStamp("CheckGeom");
1231 if (!MisalignGeometry(fLoadAlignData)) {
1232 Abort("MisalignGeometry", TSelector::kAbortProcess);
1235 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1236 AliSysInfo::AddStamp("MisalignGeom");
1239 Abort("InitGRP", TSelector::kAbortProcess);
1242 AliSysInfo::AddStamp("InitGRP");
1245 Abort("LoadCDB", TSelector::kAbortProcess);
1248 AliSysInfo::AddStamp("LoadCDB");
1250 // Read the reconstruction parameters from OCDB
1251 if (!InitRecoParams()) {
1252 AliWarning("Not all detectors have correct RecoParam objects initialized");
1254 AliSysInfo::AddStamp("InitRecoParams");
1257 if (reco) *reco = *this;
1258 fInput->Add(gGeoManager);
1260 fInput->Add(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()));
1261 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1262 AliMagF *magFieldMap = (AliMagF*)AliTracker::GetFieldMap();
1263 magFieldMap->SetName("MagneticFieldMap");
1264 fInput->Add(magFieldMap);
1269 //_____________________________________________________________________________
1270 void AliReconstruction::SlaveBegin(TTree*)
1272 // Initialization related to run-loader,
1273 // vertexer, trackers, recontructors
1274 // In proof mode it is executed on the slave
1275 AliCodeTimerAuto("");
1277 TProofOutputFile *outProofFile = NULL;
1279 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1282 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1284 AliGeomManager::SetGeometry(tgeo);
1286 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1287 Int_t runNumber = -1;
1288 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1289 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1290 man->SetCacheFlag(kTRUE);
1291 man->SetLock(kTRUE);
1295 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1296 AliTracker::SetFieldMap(map,fUniformField);
1298 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1299 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1300 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1301 fOutput->Add(outProofFile);
1303 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1306 // get the run loader
1307 if (!InitRunLoader()) {
1308 Abort("InitRunLoader", TSelector::kAbortProcess);
1311 AliSysInfo::AddStamp("LoadLoader");
1313 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1316 if (fRunVertexFinder && !CreateVertexer()) {
1317 Abort("CreateVertexer", TSelector::kAbortProcess);
1320 AliSysInfo::AddStamp("CreateVertexer");
1323 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1324 Abort("CreateTrackers", TSelector::kAbortProcess);
1327 AliSysInfo::AddStamp("CreateTrackers");
1329 // create the ESD output file and tree
1330 if (!outProofFile) {
1331 ffile = TFile::Open("AliESDs.root", "RECREATE");
1332 ffile->SetCompressionLevel(2);
1333 if (!ffile->IsOpen()) {
1334 Abort("OpenESDFile", TSelector::kAbortProcess);
1339 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1340 Abort(Form("Problems opening output PROOF file: %s/%s",
1341 outProofFile->GetDir(), outProofFile->GetFileName()),
1342 TSelector::kAbortProcess);
1347 ftree = new TTree("esdTree", "Tree with ESD objects");
1348 fesd = new AliESDEvent();
1349 fesd->CreateStdContent();
1350 if (fWriteESDfriend) {
1351 fesdf = new AliESDfriend();
1352 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1353 br->SetFile("AliESDfriends.root");
1354 fesd->AddObject(fesdf);
1356 fesd->WriteToTree(ftree);
1357 ftree->GetUserInfo()->Add(fesd);
1359 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1360 fhltesd = new AliESDEvent();
1361 fhltesd->CreateStdContent();
1363 // read the ESD template from CDB
1364 // HLT is allowed to put non-std content to its ESD, the non-std
1365 // objects need to be created before invocation of WriteToTree in
1366 // order to create all branches. Initialization is done from an
1367 // ESD layout template in CDB
1368 AliCDBManager* man = AliCDBManager::Instance();
1369 AliCDBPath hltESDConfigPath("HLT/ConfigHLT/esdLayout");
1370 AliCDBEntry* hltESDConfig=NULL;
1371 if (man->GetId(hltESDConfigPath)!=NULL &&
1372 (hltESDConfig=man->Get(hltESDConfigPath))!=NULL) {
1373 AliESDEvent* pESDLayout=dynamic_cast<AliESDEvent*>(hltESDConfig->GetObject());
1375 // init all internal variables from the list of objects
1376 pESDLayout->GetStdContent();
1378 // copy content and create non-std objects
1379 *fhltesd=*pESDLayout;
1382 AliError(Form("error setting hltEsd layout from %s: invalid object type",
1383 hltESDConfigPath.GetPath().Data()));
1387 fhltesd->WriteToTree(fhlttree);
1388 fhlttree->GetUserInfo()->Add(fhltesd);
1390 ProcInfo_t ProcInfo;
1391 gSystem->GetProcInfo(&ProcInfo);
1392 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1395 //Initialize the QA and start of cycle
1397 fQASteer = new AliQADataMakerSteer("rec") ;
1398 fQASteer->SetActiveDetectors(fQADetectors) ;
1399 for (Int_t det = 0 ; det < AliQA::kNDET ; det++) {
1400 fQASteer->SetCycleLength(AliQA::DETECTORINDEX_t(det), fQACycles[det]) ;
1401 if (fQAWriteExpert[det])
1402 fQASteer->SetWriteExpert(AliQA::DETECTORINDEX_t(det)) ;
1405 if (!fRawReader && fQATasks.Contains(AliQA::kRAWS))
1406 fQATasks.ReplaceAll(Form("%d",AliQA::kRAWS), "") ;
1407 fQASteer->SetTasks(fQATasks) ;
1408 fQASteer->InitQADataMaker(AliCDBManager::Instance()->GetRun(), fRecoParam) ;
1412 Bool_t sameCycle = kFALSE ;
1413 if (!fQASteer) fQASteer = new AliQADataMakerSteer("rec") ;
1414 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1415 AliInfo(Form("Initializing the global QA data maker"));
1416 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1417 qadm->StartOfCycle(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1418 TObjArray *arr=qadm->Init(AliQA::kRECPOINTS);
1419 AliTracker::SetResidualsArray(arr);
1422 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1423 qadm->StartOfCycle(AliQA::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1424 qadm->Init(AliQA::kESDS);
1428 //Initialize the Plane Efficiency framework
1429 if (fRunPlaneEff && !InitPlaneEff()) {
1430 Abort("InitPlaneEff", TSelector::kAbortProcess);
1434 if (strcmp(gProgName,"alieve") == 0)
1435 fRunAliEVE = InitAliEVE();
1440 //_____________________________________________________________________________
1441 Bool_t AliReconstruction::Process(Long64_t entry)
1443 // run the reconstruction over a single entry
1444 // from the chain with raw data
1445 AliCodeTimerAuto("");
1447 TTree *currTree = fChain->GetTree();
1448 AliRawEvent *event = new AliRawEvent;
1449 currTree->SetBranchAddress("rawevent",&event);
1450 currTree->GetEntry(entry);
1451 fRawReader = new AliRawReaderRoot(event);
1452 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1460 //_____________________________________________________________________________
1461 void AliReconstruction::Init(TTree *tree)
1464 AliError("The input tree is not found!");
1470 //_____________________________________________________________________________
1471 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1473 // run the reconstruction over a single event
1474 // The event loop is steered in Run method
1476 AliCodeTimerAuto("");
1478 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1479 fRunLoader->SetEventNumber(iEvent);
1480 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1482 fRunLoader->TreeE()->Fill();
1483 if (fRawReader && fRawReader->UseAutoSaveESD())
1484 fRunLoader->TreeE()->AutoSave("SaveSelf");
1487 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1491 AliInfo(Form("processing event %d", iEvent));
1493 // Fill Event-info object
1495 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo);
1497 // Set the reco-params
1499 TString detStr = fLoadCDB;
1500 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1501 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1502 AliReconstructor *reconstructor = GetReconstructor(iDet);
1503 if (reconstructor && fRecoParam.GetDetRecoParamArray(iDet)) {
1504 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
1505 reconstructor->SetRecoParam(par);
1510 fRunLoader->GetEvent(iEvent);
1514 fQASteer->RunOneEvent(fRawReader) ;
1516 // local single event reconstruction
1517 if (!fRunLocalReconstruction.IsNull()) {
1518 TString detectors=fRunLocalReconstruction;
1519 // run HLT event reconstruction first
1520 // ;-( IsSelected changes the string
1521 if (IsSelected("HLT", detectors) &&
1522 !RunLocalEventReconstruction("HLT")) {
1523 if (fStopOnError) {CleanUp(); return kFALSE;}
1525 detectors=fRunLocalReconstruction;
1526 detectors.ReplaceAll("HLT", "");
1527 if (!RunLocalEventReconstruction(detectors)) {
1528 if (fStopOnError) {CleanUp(); return kFALSE;}
1532 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1533 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1534 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1535 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1537 // Set magnetic field from the tracker
1538 fesd->SetMagneticField(AliTracker::GetBz());
1539 fhltesd->SetMagneticField(AliTracker::GetBz());
1541 // Set most probable pt, for B=0 tracking
1542 // Get the global reco-params. They are atposition 16 inside the array of detectors in fRecoParam
1543 const AliGRPRecoParam *grpRecoParam = dynamic_cast<const AliGRPRecoParam*>(fRecoParam.GetDetRecoParam(fgkNDetectors));
1544 if (grpRecoParam) AliExternalTrackParam::SetMostProbablePt(grpRecoParam->GetMostProbablePt());
1546 // Fill raw-data error log into the ESD
1547 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1550 if (fRunVertexFinder) {
1551 if (!RunVertexFinder(fesd)) {
1552 if (fStopOnError) {CleanUp(); return kFALSE;}
1557 if (!fRunTracking.IsNull()) {
1558 if (fRunMuonTracking) {
1559 if (!RunMuonTracking(fesd)) {
1560 if (fStopOnError) {CleanUp(); return kFALSE;}
1566 if (!fRunTracking.IsNull()) {
1567 if (!RunTracking(fesd)) {
1568 if (fStopOnError) {CleanUp(); return kFALSE;}
1573 if (!fFillESD.IsNull()) {
1574 TString detectors=fFillESD;
1575 // run HLT first and on hltesd
1576 // ;-( IsSelected changes the string
1577 if (IsSelected("HLT", detectors) &&
1578 !FillESD(fhltesd, "HLT")) {
1579 if (fStopOnError) {CleanUp(); return kFALSE;}
1582 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1583 if (detectors.Contains("ALL")) {
1585 for (Int_t idet=0; idet<fgkNDetectors; ++idet){
1586 detectors += fgkDetectorName[idet];
1590 detectors.ReplaceAll("HLT", "");
1591 if (!FillESD(fesd, detectors)) {
1592 if (fStopOnError) {CleanUp(); return kFALSE;}
1596 // fill Event header information from the RawEventHeader
1597 if (fRawReader){FillRawEventHeaderESD(fesd);}
1600 AliESDpid::MakePID(fesd);
1602 if (fFillTriggerESD) {
1603 if (!FillTriggerESD(fesd)) {
1604 if (fStopOnError) {CleanUp(); return kFALSE;}
1611 // Propagate track to the beam pipe (if not already done by ITS)
1613 const Int_t ntracks = fesd->GetNumberOfTracks();
1614 const Double_t kBz = fesd->GetMagneticField();
1615 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1618 UShort_t *selectedIdx=new UShort_t[ntracks];
1620 for (Int_t itrack=0; itrack<ntracks; itrack++){
1621 const Double_t kMaxStep = 5; //max step over the material
1624 AliESDtrack *track = fesd->GetTrack(itrack);
1625 if (!track) continue;
1627 AliExternalTrackParam *tpcTrack =
1628 (AliExternalTrackParam *)track->GetTPCInnerParam();
1632 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kFALSE);
1635 Int_t n=trkArray.GetEntriesFast();
1636 selectedIdx[n]=track->GetID();
1637 trkArray.AddLast(tpcTrack);
1640 //Tracks refitted by ITS should already be at the SPD vertex
1641 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1644 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kFALSE);
1645 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1650 // Improve the reconstructed primary vertex position using the tracks
1652 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1653 if(fesd->GetPrimaryVertexSPD()) {
1654 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1655 if(vtitle.Contains("cosmics")) {
1656 runVertexFinderTracks=kFALSE;
1660 if (runVertexFinderTracks) {
1661 // TPC + ITS primary vertex
1662 ftVertexer->SetITSMode();
1663 ftVertexer->SetConstraintOff();
1664 // get cuts for vertexer from AliGRPRecoParam
1666 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1667 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1668 grpRecoParam->GetVertexerTracksCutsITS(cutsVertexer);
1669 ftVertexer->SetCuts(cutsVertexer);
1670 delete [] cutsVertexer; cutsVertexer = NULL;
1671 if(fDiamondProfile && grpRecoParam->GetVertexerTracksConstraintITS())
1672 ftVertexer->SetVtxStart(fDiamondProfile);
1674 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1676 if (pvtx->GetStatus()) {
1677 fesd->SetPrimaryVertex(pvtx);
1678 for (Int_t i=0; i<ntracks; i++) {
1679 AliESDtrack *t = fesd->GetTrack(i);
1680 t->RelateToVertex(pvtx, kBz, kVeryBig);
1685 // TPC-only primary vertex
1686 ftVertexer->SetTPCMode();
1687 ftVertexer->SetConstraintOff();
1688 // get cuts for vertexer from AliGRPRecoParam
1690 Int_t nCutsVertexer = grpRecoParam->GetVertexerTracksNCuts();
1691 Double_t *cutsVertexer = new Double_t[nCutsVertexer];
1692 grpRecoParam->GetVertexerTracksCutsTPC(cutsVertexer);
1693 ftVertexer->SetCuts(cutsVertexer);
1694 delete [] cutsVertexer; cutsVertexer = NULL;
1695 if(fDiamondProfileTPC && grpRecoParam->GetVertexerTracksConstraintTPC())
1696 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1698 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1700 if (pvtx->GetStatus()) {
1701 fesd->SetPrimaryVertexTPC(pvtx);
1702 for (Int_t i=0; i<ntracks; i++) {
1703 AliESDtrack *t = fesd->GetTrack(i);
1704 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1710 delete[] selectedIdx;
1712 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1717 AliV0vertexer vtxer;
1718 vtxer.Tracks2V0vertices(fesd);
1720 if (fRunCascadeFinder) {
1722 AliCascadeVertexer cvtxer;
1723 cvtxer.V0sTracks2CascadeVertices(fesd);
1728 if (fCleanESD) CleanESD(fesd);
1731 fQASteer->RunOneEvent(fesd) ;
1734 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1735 if (qadm && fQATasks.Contains(Form("%d", AliQA::kESDS)))
1736 qadm->Exec(AliQA::kESDS, fesd);
1739 if (fWriteESDfriend) {
1740 // fesdf->~AliESDfriend();
1741 // new (fesdf) AliESDfriend(); // Reset...
1742 fesd->GetESDfriend(fesdf);
1746 // Auto-save the ESD tree in case of prompt reco @P2
1747 if (fRawReader && fRawReader->UseAutoSaveESD()) {
1748 ftree->AutoSave("SaveSelf");
1749 TFile *friendfile = (TFile *)(gROOT->GetListOfFiles()->FindObject("AliESDfriends.root"));
1750 if (friendfile) friendfile->Save();
1757 if (fRunAliEVE) RunAliEVE();
1761 if (fWriteESDfriend) {
1762 fesdf->~AliESDfriend();
1763 new (fesdf) AliESDfriend(); // Reset...
1766 ProcInfo_t ProcInfo;
1767 gSystem->GetProcInfo(&ProcInfo);
1768 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1771 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1772 if (fReconstructor[iDet])
1773 fReconstructor[iDet]->SetRecoParam(NULL);
1776 if (fRunQA || fRunGlobalQA)
1777 fQASteer->Increment() ;
1782 //_____________________________________________________________________________
1783 void AliReconstruction::SlaveTerminate()
1785 // Finalize the run on the slave side
1786 // Called after the exit
1787 // from the event loop
1788 AliCodeTimerAuto("");
1790 if (fIsNewRunLoader) { // galice.root didn't exist
1791 fRunLoader->WriteHeader("OVERWRITE");
1792 fRunLoader->CdGAFile();
1793 fRunLoader->Write(0, TObject::kOverwrite);
1796 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1797 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1799 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1800 cdbMapCopy->SetOwner(1);
1801 cdbMapCopy->SetName("cdbMap");
1802 TIter iter(cdbMap->GetTable());
1805 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1806 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1807 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1808 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1811 TList *cdbListCopy = new TList();
1812 cdbListCopy->SetOwner(1);
1813 cdbListCopy->SetName("cdbList");
1815 TIter iter2(cdbList);
1818 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1819 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1822 ftree->GetUserInfo()->Add(cdbMapCopy);
1823 ftree->GetUserInfo()->Add(cdbListCopy);
1826 if(fESDPar.Contains("ESD.par")){
1827 AliInfo("Attaching ESD.par to Tree");
1828 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1829 ftree->GetUserInfo()->Add(fn);
1835 if (fWriteESDfriend)
1836 ftree->SetBranchStatus("ESDfriend*",0);
1837 // we want to have only one tree version number
1838 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1841 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1842 if (fRunPlaneEff && !FinishPlaneEff()) {
1843 AliWarning("Finish PlaneEff evaluation failed");
1846 // End of cycle for the in-loop
1848 fQASteer->EndOfCycle() ;
1851 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1853 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1854 qadm->EndOfCycle(AliQA::kRECPOINTS);
1855 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1856 qadm->EndOfCycle(AliQA::kESDS);
1864 //_____________________________________________________________________________
1865 void AliReconstruction::Terminate()
1867 // Create tags for the events in the ESD tree (the ESD tree is always present)
1868 // In case of empty events the tags will contain dummy values
1869 AliCodeTimerAuto("");
1871 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1872 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1874 // Cleanup of CDB manager: cache and active storages!
1875 AliCDBManager::Instance()->ClearCache();
1878 //_____________________________________________________________________________
1879 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1881 // run the local reconstruction
1883 static Int_t eventNr=0;
1884 AliCodeTimerAuto("")
1886 TString detStr = detectors;
1887 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1888 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1889 AliReconstructor* reconstructor = GetReconstructor(iDet);
1890 if (!reconstructor) continue;
1891 AliLoader* loader = fLoader[iDet];
1892 // Matthias April 2008: temporary fix to run HLT reconstruction
1893 // although the HLT loader is missing
1894 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1896 reconstructor->Reconstruct(fRawReader, NULL);
1899 reconstructor->Reconstruct(dummy, NULL);
1904 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1907 // conversion of digits
1908 if (fRawReader && reconstructor->HasDigitConversion()) {
1909 AliInfo(Form("converting raw data digits into root objects for %s",
1910 fgkDetectorName[iDet]));
1911 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1912 // fgkDetectorName[iDet]));
1913 loader->LoadDigits("update");
1914 loader->CleanDigits();
1915 loader->MakeDigitsContainer();
1916 TTree* digitsTree = loader->TreeD();
1917 reconstructor->ConvertDigits(fRawReader, digitsTree);
1918 loader->WriteDigits("OVERWRITE");
1919 loader->UnloadDigits();
1921 // local reconstruction
1922 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1923 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1924 loader->LoadRecPoints("update");
1925 loader->CleanRecPoints();
1926 loader->MakeRecPointsContainer();
1927 TTree* clustersTree = loader->TreeR();
1928 if (fRawReader && !reconstructor->HasDigitConversion()) {
1929 reconstructor->Reconstruct(fRawReader, clustersTree);
1931 loader->LoadDigits("read");
1932 TTree* digitsTree = loader->TreeD();
1934 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1935 if (fStopOnError) return kFALSE;
1937 reconstructor->Reconstruct(digitsTree, clustersTree);
1939 loader->UnloadDigits();
1942 TString detQAStr(fQADetectors) ;
1944 fQASteer->RunOneEventInOneDetector(iDet, clustersTree) ;
1946 loader->WriteRecPoints("OVERWRITE");
1947 loader->UnloadRecPoints();
1948 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1950 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1951 AliError(Form("the following detectors were not found: %s",
1953 if (fStopOnError) return kFALSE;
1959 //_____________________________________________________________________________
1960 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1962 // run the barrel tracking
1964 AliCodeTimerAuto("")
1966 AliESDVertex* vertex = NULL;
1967 Double_t vtxPos[3] = {0, 0, 0};
1968 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1969 TArrayF mcVertex(3);
1970 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1971 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1972 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1976 AliInfo("running the ITS vertex finder");
1978 fLoader[0]->LoadRecPoints();
1979 TTree* cltree = fLoader[0]->TreeR();
1981 if(fDiamondProfileSPD) fVertexer->SetVtxStart(fDiamondProfileSPD);
1982 vertex = fVertexer->FindVertexForCurrentEvent(cltree);
1985 AliError("Can't get the ITS cluster tree");
1987 fLoader[0]->UnloadRecPoints();
1990 AliError("Can't get the ITS loader");
1993 AliWarning("Vertex not found");
1994 vertex = new AliESDVertex();
1995 vertex->SetName("default");
1998 vertex->SetName("reconstructed");
2002 AliInfo("getting the primary vertex from MC");
2003 vertex = new AliESDVertex(vtxPos, vtxErr);
2007 vertex->GetXYZ(vtxPos);
2008 vertex->GetSigmaXYZ(vtxErr);
2010 AliWarning("no vertex reconstructed");
2011 vertex = new AliESDVertex(vtxPos, vtxErr);
2013 esd->SetPrimaryVertexSPD(vertex);
2014 // if SPD multiplicity has been determined, it is stored in the ESD
2015 AliMultiplicity *mult = fVertexer->GetMultiplicity();
2016 if(mult)esd->SetMultiplicity(mult);
2018 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2019 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
2026 //_____________________________________________________________________________
2027 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
2029 // run the HLT barrel tracking
2031 AliCodeTimerAuto("")
2034 AliError("Missing runLoader!");
2038 AliInfo("running HLT tracking");
2040 // Get a pointer to the HLT reconstructor
2041 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
2042 if (!reconstructor) return kFALSE;
2045 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2046 TString detName = fgkDetectorName[iDet];
2047 AliDebug(1, Form("%s HLT tracking", detName.Data()));
2048 reconstructor->SetOption(detName.Data());
2049 AliTracker *tracker = reconstructor->CreateTracker();
2051 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
2052 if (fStopOnError) return kFALSE;
2056 Double_t vtxErr[3]={0.005,0.005,0.010};
2057 const AliESDVertex *vertex = esd->GetVertex();
2058 vertex->GetXYZ(vtxPos);
2059 tracker->SetVertex(vtxPos,vtxErr);
2061 fLoader[iDet]->LoadRecPoints("read");
2062 TTree* tree = fLoader[iDet]->TreeR();
2064 AliError(Form("Can't get the %s cluster tree", detName.Data()));
2067 tracker->LoadClusters(tree);
2069 if (tracker->Clusters2Tracks(esd) != 0) {
2070 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
2074 tracker->UnloadClusters();
2082 //_____________________________________________________________________________
2083 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
2085 // run the muon spectrometer tracking
2087 AliCodeTimerAuto("")
2090 AliError("Missing runLoader!");
2093 Int_t iDet = 7; // for MUON
2095 AliInfo("is running...");
2097 // Get a pointer to the MUON reconstructor
2098 AliReconstructor *reconstructor = GetReconstructor(iDet);
2099 if (!reconstructor) return kFALSE;
2102 TString detName = fgkDetectorName[iDet];
2103 AliDebug(1, Form("%s tracking", detName.Data()));
2104 AliTracker *tracker = reconstructor->CreateTracker();
2106 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2111 fLoader[iDet]->LoadRecPoints("read");
2113 tracker->LoadClusters(fLoader[iDet]->TreeR());
2115 Int_t rv = tracker->Clusters2Tracks(esd);
2119 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2123 fLoader[iDet]->UnloadRecPoints();
2125 tracker->UnloadClusters();
2133 //_____________________________________________________________________________
2134 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
2136 // run the barrel tracking
2137 static Int_t eventNr=0;
2138 AliCodeTimerAuto("")
2140 AliInfo("running tracking");
2142 //Fill the ESD with the T0 info (will be used by the TOF)
2143 if (fReconstructor[11] && fLoader[11]) {
2144 fLoader[11]->LoadRecPoints("READ");
2145 TTree *treeR = fLoader[11]->TreeR();
2146 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2149 // pass 1: TPC + ITS inwards
2150 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2151 if (!fTracker[iDet]) continue;
2152 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2155 fLoader[iDet]->LoadRecPoints("read");
2156 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2157 TTree* tree = fLoader[iDet]->TreeR();
2159 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2162 fTracker[iDet]->LoadClusters(tree);
2163 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2165 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2166 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2169 // preliminary PID in TPC needed by the ITS tracker
2171 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2172 AliESDpid::MakePID(esd);
2174 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2177 // pass 2: ALL backwards
2179 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2180 if (!fTracker[iDet]) continue;
2181 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2184 if (iDet > 1) { // all except ITS, TPC
2186 fLoader[iDet]->LoadRecPoints("read");
2187 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2188 tree = fLoader[iDet]->TreeR();
2190 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2193 fTracker[iDet]->LoadClusters(tree);
2194 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2198 if (iDet>1) // start filling residuals for the "outer" detectors
2199 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2201 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2202 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2207 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2208 fTracker[iDet]->UnloadClusters();
2209 fLoader[iDet]->UnloadRecPoints();
2211 // updated PID in TPC needed by the ITS tracker -MI
2213 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2214 AliESDpid::MakePID(esd);
2216 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2218 //stop filling residuals for the "outer" detectors
2219 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2221 // pass 3: TRD + TPC + ITS refit inwards
2223 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2224 if (!fTracker[iDet]) continue;
2225 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2228 if (iDet<2) // start filling residuals for TPC and ITS
2229 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2231 if (fTracker[iDet]->RefitInward(esd) != 0) {
2232 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2235 // run postprocessing
2236 if (fTracker[iDet]->PostProcess(esd) != 0) {
2237 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2240 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2243 // write space-points to the ESD in case alignment data output
2245 if (fWriteAlignmentData)
2246 WriteAlignmentData(esd);
2248 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2249 if (!fTracker[iDet]) continue;
2251 fTracker[iDet]->UnloadClusters();
2252 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2253 fLoader[iDet]->UnloadRecPoints();
2254 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2256 // stop filling residuals for TPC and ITS
2257 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2263 //_____________________________________________________________________________
2264 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2266 // Remove the data which are not needed for the physics analysis.
2269 Int_t nTracks=esd->GetNumberOfTracks();
2270 Int_t nV0s=esd->GetNumberOfV0s();
2272 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2274 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2275 Bool_t rc=esd->Clean(cleanPars);
2277 nTracks=esd->GetNumberOfTracks();
2278 nV0s=esd->GetNumberOfV0s();
2280 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2285 //_____________________________________________________________________________
2286 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2288 // fill the event summary data
2290 AliCodeTimerAuto("")
2291 static Int_t eventNr=0;
2292 TString detStr = detectors;
2294 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2295 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2296 AliReconstructor* reconstructor = GetReconstructor(iDet);
2297 if (!reconstructor) continue;
2298 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2299 TTree* clustersTree = NULL;
2300 if (fLoader[iDet]) {
2301 fLoader[iDet]->LoadRecPoints("read");
2302 clustersTree = fLoader[iDet]->TreeR();
2303 if (!clustersTree) {
2304 AliError(Form("Can't get the %s clusters tree",
2305 fgkDetectorName[iDet]));
2306 if (fStopOnError) return kFALSE;
2309 if (fRawReader && !reconstructor->HasDigitConversion()) {
2310 reconstructor->FillESD(fRawReader, clustersTree, esd);
2312 TTree* digitsTree = NULL;
2313 if (fLoader[iDet]) {
2314 fLoader[iDet]->LoadDigits("read");
2315 digitsTree = fLoader[iDet]->TreeD();
2317 AliError(Form("Can't get the %s digits tree",
2318 fgkDetectorName[iDet]));
2319 if (fStopOnError) return kFALSE;
2322 reconstructor->FillESD(digitsTree, clustersTree, esd);
2323 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2325 if (fLoader[iDet]) {
2326 fLoader[iDet]->UnloadRecPoints();
2330 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2331 AliError(Form("the following detectors were not found: %s",
2333 if (fStopOnError) return kFALSE;
2335 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2340 //_____________________________________________________________________________
2341 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2343 // Reads the trigger decision which is
2344 // stored in Trigger.root file and fills
2345 // the corresponding esd entries
2347 AliCodeTimerAuto("")
2349 AliInfo("Filling trigger information into the ESD");
2352 AliCTPRawStream input(fRawReader);
2353 if (!input.Next()) {
2354 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2357 if (esd->GetTriggerMask() != input.GetClassMask())
2358 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2359 input.GetClassMask(),esd->GetTriggerMask()));
2360 if (esd->GetOrbitNumber() != input.GetOrbitID())
2361 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2362 input.GetOrbitID(),esd->GetOrbitNumber()));
2363 if (esd->GetBunchCrossNumber() != input.GetBCID())
2364 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2365 input.GetBCID(),esd->GetBunchCrossNumber()));
2368 // Here one has to add the filling of trigger inputs and
2369 // interaction records
2379 //_____________________________________________________________________________
2380 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2383 // Filling information from RawReader Header
2386 if (!fRawReader) return kFALSE;
2388 AliInfo("Filling information from RawReader Header");
2390 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2391 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2392 esd->SetPeriodNumber(fRawReader->GetPeriod());
2394 esd->SetTimeStamp(fRawReader->GetTimestamp());
2395 esd->SetEventType(fRawReader->GetType());
2401 //_____________________________________________________________________________
2402 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2404 // check whether detName is contained in detectors
2405 // if yes, it is removed from detectors
2407 // check if all detectors are selected
2408 if ((detectors.CompareTo("ALL") == 0) ||
2409 detectors.BeginsWith("ALL ") ||
2410 detectors.EndsWith(" ALL") ||
2411 detectors.Contains(" ALL ")) {
2416 // search for the given detector
2417 Bool_t result = kFALSE;
2418 if ((detectors.CompareTo(detName) == 0) ||
2419 detectors.BeginsWith(detName+" ") ||
2420 detectors.EndsWith(" "+detName) ||
2421 detectors.Contains(" "+detName+" ")) {
2422 detectors.ReplaceAll(detName, "");
2426 // clean up the detectors string
2427 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2428 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2429 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2434 //_____________________________________________________________________________
2435 Bool_t AliReconstruction::InitRunLoader()
2437 // get or create the run loader
2439 if (gAlice) delete gAlice;
2442 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2443 // load all base libraries to get the loader classes
2444 TString libs = gSystem->GetLibraries();
2445 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2446 TString detName = fgkDetectorName[iDet];
2447 if (detName == "HLT") continue;
2448 if (libs.Contains("lib" + detName + "base.so")) continue;
2449 gSystem->Load("lib" + detName + "base.so");
2451 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2453 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2458 fRunLoader->CdGAFile();
2459 fRunLoader->LoadgAlice();
2461 //PH This is a temporary fix to give access to the kinematics
2462 //PH that is needed for the labels of ITS clusters
2463 fRunLoader->LoadHeader();
2464 fRunLoader->LoadKinematics();
2466 } else { // galice.root does not exist
2468 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2470 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2471 AliConfig::GetDefaultEventFolderName(),
2474 AliError(Form("could not create run loader in file %s",
2475 fGAliceFileName.Data()));
2479 fIsNewRunLoader = kTRUE;
2480 fRunLoader->MakeTree("E");
2482 if (fNumberOfEventsPerFile > 0)
2483 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2485 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2491 //_____________________________________________________________________________
2492 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2494 // get the reconstructor object and the loader for a detector
2496 if (fReconstructor[iDet]) {
2497 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2498 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2499 fReconstructor[iDet]->SetRecoParam(par);
2501 return fReconstructor[iDet];
2504 // load the reconstructor object
2505 TPluginManager* pluginManager = gROOT->GetPluginManager();
2506 TString detName = fgkDetectorName[iDet];
2507 TString recName = "Ali" + detName + "Reconstructor";
2509 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2511 AliReconstructor* reconstructor = NULL;
2512 // first check if a plugin is defined for the reconstructor
2513 TPluginHandler* pluginHandler =
2514 pluginManager->FindHandler("AliReconstructor", detName);
2515 // if not, add a plugin for it
2516 if (!pluginHandler) {
2517 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2518 TString libs = gSystem->GetLibraries();
2519 if (libs.Contains("lib" + detName + "base.so") ||
2520 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2521 pluginManager->AddHandler("AliReconstructor", detName,
2522 recName, detName + "rec", recName + "()");
2524 pluginManager->AddHandler("AliReconstructor", detName,
2525 recName, detName, recName + "()");
2527 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2529 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2530 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2532 if (reconstructor) {
2533 TObject* obj = fOptions.FindObject(detName.Data());
2534 if (obj) reconstructor->SetOption(obj->GetTitle());
2535 reconstructor->Init();
2536 fReconstructor[iDet] = reconstructor;
2539 // get or create the loader
2540 if (detName != "HLT") {
2541 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2542 if (!fLoader[iDet]) {
2543 AliConfig::Instance()
2544 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2546 // first check if a plugin is defined for the loader
2548 pluginManager->FindHandler("AliLoader", detName);
2549 // if not, add a plugin for it
2550 if (!pluginHandler) {
2551 TString loaderName = "Ali" + detName + "Loader";
2552 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2553 pluginManager->AddHandler("AliLoader", detName,
2554 loaderName, detName + "base",
2555 loaderName + "(const char*, TFolder*)");
2556 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2558 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2560 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2561 fRunLoader->GetEventFolder());
2563 if (!fLoader[iDet]) { // use default loader
2564 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2566 if (!fLoader[iDet]) {
2567 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2568 if (fStopOnError) return NULL;
2570 fRunLoader->AddLoader(fLoader[iDet]);
2571 fRunLoader->CdGAFile();
2572 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2573 fRunLoader->Write(0, TObject::kOverwrite);
2578 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2579 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2580 reconstructor->SetRecoParam(par);
2582 return reconstructor;
2585 //_____________________________________________________________________________
2586 Bool_t AliReconstruction::CreateVertexer()
2588 // create the vertexer
2591 AliReconstructor* itsReconstructor = GetReconstructor(0);
2592 if (itsReconstructor) {
2593 fVertexer = itsReconstructor->CreateVertexer();
2596 AliWarning("couldn't create a vertexer for ITS");
2597 if (fStopOnError) return kFALSE;
2603 //_____________________________________________________________________________
2604 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2606 // create the trackers
2607 AliInfo("Creating trackers");
2609 TString detStr = detectors;
2610 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2611 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2612 AliReconstructor* reconstructor = GetReconstructor(iDet);
2613 if (!reconstructor) continue;
2614 TString detName = fgkDetectorName[iDet];
2615 if (detName == "HLT") {
2616 fRunHLTTracking = kTRUE;
2619 if (detName == "MUON") {
2620 fRunMuonTracking = kTRUE;
2625 fTracker[iDet] = reconstructor->CreateTracker();
2626 if (!fTracker[iDet] && (iDet < 7)) {
2627 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2628 if (fStopOnError) return kFALSE;
2630 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2636 //_____________________________________________________________________________
2637 void AliReconstruction::CleanUp()
2639 // delete trackers and the run loader and close and delete the file
2641 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2642 delete fReconstructor[iDet];
2643 fReconstructor[iDet] = NULL;
2644 fLoader[iDet] = NULL;
2645 delete fTracker[iDet];
2646 fTracker[iDet] = NULL;
2657 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2658 delete fDiamondProfileSPD;
2659 fDiamondProfileSPD = NULL;
2660 delete fDiamondProfile;
2661 fDiamondProfile = NULL;
2662 delete fDiamondProfileTPC;
2663 fDiamondProfileTPC = NULL;
2669 delete fParentRawReader;
2670 fParentRawReader=NULL;
2679 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2681 // Write space-points which are then used in the alignment procedures
2682 // For the moment only ITS, TPC, TRD and TOF
2684 Int_t ntracks = esd->GetNumberOfTracks();
2685 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2687 AliESDtrack *track = esd->GetTrack(itrack);
2690 for (Int_t iDet = 5; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2691 nsp += track->GetNcls(iDet);
2693 if (iDet==0) { // ITS "extra" clusters
2694 track->GetClusters(iDet,idx);
2695 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2700 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2701 track->SetTrackPointArray(sp);
2703 for (Int_t iDet = 5; iDet >= 0; iDet--) {
2704 AliTracker *tracker = fTracker[iDet];
2705 if (!tracker) continue;
2706 Int_t nspdet = track->GetClusters(iDet,idx);
2708 if (iDet==0) // ITS "extra" clusters
2709 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2711 if (nspdet <= 0) continue;
2715 while (isp2 < nspdet) {
2716 Bool_t isvalid=kTRUE;
2718 Int_t index=idx[isp++];
2719 if (index < 0) continue;
2721 TString dets = fgkDetectorName[iDet];
2722 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2723 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2724 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2725 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2726 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2728 isvalid = tracker->GetTrackPoint(index,p);
2731 if (!isvalid) continue;
2732 sp->AddPoint(isptrack,&p); isptrack++;
2739 //_____________________________________________________________________________
2740 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2742 // The method reads the raw-data error log
2743 // accumulated within the rawReader.
2744 // It extracts the raw-data errors related to
2745 // the current event and stores them into
2746 // a TClonesArray inside the esd object.
2748 if (!fRawReader) return;
2750 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2752 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2754 if (iEvent != log->GetEventNumber()) continue;
2756 esd->AddRawDataErrorLog(log);
2761 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2762 // Dump a file content into a char in TNamed
2764 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2765 Int_t kBytes = (Int_t)in.tellg();
2766 printf("Size: %d \n",kBytes);
2769 char* memblock = new char [kBytes];
2770 in.seekg (0, ios::beg);
2771 in.read (memblock, kBytes);
2773 TString fData(memblock,kBytes);
2774 fn = new TNamed(pName,fData);
2775 printf("fData Size: %d \n",fData.Sizeof());
2776 printf("pName Size: %d \n",pName.Sizeof());
2777 printf("fn Size: %d \n",fn->Sizeof());
2781 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2787 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2788 // This is not really needed in AliReconstruction at the moment
2789 // but can serve as a template
2791 TList *fList = fTree->GetUserInfo();
2792 TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2793 printf("fn Size: %d \n",fn->Sizeof());
2795 TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2796 const char* cdata = fn->GetTitle();
2797 printf("fTmp Size %d\n",fTmp.Sizeof());
2799 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2800 printf("calculated size %d\n",size);
2801 ofstream out(pName.Data(),ios::out | ios::binary);
2802 out.write(cdata,size);
2807 //_____________________________________________________________________________
2808 void AliReconstruction::CheckQA()
2810 // check the QA of SIM for this run and remove the detectors
2811 // with status Fatal
2813 TString newRunLocalReconstruction ;
2814 TString newRunTracking ;
2815 TString newFillESD ;
2817 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2818 TString detName(AliQA::GetDetName(iDet)) ;
2819 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2820 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2821 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2823 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2824 fRunLocalReconstruction.Contains("ALL") ) {
2825 newRunLocalReconstruction += detName ;
2826 newRunLocalReconstruction += " " ;
2828 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2829 fRunTracking.Contains("ALL") ) {
2830 newRunTracking += detName ;
2831 newRunTracking += " " ;
2833 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2834 fFillESD.Contains("ALL") ) {
2835 newFillESD += detName ;
2840 fRunLocalReconstruction = newRunLocalReconstruction ;
2841 fRunTracking = newRunTracking ;
2842 fFillESD = newFillESD ;
2845 //_____________________________________________________________________________
2846 Int_t AliReconstruction::GetDetIndex(const char* detector)
2848 // return the detector index corresponding to detector
2850 for (index = 0; index < fgkNDetectors ; index++) {
2851 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2856 //_____________________________________________________________________________
2857 Bool_t AliReconstruction::FinishPlaneEff() {
2859 // Here execute all the necessary operationis, at the end of the tracking phase,
2860 // in case that evaluation of PlaneEfficiencies was required for some detector.
2861 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2863 // This Preliminary version works only FOR ITS !!!!!
2864 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2867 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2870 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2871 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2872 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2873 if(fTracker[iDet]) {
2874 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2875 TString name=planeeff->GetName();
2877 TFile* pefile = TFile::Open(name, "RECREATE");
2878 ret=(Bool_t)planeeff->Write();
2880 if(planeeff->GetCreateHistos()) {
2881 TString hname=planeeff->GetName();
2882 hname+="Histo.root";
2883 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
2889 //_____________________________________________________________________________
2890 Bool_t AliReconstruction::InitPlaneEff() {
2892 // Here execute all the necessary operations, before of the tracking phase,
2893 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2894 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2895 // which should be updated/recalculated.
2897 // This Preliminary version will work only FOR ITS !!!!!
2898 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2901 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2903 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2907 //_____________________________________________________________________________
2908 Bool_t AliReconstruction::InitAliEVE()
2910 // This method should be called only in case
2911 // AliReconstruction is run
2912 // within the alieve environment.
2913 // It will initialize AliEVE in a way
2914 // so that it can visualize event processed
2915 // by AliReconstruction.
2916 // The return flag shows whenever the
2917 // AliEVE initialization was successful or not.
2920 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2921 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2922 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2924 gROOT->ProcessLine("if (!AliEveEventManager::GetMaster()){new AliEveEventManager();AliEveEventManager::GetMaster()->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(AliEveEventManager::GetMaster());};");
2925 gROOT->ProcessLine("alieve_online_init()");
2930 //_____________________________________________________________________________
2931 void AliReconstruction::RunAliEVE()
2933 // Runs AliEVE visualisation of
2934 // the current event.
2935 // Should be executed only after
2936 // successful initialization of AliEVE.
2938 AliInfo("Running AliEVE...");
2939 gROOT->ProcessLine(Form("AliEveEventManager::GetMaster()->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2943 //_____________________________________________________________________________
2944 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2946 // Allows to run QA for a selected set of detectors
2947 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2948 // all selected detectors run the same selected tasks
2950 if (!detAndAction.Contains(":")) {
2951 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2955 Int_t colon = detAndAction.Index(":") ;
2956 fQADetectors = detAndAction(0, colon) ;
2957 if (fQADetectors.Contains("ALL") )
2958 fQADetectors = fFillESD ;
2959 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2960 if (fQATasks.Contains("ALL") ) {
2961 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2963 fQATasks.ToUpper() ;
2965 if ( fQATasks.Contains("RAW") )
2966 tempo = Form("%d ", AliQA::kRAWS) ;
2967 if ( fQATasks.Contains("RECPOINT") )
2968 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2969 if ( fQATasks.Contains("ESD") )
2970 tempo += Form("%d ", AliQA::kESDS) ;
2972 if (fQATasks.IsNull()) {
2973 AliInfo("No QA requested\n") ;
2978 TString tempo(fQATasks) ;
2979 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2980 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2981 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2982 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2987 //_____________________________________________________________________________
2988 Bool_t AliReconstruction::InitRecoParams()
2990 // The method accesses OCDB and retrieves all
2991 // the available reco-param objects from there.
2993 Bool_t isOK = kTRUE;
2995 TString detStr = fLoadCDB;
2996 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2998 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
3000 if (fRecoParam.GetDetRecoParamArray(iDet)) {
3001 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
3005 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
3007 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
3008 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
3010 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
3014 TObject *recoParamObj = entry->GetObject();
3015 if (dynamic_cast<TObjArray*>(recoParamObj)) {
3016 // The detector has a normal TobjArray of AliDetectorRecoParam objects
3017 // Registering them in AliRecoParam
3018 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
3020 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
3021 // The detector has only onse set of reco parameters
3022 // Registering it in AliRecoParam
3023 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
3024 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
3025 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
3028 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
3032 AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
3036 if (AliDebugLevel() > 0) fRecoParam.Print();
3041 //_____________________________________________________________________________
3042 Bool_t AliReconstruction::GetEventInfo()
3044 // Fill the event info object
3046 AliCodeTimerAuto("")
3048 AliCentralTrigger *aCTP = NULL;
3050 fEventInfo.SetEventType(fRawReader->GetType());
3052 ULong64_t mask = fRawReader->GetClassMask();
3053 fEventInfo.SetTriggerMask(mask);
3054 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
3055 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
3057 aCTP = new AliCentralTrigger();
3058 TString configstr("");
3059 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
3060 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
3064 aCTP->SetClassMask(mask);
3065 aCTP->SetClusterMask(clmask);
3068 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
3070 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
3071 aCTP = fRunLoader->GetTrigger();
3072 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
3073 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
3076 AliWarning("No trigger can be loaded! The trigger information will not be used!");
3081 AliTriggerConfiguration *config = aCTP->GetConfiguration();
3083 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
3084 if (fRawReader) delete aCTP;
3088 UChar_t clustmask = 0;
3090 ULong64_t trmask = fEventInfo.GetTriggerMask();
3091 const TObjArray& classesArray = config->GetClasses();
3092 Int_t nclasses = classesArray.GetEntriesFast();
3093 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
3094 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
3096 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
3097 fesd->SetTriggerClass(trclass->GetName(),trindex);
3098 if (trmask & (1 << trindex)) {
3100 trclasses += trclass->GetName();
3102 clustmask |= trclass->GetCluster()->GetClusterMask();
3106 fEventInfo.SetTriggerClasses(trclasses);
3108 // Set the information in ESD
3109 fesd->SetTriggerMask(trmask);
3110 fesd->SetTriggerCluster(clustmask);
3112 if (!aCTP->CheckTriggeredDetectors()) {
3113 if (fRawReader) delete aCTP;
3117 if (fRawReader) delete aCTP;
3119 // We have to fill also the HLT decision here!!
3125 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
3127 // Match the detector list found in the rec.C or the default 'ALL'
3128 // to the list found in the GRP (stored there by the shuttle PP which
3129 // gets the information from ECS)
3130 static TString resultList;
3131 TString detList = detectorList;
3135 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
3136 if ((detectorMask >> iDet) & 0x1) {
3137 TString det = AliDAQ::OfflineModuleName(iDet);
3138 if ((detList.CompareTo("ALL") == 0) ||
3139 detList.BeginsWith("ALL ") ||
3140 detList.EndsWith(" ALL") ||
3141 detList.Contains(" ALL ") ||
3142 (detList.CompareTo(det) == 0) ||
3143 detList.BeginsWith(det) ||
3144 detList.EndsWith(det) ||
3145 detList.Contains( " "+det+" " )) {
3146 if (!resultList.EndsWith(det + " ")) {
3155 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3156 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3157 if ((detList.CompareTo("ALL") == 0) ||
3158 detList.BeginsWith("ALL ") ||
3159 detList.EndsWith(" ALL") ||
3160 detList.Contains(" ALL ") ||
3161 (detList.CompareTo(hltDet) == 0) ||
3162 detList.BeginsWith(hltDet) ||
3163 detList.EndsWith(hltDet) ||
3164 detList.Contains( " "+hltDet+" " )) {
3165 resultList += hltDet;
3169 return resultList.Data();
3173 //______________________________________________________________________________
3174 void AliReconstruction::Abort(const char *method, EAbort what)
3176 // Abort processing. If what = kAbortProcess, the Process() loop will be
3177 // aborted. If what = kAbortFile, the current file in a chain will be
3178 // aborted and the processing will continue with the next file, if there
3179 // is no next file then Process() will be aborted. Abort() can also be
3180 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3181 // the SlaveTerminate() and Terminate() are always called. The abort flag
3182 // can be checked in these methods using GetAbort().
3184 // The method is overwritten in AliReconstruction for better handling of
3185 // reco specific errors
3187 if (!fStopOnError) return;
3191 TString whyMess = method;
3192 whyMess += " failed! Aborting...";
3194 AliError(whyMess.Data());
3197 TString mess = "Abort";
3198 if (fAbort == kAbortProcess)
3199 mess = "AbortProcess";
3200 else if (fAbort == kAbortFile)
3203 Info(mess, whyMess.Data());