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 "AliRunInfo.h"
190 #include "AliEventInfo.h"
194 ClassImp(AliReconstruction)
196 //_____________________________________________________________________________
197 const char* AliReconstruction::fgkDetectorName[AliReconstruction::fgkNDetectors] = {"ITS", "TPC", "TRD", "TOF", "PHOS", "HMPID", "EMCAL", "MUON", "FMD", "ZDC", "PMD", "T0", "VZERO", "ACORDE", "HLT"};
199 //_____________________________________________________________________________
200 AliReconstruction::AliReconstruction(const char* gAliceFilename) :
202 fUniformField(kFALSE),
203 fForcedFieldMap(NULL),
204 fRunVertexFinder(kTRUE),
205 fRunVertexFinderTracks(kTRUE),
206 fRunHLTTracking(kFALSE),
207 fRunMuonTracking(kFALSE),
209 fRunCascadeFinder(kTRUE),
210 fStopOnError(kFALSE),
211 fWriteAlignmentData(kFALSE),
212 fWriteESDfriend(kFALSE),
213 fFillTriggerESD(kTRUE),
221 fRunLocalReconstruction("ALL"),
225 fUseTrackingErrorsForAlignment(""),
226 fGAliceFileName(gAliceFilename),
231 fNumberOfEventsPerFile(1),
233 fLoadAlignFromCDB(kTRUE),
234 fLoadAlignData("ALL"),
242 fParentRawReader(NULL),
247 fDiamondProfile(NULL),
248 fDiamondProfileTPC(NULL),
249 fMeanVertexConstraint(kTRUE),
253 fAlignObjArray(NULL),
256 fInitCDBCalled(kFALSE),
257 fSetRunNumberFromDataCalled(kFALSE),
263 fSameQACycle(kFALSE),
265 fRunPlaneEff(kFALSE),
274 fIsNewRunLoader(kFALSE),
278 // create reconstruction object with default parameters
281 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
282 fReconstructor[iDet] = NULL;
283 fLoader[iDet] = NULL;
284 fTracker[iDet] = NULL;
285 fQACycles[iDet] = 999999;
290 //_____________________________________________________________________________
291 AliReconstruction::AliReconstruction(const AliReconstruction& rec) :
293 fUniformField(rec.fUniformField),
294 fForcedFieldMap(NULL),
295 fRunVertexFinder(rec.fRunVertexFinder),
296 fRunVertexFinderTracks(rec.fRunVertexFinderTracks),
297 fRunHLTTracking(rec.fRunHLTTracking),
298 fRunMuonTracking(rec.fRunMuonTracking),
299 fRunV0Finder(rec.fRunV0Finder),
300 fRunCascadeFinder(rec.fRunCascadeFinder),
301 fStopOnError(rec.fStopOnError),
302 fWriteAlignmentData(rec.fWriteAlignmentData),
303 fWriteESDfriend(rec.fWriteESDfriend),
304 fFillTriggerESD(rec.fFillTriggerESD),
306 fCleanESD(rec.fCleanESD),
307 fV0DCAmax(rec.fV0DCAmax),
308 fV0CsPmin(rec.fV0CsPmin),
312 fRunLocalReconstruction(rec.fRunLocalReconstruction),
313 fRunTracking(rec.fRunTracking),
314 fFillESD(rec.fFillESD),
315 fLoadCDB(rec.fLoadCDB),
316 fUseTrackingErrorsForAlignment(rec.fUseTrackingErrorsForAlignment),
317 fGAliceFileName(rec.fGAliceFileName),
318 fRawInput(rec.fRawInput),
319 fEquipIdMap(rec.fEquipIdMap),
320 fFirstEvent(rec.fFirstEvent),
321 fLastEvent(rec.fLastEvent),
322 fNumberOfEventsPerFile(rec.fNumberOfEventsPerFile),
324 fLoadAlignFromCDB(rec.fLoadAlignFromCDB),
325 fLoadAlignData(rec.fLoadAlignData),
326 fESDPar(rec.fESDPar),
327 fUseHLTData(rec.fUseHLTData),
333 fParentRawReader(NULL),
335 fRecoParam(rec.fRecoParam),
338 fDiamondProfile(rec.fDiamondProfile),
339 fDiamondProfileTPC(rec.fDiamondProfileTPC),
340 fMeanVertexConstraint(rec.fMeanVertexConstraint),
344 fAlignObjArray(rec.fAlignObjArray),
345 fCDBUri(rec.fCDBUri),
347 fInitCDBCalled(rec.fInitCDBCalled),
348 fSetRunNumberFromDataCalled(rec.fSetRunNumberFromDataCalled),
349 fQADetectors(rec.fQADetectors),
351 fQATasks(rec.fQATasks),
353 fRunGlobalQA(rec.fRunGlobalQA),
354 fSameQACycle(rec.fSameQACycle),
355 fRunPlaneEff(rec.fRunPlaneEff),
364 fIsNewRunLoader(rec.fIsNewRunLoader),
370 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
371 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
373 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
374 fReconstructor[iDet] = NULL;
375 fLoader[iDet] = NULL;
376 fTracker[iDet] = NULL;
377 fQACycles[iDet] = rec.fQACycles[iDet];
379 for (Int_t i = 0; i < rec.fSpecCDBUri.GetEntriesFast(); i++) {
380 if (rec.fSpecCDBUri[i]) fSpecCDBUri.Add(rec.fSpecCDBUri[i]->Clone());
384 //_____________________________________________________________________________
385 AliReconstruction& AliReconstruction::operator = (const AliReconstruction& rec)
387 // assignment operator
388 // Used in PROOF mode
389 // Be very careful while modifing it!
390 // Simple rules to follow:
391 // for persistent data members - use their assignment operators
392 // for non-persistent ones - do nothing or take the default values from constructor
393 // TSelector members should not be touched
394 if(&rec == this) return *this;
396 fUniformField = rec.fUniformField;
397 fForcedFieldMap = NULL;
398 fRunVertexFinder = rec.fRunVertexFinder;
399 fRunVertexFinderTracks = rec.fRunVertexFinderTracks;
400 fRunHLTTracking = rec.fRunHLTTracking;
401 fRunMuonTracking = rec.fRunMuonTracking;
402 fRunV0Finder = rec.fRunV0Finder;
403 fRunCascadeFinder = rec.fRunCascadeFinder;
404 fStopOnError = rec.fStopOnError;
405 fWriteAlignmentData = rec.fWriteAlignmentData;
406 fWriteESDfriend = rec.fWriteESDfriend;
407 fFillTriggerESD = rec.fFillTriggerESD;
409 fCleanESD = rec.fCleanESD;
410 fV0DCAmax = rec.fV0DCAmax;
411 fV0CsPmin = rec.fV0CsPmin;
415 fRunLocalReconstruction = rec.fRunLocalReconstruction;
416 fRunTracking = rec.fRunTracking;
417 fFillESD = rec.fFillESD;
418 fLoadCDB = rec.fLoadCDB;
419 fUseTrackingErrorsForAlignment = rec.fUseTrackingErrorsForAlignment;
420 fGAliceFileName = rec.fGAliceFileName;
421 fRawInput = rec.fRawInput;
422 fEquipIdMap = rec.fEquipIdMap;
423 fFirstEvent = rec.fFirstEvent;
424 fLastEvent = rec.fLastEvent;
425 fNumberOfEventsPerFile = rec.fNumberOfEventsPerFile;
427 for (Int_t i = 0; i < rec.fOptions.GetEntriesFast(); i++) {
428 if (rec.fOptions[i]) fOptions.Add(rec.fOptions[i]->Clone());
431 fLoadAlignFromCDB = rec.fLoadAlignFromCDB;
432 fLoadAlignData = rec.fLoadAlignData;
433 fESDPar = rec.fESDPar;
434 fUseHLTData = rec.fUseHLTData;
436 delete fRunInfo; fRunInfo = NULL;
437 if (rec.fRunInfo) fRunInfo = new AliRunInfo(*rec.fRunInfo);
439 fEventInfo = rec.fEventInfo;
443 fParentRawReader = NULL;
445 fRecoParam = rec.fRecoParam;
447 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
448 delete fReconstructor[iDet]; fReconstructor[iDet] = NULL;
449 delete fLoader[iDet]; fLoader[iDet] = NULL;
450 delete fTracker[iDet]; fTracker[iDet] = NULL;
451 fQACycles[iDet] = rec.fQACycles[iDet];
455 delete fDiamondProfile; fDiamondProfile = NULL;
456 if (rec.fDiamondProfile) fDiamondProfile = new AliESDVertex(*rec.fDiamondProfile);
457 delete fDiamondProfileTPC; fDiamondProfileTPC = NULL;
458 if (rec.fDiamondProfileTPC) fDiamondProfileTPC = new AliESDVertex(*rec.fDiamondProfileTPC);
459 fMeanVertexConstraint = rec.fMeanVertexConstraint;
461 delete fGRPData; fGRPData = NULL;
462 if (rec.fGRPData) fGRPData = (TMap*)((rec.fGRPData)->Clone());
464 delete fAlignObjArray; fAlignObjArray = NULL;
467 fSpecCDBUri.Delete();
468 fInitCDBCalled = rec.fInitCDBCalled;
469 fSetRunNumberFromDataCalled = rec.fSetRunNumberFromDataCalled;
470 fQADetectors = rec.fQADetectors;
472 fQATasks = rec.fQATasks;
474 fRunGlobalQA = rec.fRunGlobalQA;
475 fSameQACycle = rec.fSameQACycle;
476 fRunPlaneEff = rec.fRunPlaneEff;
485 fIsNewRunLoader = rec.fIsNewRunLoader;
492 //_____________________________________________________________________________
493 AliReconstruction::~AliReconstruction()
498 delete fForcedFieldMap;
500 if (fAlignObjArray) {
501 fAlignObjArray->Delete();
502 delete fAlignObjArray;
504 fSpecCDBUri.Delete();
506 AliCodeTimer::Instance()->Print();
509 //_____________________________________________________________________________
510 void AliReconstruction::InitCDB()
512 // activate a default CDB storage
513 // First check if we have any CDB storage set, because it is used
514 // to retrieve the calibration and alignment constants
515 AliCodeTimerAuto("");
517 if (fInitCDBCalled) return;
518 fInitCDBCalled = kTRUE;
520 AliCDBManager* man = AliCDBManager::Instance();
521 if (man->IsDefaultStorageSet())
523 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
524 AliWarning("Default CDB storage has been already set !");
525 AliWarning(Form("Ignoring the default storage declared in AliReconstruction: %s",fCDBUri.Data()));
526 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
527 fCDBUri = man->GetDefaultStorage()->GetURI();
530 if (fCDBUri.Length() > 0)
532 AliDebug(2,"!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
533 AliDebug(2, Form("Default CDB storage is set to: %s", fCDBUri.Data()));
534 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
536 fCDBUri="local://$ALICE_ROOT";
537 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
538 AliWarning("Default CDB storage not yet set !!!!");
539 AliWarning(Form("Setting it now to: %s", fCDBUri.Data()));
540 AliWarning("!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
543 man->SetDefaultStorage(fCDBUri);
546 // Now activate the detector specific CDB storage locations
547 for (Int_t i = 0; i < fSpecCDBUri.GetEntriesFast(); i++) {
548 TObject* obj = fSpecCDBUri[i];
550 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
551 AliDebug(2, Form("Specific CDB storage for %s is set to: %s",obj->GetName(),obj->GetTitle()));
552 AliDebug(2, "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!");
553 man->SetSpecificStorage(obj->GetName(), obj->GetTitle());
555 AliSysInfo::AddStamp("InitCDB");
558 //_____________________________________________________________________________
559 void AliReconstruction::SetDefaultStorage(const char* uri) {
560 // Store the desired default CDB storage location
561 // Activate it later within the Run() method
567 //_____________________________________________________________________________
568 void AliReconstruction::SetSpecificStorage(const char* calibType, const char* uri) {
569 // Store a detector-specific CDB storage location
570 // Activate it later within the Run() method
572 AliCDBPath aPath(calibType);
573 if(!aPath.IsValid()){
574 // if calibType is not wildcard but it is a valid detector, add "/*" to make it a valid path
575 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
576 if(!strcmp(calibType, fgkDetectorName[iDet])) {
577 aPath.SetPath(Form("%s/*", calibType));
578 AliInfo(Form("Path for specific storage set to %s", aPath.GetPath().Data()));
582 if(!aPath.IsValid()){
583 AliError(Form("Not a valid path or detector: %s", calibType));
588 // // check that calibType refers to a "valid" detector name
589 // Bool_t isDetector = kFALSE;
590 // for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
591 // TString detName = fgkDetectorName[iDet];
592 // if(aPath.GetLevel0() == detName) {
593 // isDetector = kTRUE;
599 // AliError(Form("Not a valid detector: %s", aPath.GetLevel0().Data()));
603 TObject* obj = fSpecCDBUri.FindObject(aPath.GetPath().Data());
604 if (obj) fSpecCDBUri.Remove(obj);
605 fSpecCDBUri.Add(new TNamed(aPath.GetPath().Data(), uri));
609 //_____________________________________________________________________________
610 Bool_t AliReconstruction::SetRunNumberFromData()
612 // The method is called in Run() in order
613 // to set a correct run number.
614 // In case of raw data reconstruction the
615 // run number is taken from the raw data header
617 if (fSetRunNumberFromDataCalled) return kTRUE;
618 fSetRunNumberFromDataCalled = kTRUE;
620 AliCDBManager* man = AliCDBManager::Instance();
623 if(fRawReader->NextEvent()) {
624 if(man->GetRun() > 0) {
625 AliWarning("Run number is taken from raw-event header! Ignoring settings in AliCDBManager!");
627 man->SetRun(fRawReader->GetRunNumber());
628 fRawReader->RewindEvents();
631 if(man->GetRun() > 0) {
632 AliWarning("No raw-data events are found ! Using settings in AliCDBManager !");
635 AliWarning("Neither raw events nor settings in AliCDBManager are found !");
641 AliRunLoader *rl = AliRunLoader::Open(fGAliceFileName.Data());
643 AliError(Form("No run loader found in file %s", fGAliceFileName.Data()));
648 // read run number from gAlice
649 if(rl->GetHeader()) {
650 man->SetRun(rl->GetHeader()->GetRun());
655 AliError("Neither run-loader header nor RawReader objects are found !");
667 //_____________________________________________________________________________
668 void AliReconstruction::SetCDBLock() {
669 // Set CDB lock: from now on it is forbidden to reset the run number
670 // or the default storage or to activate any further storage!
672 AliCDBManager::Instance()->SetLock(1);
675 //_____________________________________________________________________________
676 Bool_t AliReconstruction::MisalignGeometry(const TString& detectors)
678 // Read the alignment objects from CDB.
679 // Each detector is supposed to have the
680 // alignment objects in DET/Align/Data CDB path.
681 // All the detector objects are then collected,
682 // sorted by geometry level (starting from ALIC) and
683 // then applied to the TGeo geometry.
684 // Finally an overlaps check is performed.
686 // Load alignment data from CDB and fill fAlignObjArray
687 if(fLoadAlignFromCDB){
689 TString detStr = detectors;
690 TString loadAlObjsListOfDets = "";
692 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
693 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
694 loadAlObjsListOfDets += fgkDetectorName[iDet];
695 loadAlObjsListOfDets += " ";
696 } // end loop over detectors
697 loadAlObjsListOfDets.Prepend("GRP "); //add alignment objects for non-sensitive modules
698 AliGeomManager::ApplyAlignObjsFromCDB(loadAlObjsListOfDets.Data());
699 AliCDBManager::Instance()->UnloadFromCache("*/Align/*");
701 // Check if the array with alignment objects was
702 // provided by the user. If yes, apply the objects
703 // to the present TGeo geometry
704 if (fAlignObjArray) {
705 if (gGeoManager && gGeoManager->IsClosed()) {
706 if (AliGeomManager::ApplyAlignObjsToGeom(*fAlignObjArray) == kFALSE) {
707 AliError("The misalignment of one or more volumes failed!"
708 "Compare the list of simulated detectors and the list of detector alignment data!");
713 AliError("Can't apply the misalignment! gGeoManager doesn't exist or it is still opened!");
719 if (fAlignObjArray) {
720 fAlignObjArray->Delete();
721 delete fAlignObjArray; fAlignObjArray=NULL;
727 //_____________________________________________________________________________
728 void AliReconstruction::SetGAliceFile(const char* fileName)
730 // set the name of the galice file
732 fGAliceFileName = fileName;
735 //_____________________________________________________________________________
736 void AliReconstruction::SetInput(const char* input)
738 // In case the input string starts with 'mem://', we run in an online mode
739 // and AliRawReaderDateOnline object is created. In all other cases a raw-data
740 // file is assumed. One can give as an input:
741 // mem://: - events taken from DAQ monitoring libs online
743 // mem://<filename> - emulation of the above mode (via DATE monitoring libs)
744 if (input) fRawInput = input;
747 //_____________________________________________________________________________
748 void AliReconstruction::SetOption(const char* detector, const char* option)
750 // set options for the reconstruction of a detector
752 TObject* obj = fOptions.FindObject(detector);
753 if (obj) fOptions.Remove(obj);
754 fOptions.Add(new TNamed(detector, option));
757 //_____________________________________________________________________________
758 void AliReconstruction::SetRecoParam(const char* detector, AliDetectorRecoParam *par)
760 // Set custom reconstruction parameters for a given detector
761 // Single set of parameters for all the events
762 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
763 if(!strcmp(detector, fgkDetectorName[iDet])) {
765 fRecoParam.AddDetRecoParam(iDet,par);
772 //_____________________________________________________________________________
773 Bool_t AliReconstruction::SetFieldMap(Float_t l3Current, Float_t diCurrent, Float_t factor, const char *path) {
774 //------------------------------------------------
775 // The magnetic field map, defined externally...
776 // L3 current 30000 A -> 0.5 T
777 // L3 current 12000 A -> 0.2 T
778 // dipole current 6000 A
779 // The polarities must be the same
780 //------------------------------------------------
781 const Float_t l3NominalCurrent1=30000.; // (A)
782 const Float_t l3NominalCurrent2=12000.; // (A)
783 const Float_t diNominalCurrent =6000. ; // (A)
785 const Float_t tolerance=0.03; // relative current tolerance
786 const Float_t zero=77.; // "zero" current (A)
789 Bool_t dipoleON=kFALSE;
791 TString s=(factor < 0) ? "L3: -" : "L3: +";
793 l3Current = TMath::Abs(l3Current);
794 if (TMath::Abs(l3Current-l3NominalCurrent1)/l3NominalCurrent1 < tolerance) {
795 map=AliMagWrapCheb::k5kG;
798 if (TMath::Abs(l3Current-l3NominalCurrent2)/l3NominalCurrent2 < tolerance) {
799 map=AliMagWrapCheb::k2kG;
802 if (l3Current < zero) {
803 map=AliMagWrapCheb::k2kG;
805 factor=0.; // in fact, this is a global factor...
806 fUniformField=kTRUE; // track with the uniform (zero) B field
808 AliError(Form("Wrong L3 current (%f A)!",l3Current));
812 diCurrent = TMath::Abs(diCurrent);
813 if (TMath::Abs(diCurrent-diNominalCurrent)/diNominalCurrent < tolerance) {
814 // 3% current tolerance...
818 if (diCurrent < zero) { // some small current..
822 AliError(Form("Wrong dipole current (%f A)!",diCurrent));
826 delete fForcedFieldMap;
828 new AliMagWrapCheb("B field map ",s,2,factor,10.,map,dipoleON,path);
830 fForcedFieldMap->Print();
832 AliTracker::SetFieldMap(fForcedFieldMap,fUniformField);
838 Bool_t AliReconstruction::InitGRP() {
839 //------------------------------------
840 // Initialization of the GRP entry
841 //------------------------------------
842 AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/GRP/Data");
845 fGRPData = dynamic_cast<TMap*>(entry->GetObject());
847 AliCDBManager::Instance()->UnloadFromCache("GRP/GRP/Data");
851 AliError("No GRP entry found in OCDB!");
855 TObjString *lhcState=
856 dynamic_cast<TObjString*>(fGRPData->GetValue("fLHCState"));
858 AliError("GRP/GRP/Data entry: missing value for the LHC state ! Using UNKNOWN");
861 TObjString *beamType=
862 dynamic_cast<TObjString*>(fGRPData->GetValue("fAliceBeamType"));
864 AliError("GRP/GRP/Data entry: missing value for the beam type ! Using UNKNOWN");
867 TObjString *beamEnergyStr=
868 dynamic_cast<TObjString*>(fGRPData->GetValue("fAliceBeamEnergy"));
869 if (!beamEnergyStr) {
870 AliError("GRP/GRP/Data entry: missing value for the beam energy ! Using 0");
874 dynamic_cast<TObjString*>(fGRPData->GetValue("fRunType"));
876 AliError("GRP/GRP/Data entry: missing value for the run type ! Using UNKNOWN");
879 TObjString *activeDetectors=
880 dynamic_cast<TObjString*>(fGRPData->GetValue("fDetectorMask"));
881 if (!activeDetectors) {
882 AliError("GRP/GRP/Data entry: missing value for the detector mask ! Using 1074790399");
885 fRunInfo = new AliRunInfo(lhcState ? lhcState->GetString().Data() : "UNKNOWN",
886 beamType ? beamType->GetString().Data() : "UNKNOWN",
887 beamEnergyStr ? beamEnergyStr->GetString().Atof() : 0,
888 runType ? runType->GetString().Data() : "UNKNOWN",
889 activeDetectors ? activeDetectors->GetString().Atoi() : 1074790399);
891 // Process the list of active detectors
892 if (activeDetectors && activeDetectors->GetString().IsDigit()) {
893 UInt_t detMask = activeDetectors->GetString().Atoi();
894 fLoadCDB.Form("%s %s %s %s",
895 fRunLocalReconstruction.Data(),
898 fQADetectors.Data());
899 fRunLocalReconstruction = MatchDetectorList(fRunLocalReconstruction,detMask);
900 fRunTracking = MatchDetectorList(fRunTracking,detMask);
901 fFillESD = MatchDetectorList(fFillESD,detMask);
902 fQADetectors = MatchDetectorList(fQADetectors,detMask);
903 fLoadCDB = MatchDetectorList(fLoadCDB,detMask);
906 AliInfo("===================================================================================");
907 AliInfo(Form("Running local reconstruction for detectors: %s",fRunLocalReconstruction.Data()));
908 AliInfo(Form("Running tracking for detectors: %s",fRunTracking.Data()));
909 AliInfo(Form("Filling ESD for detectors: %s",fFillESD.Data()));
910 AliInfo(Form("Quality assurance is active for detectors: %s",fQADetectors.Data()));
911 AliInfo(Form("CDB and reconstruction parameters are loaded for detectors: %s",fLoadCDB.Data()));
912 AliInfo("===================================================================================");
914 //*** Dealing with the magnetic field map
915 if (AliTracker::GetFieldMap()) {
916 AliInfo("Running with the externally set B field !");
918 // Construct the field map out of the information retrieved from GRP.
923 TObjString *l3Current=
924 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Current"));
926 AliError("GRP/GRP/Data entry: missing value for the L3 current !");
929 TObjString *l3Polarity=
930 dynamic_cast<TObjString*>(fGRPData->GetValue("fL3Polarity"));
932 AliError("GRP/GRP/Data entry: missing value for the L3 polarity !");
937 TObjString *diCurrent=
938 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipoleCurrent"));
940 AliError("GRP/GRP/Data entry: missing value for the dipole current !");
943 TObjString *diPolarity=
944 dynamic_cast<TObjString*>(fGRPData->GetValue("fDipolePolarity"));
946 AliError("GRP/GRP/Data entry: missing value for the dipole polarity !");
951 Float_t l3Cur=TMath::Abs(atof(l3Current->GetName()));
952 Float_t diCur=TMath::Abs(atof(diCurrent->GetName()));
953 Float_t l3Pol=atof(l3Polarity->GetName());
955 if (l3Pol != 0.) factor=-1.;
958 if (!SetFieldMap(l3Cur, diCur, factor)) {
959 AliFatal("Failed to creat a B field map ! Exiting...");
961 AliInfo("Running with the B field constructed out of GRP !");
964 AliFatal("B field is neither set nor constructed from GRP ! Exitig...");
969 //*** Get the diamond profile from OCDB
970 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertex");
972 if (fMeanVertexConstraint)
973 fDiamondProfile = dynamic_cast<AliESDVertex*> (entry->GetObject());
975 AliError("No diamond profile found in OCDB!");
978 entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexTPC");
980 if (fMeanVertexConstraint)
981 fDiamondProfileTPC = dynamic_cast<AliESDVertex*> (entry->GetObject());
983 AliError("No diamond profile found in OCDB!");
989 //_____________________________________________________________________________
990 Bool_t AliReconstruction::LoadCDB()
992 AliCodeTimerAuto("");
994 AliCDBManager::Instance()->Get("GRP/CTP/Config");
996 TString detStr = fLoadCDB;
997 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
998 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
999 AliCDBManager::Instance()->GetAll(Form("%s/Calib/*",fgkDetectorName[iDet]));
1004 //_____________________________________________________________________________
1005 Bool_t AliReconstruction::Run(const char* input)
1008 AliCodeTimerAuto("");
1011 if (GetAbort() != TSelector::kContinue) return kFALSE;
1013 TChain *chain = NULL;
1014 if (fRawReader && (chain = fRawReader->GetChain())) {
1017 gProof->AddInput(this);
1019 outputFile.SetProtocol("root",kTRUE);
1020 outputFile.SetHost(gSystem->HostName());
1021 outputFile.SetFile(Form("%s/AliESDs.root",gSystem->pwd()));
1022 AliInfo(Form("Output file with ESDs is %s",outputFile.GetUrl()));
1023 gProof->AddInput(new TNamed("PROOF_OUTPUTFILE",outputFile.GetUrl()));
1025 chain->Process("AliReconstruction");
1028 chain->Process(this);
1033 if (GetAbort() != TSelector::kContinue) return kFALSE;
1035 if (GetAbort() != TSelector::kContinue) return kFALSE;
1036 //******* The loop over events
1038 while ((iEvent < fRunLoader->GetNumberOfEvents()) ||
1039 (fRawReader && fRawReader->NextEvent())) {
1040 if (!ProcessEvent(iEvent)) {
1041 Abort("ProcessEvent",TSelector::kAbortFile);
1047 if (GetAbort() != TSelector::kContinue) return kFALSE;
1049 if (GetAbort() != TSelector::kContinue) return kFALSE;
1055 //_____________________________________________________________________________
1056 void AliReconstruction::InitRawReader(const char* input)
1058 AliCodeTimerAuto("");
1060 // Init raw-reader and
1061 // set the input in case of raw data
1062 if (input) fRawInput = input;
1063 fRawReader = AliRawReader::Create(fRawInput.Data());
1065 AliInfo("Reconstruction will run over digits");
1067 if (!fEquipIdMap.IsNull() && fRawReader)
1068 fRawReader->LoadEquipmentIdsMap(fEquipIdMap);
1070 if (!fUseHLTData.IsNull()) {
1071 // create the RawReaderHLT which performs redirection of HLT input data for
1072 // the specified detectors
1073 AliRawReader* pRawReader=AliRawHLTManager::CreateRawReaderHLT(fRawReader, fUseHLTData.Data());
1075 fParentRawReader=fRawReader;
1076 fRawReader=pRawReader;
1078 AliError(Form("can not create Raw Reader for HLT input %s", fUseHLTData.Data()));
1081 AliSysInfo::AddStamp("CreateRawReader");
1084 //_____________________________________________________________________________
1085 void AliReconstruction::InitRun(const char* input)
1087 // Initialization of raw-reader,
1088 // run number, CDB etc.
1089 AliCodeTimerAuto("");
1090 AliSysInfo::AddStamp("Start");
1092 // Initialize raw-reader if any
1093 InitRawReader(input);
1095 // Initialize the CDB storage
1098 // Set run number in CDBManager (if it is not already set by the user)
1099 if (!SetRunNumberFromData()) {
1100 Abort("SetRunNumberFromData", TSelector::kAbortProcess);
1104 // Set CDB lock: from now on it is forbidden to reset the run number
1105 // or the default storage or to activate any further storage!
1110 //_____________________________________________________________________________
1111 void AliReconstruction::Begin(TTree *)
1113 // Initialize AlReconstruction before
1114 // going into the event loop
1115 // Should follow the TSelector convention
1116 // i.e. initialize only the object on the client side
1117 AliCodeTimerAuto("");
1119 AliReconstruction *reco = NULL;
1121 if ((reco = (AliReconstruction*)fInput->FindObject("AliReconstruction"))) {
1124 AliSysInfo::AddStamp("ReadInputInBegin");
1127 // Import ideal TGeo geometry and apply misalignment
1129 TString geom(gSystem->DirName(fGAliceFileName));
1130 geom += "/geometry.root";
1131 AliGeomManager::LoadGeometry(geom.Data());
1133 Abort("LoadGeometry", TSelector::kAbortProcess);
1136 AliSysInfo::AddStamp("LoadGeom");
1137 TString detsToCheck=fRunLocalReconstruction;
1138 if(!AliGeomManager::CheckSymNamesLUT(detsToCheck.Data())) {
1139 Abort("CheckSymNamesLUT", TSelector::kAbortProcess);
1142 AliSysInfo::AddStamp("CheckGeom");
1145 if (!MisalignGeometry(fLoadAlignData)) {
1146 Abort("MisalignGeometry", TSelector::kAbortProcess);
1149 AliCDBManager::Instance()->UnloadFromCache("GRP/Geometry/Data");
1150 AliSysInfo::AddStamp("MisalignGeom");
1153 Abort("InitGRP", TSelector::kAbortProcess);
1156 AliSysInfo::AddStamp("InitGRP");
1159 Abort("LoadCDB", TSelector::kAbortProcess);
1162 AliSysInfo::AddStamp("LoadCDB");
1164 // Read the reconstruction parameters from OCDB
1165 if (!InitRecoParams()) {
1166 AliWarning("Not all detectors have correct RecoParam objects initialized");
1168 AliSysInfo::AddStamp("InitRecoParams");
1171 if (reco) *reco = *this;
1172 fInput->Add(gGeoManager);
1174 fInput->Add(const_cast<TMap*>(AliCDBManager::Instance()->GetEntryCache()));
1175 fInput->Add(new TParameter<Int_t>("RunNumber",AliCDBManager::Instance()->GetRun()));
1176 AliMagF *magFieldMap = (AliMagF*)AliTracker::GetFieldMap();
1177 magFieldMap->SetName("MagneticFieldMap");
1178 fInput->Add(magFieldMap);
1183 //_____________________________________________________________________________
1184 void AliReconstruction::SlaveBegin(TTree*)
1186 // Initialization related to run-loader,
1187 // vertexer, trackers, recontructors
1188 // In proof mode it is executed on the slave
1189 AliCodeTimerAuto("");
1191 TProofOutputFile *outProofFile = NULL;
1193 if (AliReconstruction *reco = (AliReconstruction*)fInput->FindObject("AliReconstruction")) {
1196 if (TGeoManager *tgeo = (TGeoManager*)fInput->FindObject("Geometry")) {
1198 AliGeomManager::SetGeometry(tgeo);
1200 if (TMap *entryCache = (TMap*)fInput->FindObject("CDBEntryCache")) {
1201 Int_t runNumber = -1;
1202 if (TProof::GetParameter(fInput,"RunNumber",runNumber) == 0) {
1203 AliCDBManager *man = AliCDBManager::Instance(entryCache,runNumber);
1204 man->SetCacheFlag(kTRUE);
1205 man->SetLock(kTRUE);
1209 if (AliMagF *map = (AliMagF*)fInput->FindObject("MagneticFieldMap")) {
1210 AliTracker::SetFieldMap(map,fUniformField);
1212 if (TNamed *outputFileName = (TNamed *) fInput->FindObject("PROOF_OUTPUTFILE")) {
1213 outProofFile = new TProofOutputFile(gSystem->BaseName(TUrl(outputFileName->GetTitle()).GetFile()));
1214 outProofFile->SetOutputFileName(outputFileName->GetTitle());
1215 fOutput->Add(outProofFile);
1217 AliSysInfo::AddStamp("ReadInputInSlaveBegin");
1220 // get the run loader
1221 if (!InitRunLoader()) {
1222 Abort("InitRunLoader", TSelector::kAbortProcess);
1225 AliSysInfo::AddStamp("LoadLoader");
1227 ftVertexer = new AliVertexerTracks(AliTracker::GetBz());
1228 if(fDiamondProfile && fMeanVertexConstraint) ftVertexer->SetVtxStart(fDiamondProfile);
1231 if (fRunVertexFinder && !CreateVertexer()) {
1232 Abort("CreateVertexer", TSelector::kAbortProcess);
1235 AliSysInfo::AddStamp("CreateVertexer");
1238 if (!fRunTracking.IsNull() && !CreateTrackers(fRunTracking)) {
1239 Abort("CreateTrackers", TSelector::kAbortProcess);
1242 AliSysInfo::AddStamp("CreateTrackers");
1244 // create the ESD output file and tree
1245 if (!outProofFile) {
1246 ffile = TFile::Open("AliESDs.root", "RECREATE");
1247 ffile->SetCompressionLevel(2);
1248 if (!ffile->IsOpen()) {
1249 Abort("OpenESDFile", TSelector::kAbortProcess);
1254 if (!(ffile = outProofFile->OpenFile("RECREATE"))) {
1255 Abort(Form("Problems opening output PROOF file: %s/%s",
1256 outProofFile->GetDir(), outProofFile->GetFileName()),
1257 TSelector::kAbortProcess);
1262 ftree = new TTree("esdTree", "Tree with ESD objects");
1263 fesd = new AliESDEvent();
1264 fesd->CreateStdContent();
1265 fesd->WriteToTree(ftree);
1267 fhlttree = new TTree("HLTesdTree", "Tree with HLT ESD objects");
1268 fhltesd = new AliESDEvent();
1269 fhltesd->CreateStdContent();
1270 fhltesd->WriteToTree(fhlttree);
1273 if (fWriteESDfriend) {
1274 fesdf = new AliESDfriend();
1275 TBranch *br=ftree->Branch("ESDfriend.","AliESDfriend", &fesdf);
1276 br->SetFile("AliESDfriends.root");
1277 fesd->AddObject(fesdf);
1280 ProcInfo_t ProcInfo;
1281 gSystem->GetProcInfo(&ProcInfo);
1282 AliInfo(Form("Current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1285 //Initialize the QA and start of cycle
1287 fQASteer = new AliQADataMakerSteer("rec") ;
1288 fQASteer->SetActiveDetectors(fQADetectors) ;
1289 for (Int_t det = 0 ; det < fgkNDetectors ; det++)
1290 fQASteer->SetCycleLength(AliQA::DETECTORINDEX_t(det), fQACycles[det]) ;
1291 if (!fRawReader && fQATasks.Contains(AliQA::kRAWS))
1292 fQATasks.ReplaceAll(Form("%d",AliQA::kRAWS), "") ;
1293 fQASteer->SetTasks(fQATasks) ;
1294 fQASteer->InitQADataMaker(AliCDBManager::Instance()->GetRun(), fRecoParam) ;
1298 Bool_t sameCycle = kFALSE ;
1299 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1300 AliInfo(Form("Initializing the global QA data maker"));
1301 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS))) {
1302 qadm->StartOfCycle(AliQA::kRECPOINTS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1303 TObjArray *arr=qadm->Init(AliQA::kRECPOINTS);
1304 AliTracker::SetResidualsArray(arr);
1307 if (fQATasks.Contains(Form("%d", AliQA::kESDS))) {
1308 qadm->StartOfCycle(AliQA::kESDS, AliCDBManager::Instance()->GetRun(), sameCycle) ;
1309 qadm->Init(AliQA::kESDS);
1313 //Initialize the Plane Efficiency framework
1314 if (fRunPlaneEff && !InitPlaneEff()) {
1315 Abort("InitPlaneEff", TSelector::kAbortProcess);
1319 if (strcmp(gProgName,"alieve") == 0)
1320 fRunAliEVE = InitAliEVE();
1325 //_____________________________________________________________________________
1326 Bool_t AliReconstruction::Process(Long64_t entry)
1328 // run the reconstruction over a single entry
1329 // from the chain with raw data
1330 AliCodeTimerAuto("");
1332 TTree *currTree = fChain->GetTree();
1333 AliRawEvent *event = new AliRawEvent;
1334 currTree->SetBranchAddress("rawevent",&event);
1335 currTree->GetEntry(entry);
1336 fRawReader = new AliRawReaderRoot(event);
1337 fStatus = ProcessEvent(fRunLoader->GetNumberOfEvents());
1345 //_____________________________________________________________________________
1346 void AliReconstruction::Init(TTree *tree)
1349 AliError("The input tree is not found!");
1355 //_____________________________________________________________________________
1356 Bool_t AliReconstruction::ProcessEvent(Int_t iEvent)
1358 // run the reconstruction over a single event
1359 // The event loop is steered in Run method
1361 AliCodeTimerAuto("");
1363 if (iEvent >= fRunLoader->GetNumberOfEvents()) {
1364 fRunLoader->SetEventNumber(iEvent);
1365 fRunLoader->GetHeader()->Reset(fRawReader->GetRunNumber(),
1367 fRunLoader->TreeE()->Fill();
1368 if (fRawReader && fRawReader->UseAutoSaveESD())
1369 fRunLoader->TreeE()->AutoSave("SaveSelf");
1372 if ((iEvent < fFirstEvent) || ((fLastEvent >= 0) && (iEvent > fLastEvent))) {
1376 AliInfo(Form("processing event %d", iEvent));
1378 // Fill Event-info object
1380 fRecoParam.SetEventSpecie(fRunInfo,fEventInfo);
1382 fRunLoader->GetEvent(iEvent);
1386 fQASteer->RunOneEvent(fRawReader) ;
1388 // local single event reconstruction
1389 if (!fRunLocalReconstruction.IsNull()) {
1390 TString detectors=fRunLocalReconstruction;
1391 // run HLT event reconstruction first
1392 // ;-( IsSelected changes the string
1393 if (IsSelected("HLT", detectors) &&
1394 !RunLocalEventReconstruction("HLT")) {
1395 if (fStopOnError) {CleanUp(); return kFALSE;}
1397 detectors=fRunLocalReconstruction;
1398 detectors.ReplaceAll("HLT", "");
1399 if (!RunLocalEventReconstruction(detectors)) {
1400 if (fStopOnError) {CleanUp(); return kFALSE;}
1404 fesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1405 fhltesd->SetRunNumber(fRunLoader->GetHeader()->GetRun());
1406 fesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1407 fhltesd->SetEventNumberInFile(fRunLoader->GetHeader()->GetEventNrInRun());
1409 // Set magnetic field from the tracker
1410 fesd->SetMagneticField(AliTracker::GetBz());
1411 fhltesd->SetMagneticField(AliTracker::GetBz());
1415 // Fill raw-data error log into the ESD
1416 if (fRawReader) FillRawDataErrorLog(iEvent,fesd);
1419 if (fRunVertexFinder) {
1420 if (!RunVertexFinder(fesd)) {
1421 if (fStopOnError) {CleanUp(); return kFALSE;}
1426 if (!fRunTracking.IsNull()) {
1427 if (fRunMuonTracking) {
1428 if (!RunMuonTracking(fesd)) {
1429 if (fStopOnError) {CleanUp(); return kFALSE;}
1435 if (!fRunTracking.IsNull()) {
1436 if (!RunTracking(fesd)) {
1437 if (fStopOnError) {CleanUp(); return kFALSE;}
1442 if (!fFillESD.IsNull()) {
1443 TString detectors=fFillESD;
1444 // run HLT first and on hltesd
1445 // ;-( IsSelected changes the string
1446 if (IsSelected("HLT", detectors) &&
1447 !FillESD(fhltesd, "HLT")) {
1448 if (fStopOnError) {CleanUp(); return kFALSE;}
1451 // Temporary fix to avoid problems with HLT that overwrites the offline ESDs
1452 if (detectors.Contains("ALL")) {
1454 for (Int_t idet=0; idet<fgkNDetectors; ++idet){
1455 detectors += fgkDetectorName[idet];
1459 detectors.ReplaceAll("HLT", "");
1460 if (!FillESD(fesd, detectors)) {
1461 if (fStopOnError) {CleanUp(); return kFALSE;}
1465 // fill Event header information from the RawEventHeader
1466 if (fRawReader){FillRawEventHeaderESD(fesd);}
1469 AliESDpid::MakePID(fesd);
1471 if (fFillTriggerESD) {
1472 if (!FillTriggerESD(fesd)) {
1473 if (fStopOnError) {CleanUp(); return kFALSE;}
1480 // Propagate track to the beam pipe (if not already done by ITS)
1482 const Int_t ntracks = fesd->GetNumberOfTracks();
1483 const Double_t kBz = fesd->GetMagneticField();
1484 const Double_t kRadius = 2.8; //something less than the beam pipe radius
1487 UShort_t *selectedIdx=new UShort_t[ntracks];
1489 for (Int_t itrack=0; itrack<ntracks; itrack++){
1490 const Double_t kMaxStep = 5; //max step over the material
1493 AliESDtrack *track = fesd->GetTrack(itrack);
1494 if (!track) continue;
1496 AliExternalTrackParam *tpcTrack =
1497 (AliExternalTrackParam *)track->GetTPCInnerParam();
1501 PropagateTrackTo(tpcTrack,kRadius,track->GetMass(),kMaxStep,kTRUE);
1504 Int_t n=trkArray.GetEntriesFast();
1505 selectedIdx[n]=track->GetID();
1506 trkArray.AddLast(tpcTrack);
1509 //Tracks refitted by ITS should already be at the SPD vertex
1510 if (track->IsOn(AliESDtrack::kITSrefit)) continue;
1513 PropagateTrackTo(track,kRadius,track->GetMass(),kMaxStep,kTRUE);
1514 track->RelateToVertex(fesd->GetPrimaryVertexSPD(), kBz, kVeryBig);
1519 // Improve the reconstructed primary vertex position using the tracks
1521 Bool_t runVertexFinderTracks = fRunVertexFinderTracks;
1522 if(fesd->GetPrimaryVertexSPD()) {
1523 TString vtitle = fesd->GetPrimaryVertexSPD()->GetTitle();
1524 if(vtitle.Contains("cosmics")) {
1525 runVertexFinderTracks=kFALSE;
1528 if (runVertexFinderTracks) {
1529 // TPC + ITS primary vertex
1530 ftVertexer->SetITSMode();
1531 if(fDiamondProfile && fMeanVertexConstraint) {
1532 ftVertexer->SetVtxStart(fDiamondProfile);
1534 ftVertexer->SetConstraintOff();
1536 AliESDVertex *pvtx=ftVertexer->FindPrimaryVertex(fesd);
1538 if (pvtx->GetStatus()) {
1539 fesd->SetPrimaryVertex(pvtx);
1540 for (Int_t i=0; i<ntracks; i++) {
1541 AliESDtrack *t = fesd->GetTrack(i);
1542 t->RelateToVertex(pvtx, kBz, kVeryBig);
1547 // TPC-only primary vertex
1548 ftVertexer->SetTPCMode();
1549 if(fDiamondProfileTPC && fMeanVertexConstraint) {
1550 ftVertexer->SetVtxStart(fDiamondProfileTPC);
1552 ftVertexer->SetConstraintOff();
1554 pvtx=ftVertexer->FindPrimaryVertex(&trkArray,selectedIdx);
1556 if (pvtx->GetStatus()) {
1557 fesd->SetPrimaryVertexTPC(pvtx);
1558 for (Int_t i=0; i<ntracks; i++) {
1559 AliESDtrack *t = fesd->GetTrack(i);
1560 t->RelateToVertexTPC(pvtx, kBz, kVeryBig);
1566 delete[] selectedIdx;
1568 if(fDiamondProfile) fesd->SetDiamond(fDiamondProfile);
1573 AliV0vertexer vtxer;
1574 vtxer.Tracks2V0vertices(fesd);
1576 if (fRunCascadeFinder) {
1578 AliCascadeVertexer cvtxer;
1579 cvtxer.V0sTracks2CascadeVertices(fesd);
1584 if (fCleanESD) CleanESD(fesd);
1587 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1588 if (qadm && fQATasks.Contains(Form("%d", AliQA::kESDS)))
1589 qadm->Exec(AliQA::kESDS, fesd);
1592 if (fWriteESDfriend) {
1593 fesdf->~AliESDfriend();
1594 new (fesdf) AliESDfriend(); // Reset...
1595 fesd->GetESDfriend(fesdf);
1599 // Auto-save the ESD tree in case of prompt reco @P2
1600 if (fRawReader && fRawReader->UseAutoSaveESD())
1601 ftree->AutoSave("SaveSelf");
1607 if (fRunAliEVE) RunAliEVE();
1611 if (fWriteESDfriend) {
1612 fesdf->~AliESDfriend();
1613 new (fesdf) AliESDfriend(); // Reset...
1616 ProcInfo_t ProcInfo;
1617 gSystem->GetProcInfo(&ProcInfo);
1618 AliInfo(Form("Event %d -> Current memory usage %d %d",iEvent, ProcInfo.fMemResident, ProcInfo.fMemVirtual));
1621 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1622 if (fReconstructor[iDet])
1623 fReconstructor[iDet]->SetRecoParam(NULL);
1626 fQASteer->Increment() ;
1630 //_____________________________________________________________________________
1631 void AliReconstruction::SlaveTerminate()
1633 // Finalize the run on the slave side
1634 // Called after the exit
1635 // from the event loop
1636 AliCodeTimerAuto("");
1638 if (fIsNewRunLoader) { // galice.root didn't exist
1639 fRunLoader->WriteHeader("OVERWRITE");
1640 fRunLoader->CdGAFile();
1641 fRunLoader->Write(0, TObject::kOverwrite);
1644 ftree->GetUserInfo()->Add(fesd);
1645 fhlttree->GetUserInfo()->Add(fhltesd);
1647 const TMap *cdbMap = AliCDBManager::Instance()->GetStorageMap();
1648 const TList *cdbList = AliCDBManager::Instance()->GetRetrievedIds();
1650 TMap *cdbMapCopy = new TMap(cdbMap->GetEntries());
1651 cdbMapCopy->SetOwner(1);
1652 cdbMapCopy->SetName("cdbMap");
1653 TIter iter(cdbMap->GetTable());
1656 while((pair = dynamic_cast<TPair*> (iter.Next()))){
1657 TObjString* keyStr = dynamic_cast<TObjString*> (pair->Key());
1658 TObjString* valStr = dynamic_cast<TObjString*> (pair->Value());
1659 cdbMapCopy->Add(new TObjString(keyStr->GetName()), new TObjString(valStr->GetName()));
1662 TList *cdbListCopy = new TList();
1663 cdbListCopy->SetOwner(1);
1664 cdbListCopy->SetName("cdbList");
1666 TIter iter2(cdbList);
1669 while((id = dynamic_cast<AliCDBId*> (iter2.Next()))){
1670 cdbListCopy->Add(new TObjString(id->ToString().Data()));
1673 ftree->GetUserInfo()->Add(cdbMapCopy);
1674 ftree->GetUserInfo()->Add(cdbListCopy);
1677 if(fESDPar.Contains("ESD.par")){
1678 AliInfo("Attaching ESD.par to Tree");
1679 TNamed *fn = CopyFileToTNamed(fESDPar.Data(),"ESD.par");
1680 ftree->GetUserInfo()->Add(fn);
1686 if (fWriteESDfriend)
1687 ftree->SetBranchStatus("ESDfriend*",0);
1688 // we want to have only one tree version number
1689 ftree->Write(ftree->GetName(),TObject::kOverwrite);
1692 // Finish with Plane Efficiency evaluation: before of CleanUp !!!
1693 if (fRunPlaneEff && !FinishPlaneEff()) {
1694 AliWarning("Finish PlaneEff evaluation failed");
1697 // End of cycle for the in-loop
1699 fQASteer->RunOneEvent(fesd) ;
1700 fQASteer->EndOfCycle() ;
1703 AliQADataMaker *qadm = fQASteer->GetQADataMaker(AliQA::kGLOBAL);
1705 if (fQATasks.Contains(Form("%d", AliQA::kRECPOINTS)))
1706 qadm->EndOfCycle(AliQA::kRECPOINTS);
1707 if (fQATasks.Contains(Form("%d", AliQA::kESDS)))
1708 qadm->EndOfCycle(AliQA::kESDS);
1716 //_____________________________________________________________________________
1717 void AliReconstruction::Terminate()
1719 // Create tags for the events in the ESD tree (the ESD tree is always present)
1720 // In case of empty events the tags will contain dummy values
1721 AliCodeTimerAuto("");
1723 AliESDTagCreator *esdtagCreator = new AliESDTagCreator();
1724 esdtagCreator->CreateESDTags(fFirstEvent,fLastEvent,fGRPData);
1726 // Cleanup of CDB manager: cache and active storages!
1727 AliCDBManager::Instance()->ClearCache();
1730 //_____________________________________________________________________________
1731 Bool_t AliReconstruction::RunLocalEventReconstruction(const TString& detectors)
1733 // run the local reconstruction
1735 static Int_t eventNr=0;
1736 AliCodeTimerAuto("")
1738 TString detStr = detectors;
1739 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1740 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
1741 AliReconstructor* reconstructor = GetReconstructor(iDet);
1742 if (!reconstructor) continue;
1743 AliLoader* loader = fLoader[iDet];
1744 // Matthias April 2008: temporary fix to run HLT reconstruction
1745 // although the HLT loader is missing
1746 if (strcmp(fgkDetectorName[iDet], "HLT")==0) {
1748 reconstructor->Reconstruct(fRawReader, NULL);
1751 reconstructor->Reconstruct(dummy, NULL);
1756 AliWarning(Form("No loader is defined for %s!",fgkDetectorName[iDet]));
1759 // conversion of digits
1760 if (fRawReader && reconstructor->HasDigitConversion()) {
1761 AliInfo(Form("converting raw data digits into root objects for %s",
1762 fgkDetectorName[iDet]));
1763 // AliCodeTimerAuto(Form("converting raw data digits into root objects for %s",
1764 // fgkDetectorName[iDet]));
1765 loader->LoadDigits("update");
1766 loader->CleanDigits();
1767 loader->MakeDigitsContainer();
1768 TTree* digitsTree = loader->TreeD();
1769 reconstructor->ConvertDigits(fRawReader, digitsTree);
1770 loader->WriteDigits("OVERWRITE");
1771 loader->UnloadDigits();
1773 // local reconstruction
1774 AliInfo(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1775 //AliCodeTimerAuto(Form("running reconstruction for %s", fgkDetectorName[iDet]));
1776 loader->LoadRecPoints("update");
1777 loader->CleanRecPoints();
1778 loader->MakeRecPointsContainer();
1779 TTree* clustersTree = loader->TreeR();
1780 if (fRawReader && !reconstructor->HasDigitConversion()) {
1781 reconstructor->Reconstruct(fRawReader, clustersTree);
1783 loader->LoadDigits("read");
1784 TTree* digitsTree = loader->TreeD();
1786 AliError(Form("Can't get the %s digits tree", fgkDetectorName[iDet]));
1787 if (fStopOnError) return kFALSE;
1789 reconstructor->Reconstruct(digitsTree, clustersTree);
1791 loader->UnloadDigits();
1794 TString detQAStr(fQADetectors) ;
1796 fQASteer->RunOneEventInOneDetector(iDet, clustersTree) ;
1798 loader->WriteRecPoints("OVERWRITE");
1799 loader->UnloadRecPoints();
1800 AliSysInfo::AddStamp(Form("LRec%s_%d",fgkDetectorName[iDet],eventNr), iDet,1,eventNr);
1802 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
1803 AliError(Form("the following detectors were not found: %s",
1805 if (fStopOnError) return kFALSE;
1811 //_____________________________________________________________________________
1812 Bool_t AliReconstruction::RunVertexFinder(AliESDEvent*& esd)
1814 // run the barrel tracking
1816 AliCodeTimerAuto("")
1818 AliESDVertex* vertex = NULL;
1819 Double_t vtxPos[3] = {0, 0, 0};
1820 Double_t vtxErr[3] = {0.07, 0.07, 0.1};
1821 TArrayF mcVertex(3);
1822 if (fRunLoader->GetHeader() && fRunLoader->GetHeader()->GenEventHeader()) {
1823 fRunLoader->GetHeader()->GenEventHeader()->PrimaryVertex(mcVertex);
1824 for (Int_t i = 0; i < 3; i++) vtxPos[i] = mcVertex[i];
1828 AliInfo("running the ITS vertex finder");
1830 fLoader[0]->LoadRecPoints();
1831 TTree* cltree = fLoader[0]->TreeR();
1833 if(fDiamondProfile) fVertexer->SetVtxStart(fDiamondProfile);
1834 vertex = fVertexer->FindVertexForCurrentEvent(cltree);
1837 AliError("Can't get the ITS cluster tree");
1839 fLoader[0]->UnloadRecPoints();
1842 AliError("Can't get the ITS loader");
1845 AliWarning("Vertex not found");
1846 vertex = new AliESDVertex();
1847 vertex->SetName("default");
1850 vertex->SetName("reconstructed");
1854 AliInfo("getting the primary vertex from MC");
1855 vertex = new AliESDVertex(vtxPos, vtxErr);
1859 vertex->GetXYZ(vtxPos);
1860 vertex->GetSigmaXYZ(vtxErr);
1862 AliWarning("no vertex reconstructed");
1863 vertex = new AliESDVertex(vtxPos, vtxErr);
1865 esd->SetPrimaryVertexSPD(vertex);
1866 // if SPD multiplicity has been determined, it is stored in the ESD
1867 AliMultiplicity *mult = fVertexer->GetMultiplicity();
1868 if(mult)esd->SetMultiplicity(mult);
1870 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
1871 if (fTracker[iDet]) fTracker[iDet]->SetVertex(vtxPos, vtxErr);
1878 //_____________________________________________________________________________
1879 Bool_t AliReconstruction::RunHLTTracking(AliESDEvent*& esd)
1881 // run the HLT barrel tracking
1883 AliCodeTimerAuto("")
1886 AliError("Missing runLoader!");
1890 AliInfo("running HLT tracking");
1892 // Get a pointer to the HLT reconstructor
1893 AliReconstructor *reconstructor = GetReconstructor(fgkNDetectors-1);
1894 if (!reconstructor) return kFALSE;
1897 for (Int_t iDet = 1; iDet >= 0; iDet--) {
1898 TString detName = fgkDetectorName[iDet];
1899 AliDebug(1, Form("%s HLT tracking", detName.Data()));
1900 reconstructor->SetOption(detName.Data());
1901 AliTracker *tracker = reconstructor->CreateTracker();
1903 AliWarning(Form("couldn't create a HLT tracker for %s", detName.Data()));
1904 if (fStopOnError) return kFALSE;
1908 Double_t vtxErr[3]={0.005,0.005,0.010};
1909 const AliESDVertex *vertex = esd->GetVertex();
1910 vertex->GetXYZ(vtxPos);
1911 tracker->SetVertex(vtxPos,vtxErr);
1913 fLoader[iDet]->LoadRecPoints("read");
1914 TTree* tree = fLoader[iDet]->TreeR();
1916 AliError(Form("Can't get the %s cluster tree", detName.Data()));
1919 tracker->LoadClusters(tree);
1921 if (tracker->Clusters2Tracks(esd) != 0) {
1922 AliError(Form("HLT %s Clusters2Tracks failed", fgkDetectorName[iDet]));
1926 tracker->UnloadClusters();
1934 //_____________________________________________________________________________
1935 Bool_t AliReconstruction::RunMuonTracking(AliESDEvent*& esd)
1937 // run the muon spectrometer tracking
1939 AliCodeTimerAuto("")
1942 AliError("Missing runLoader!");
1945 Int_t iDet = 7; // for MUON
1947 AliInfo("is running...");
1949 // Get a pointer to the MUON reconstructor
1950 AliReconstructor *reconstructor = GetReconstructor(iDet);
1951 if (!reconstructor) return kFALSE;
1954 TString detName = fgkDetectorName[iDet];
1955 AliDebug(1, Form("%s tracking", detName.Data()));
1956 AliTracker *tracker = reconstructor->CreateTracker();
1958 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
1963 fLoader[iDet]->LoadRecPoints("read");
1965 tracker->LoadClusters(fLoader[iDet]->TreeR());
1967 Int_t rv = tracker->Clusters2Tracks(esd);
1971 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
1975 fLoader[iDet]->UnloadRecPoints();
1977 tracker->UnloadClusters();
1985 //_____________________________________________________________________________
1986 Bool_t AliReconstruction::RunTracking(AliESDEvent*& esd)
1988 // run the barrel tracking
1989 static Int_t eventNr=0;
1990 AliCodeTimerAuto("")
1992 AliInfo("running tracking");
1994 //Fill the ESD with the T0 info (will be used by the TOF)
1995 if (fReconstructor[11] && fLoader[11]) {
1996 fLoader[11]->LoadRecPoints("READ");
1997 TTree *treeR = fLoader[11]->TreeR();
1998 GetReconstructor(11)->FillESD((TTree *)NULL,treeR,esd);
2001 // pass 1: TPC + ITS inwards
2002 for (Int_t iDet = 1; iDet >= 0; iDet--) {
2003 if (!fTracker[iDet]) continue;
2004 AliDebug(1, Form("%s tracking", fgkDetectorName[iDet]));
2007 fLoader[iDet]->LoadRecPoints("read");
2008 AliSysInfo::AddStamp(Form("RLoadCluster%s_%d",fgkDetectorName[iDet],eventNr),iDet,1, eventNr);
2009 TTree* tree = fLoader[iDet]->TreeR();
2011 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2014 fTracker[iDet]->LoadClusters(tree);
2015 AliSysInfo::AddStamp(Form("TLoadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2017 if (fTracker[iDet]->Clusters2Tracks(esd) != 0) {
2018 AliError(Form("%s Clusters2Tracks failed", fgkDetectorName[iDet]));
2021 // preliminary PID in TPC needed by the ITS tracker
2023 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2024 AliESDpid::MakePID(esd);
2026 AliSysInfo::AddStamp(Form("Tracking0%s_%d",fgkDetectorName[iDet],eventNr), iDet,3,eventNr);
2029 // pass 2: ALL backwards
2031 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2032 if (!fTracker[iDet]) continue;
2033 AliDebug(1, Form("%s back propagation", fgkDetectorName[iDet]));
2036 if (iDet > 1) { // all except ITS, TPC
2038 fLoader[iDet]->LoadRecPoints("read");
2039 AliSysInfo::AddStamp(Form("RLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,1, eventNr);
2040 tree = fLoader[iDet]->TreeR();
2042 AliError(Form("Can't get the %s cluster tree", fgkDetectorName[iDet]));
2045 fTracker[iDet]->LoadClusters(tree);
2046 AliSysInfo::AddStamp(Form("TLoadCluster0%s_%d",fgkDetectorName[iDet],eventNr), iDet,2, eventNr);
2050 if (iDet>1) // start filling residuals for the "outer" detectors
2051 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2053 if (fTracker[iDet]->PropagateBack(esd) != 0) {
2054 AliError(Form("%s backward propagation failed", fgkDetectorName[iDet]));
2059 if (iDet > 3) { // all except ITS, TPC, TRD and TOF
2060 fTracker[iDet]->UnloadClusters();
2061 fLoader[iDet]->UnloadRecPoints();
2063 // updated PID in TPC needed by the ITS tracker -MI
2065 GetReconstructor(1)->FillESD((TTree*)NULL, (TTree*)NULL, esd);
2066 AliESDpid::MakePID(esd);
2068 AliSysInfo::AddStamp(Form("Tracking1%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2070 //stop filling residuals for the "outer" detectors
2071 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2073 // pass 3: TRD + TPC + ITS refit inwards
2075 for (Int_t iDet = 2; iDet >= 0; iDet--) {
2076 if (!fTracker[iDet]) continue;
2077 AliDebug(1, Form("%s inward refit", fgkDetectorName[iDet]));
2080 if (iDet<2) // start filling residuals for TPC and ITS
2081 if (fRunGlobalQA) AliTracker::SetFillResiduals(kTRUE);
2083 if (fTracker[iDet]->RefitInward(esd) != 0) {
2084 AliError(Form("%s inward refit failed", fgkDetectorName[iDet]));
2087 // run postprocessing
2088 if (fTracker[iDet]->PostProcess(esd) != 0) {
2089 AliError(Form("%s postprocessing failed", fgkDetectorName[iDet]));
2092 AliSysInfo::AddStamp(Form("Tracking2%s_%d",fgkDetectorName[iDet],eventNr), iDet,3, eventNr);
2095 // write space-points to the ESD in case alignment data output
2097 if (fWriteAlignmentData)
2098 WriteAlignmentData(esd);
2100 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2101 if (!fTracker[iDet]) continue;
2103 fTracker[iDet]->UnloadClusters();
2104 AliSysInfo::AddStamp(Form("TUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,4, eventNr);
2105 fLoader[iDet]->UnloadRecPoints();
2106 AliSysInfo::AddStamp(Form("RUnloadCluster%s_%d",fgkDetectorName[iDet],eventNr), iDet,5, eventNr);
2108 // stop filling residuals for TPC and ITS
2109 if (fRunGlobalQA) AliTracker::SetFillResiduals(kFALSE);
2115 //_____________________________________________________________________________
2116 Bool_t AliReconstruction::CleanESD(AliESDEvent *esd){
2118 // Remove the data which are not needed for the physics analysis.
2121 Int_t nTracks=esd->GetNumberOfTracks();
2122 Int_t nV0s=esd->GetNumberOfV0s();
2124 (Form("Number of ESD tracks and V0s before cleaning: %d %d",nTracks,nV0s));
2126 Float_t cleanPars[]={fV0DCAmax,fV0CsPmin,fDmax,fZmax};
2127 Bool_t rc=esd->Clean(cleanPars);
2129 nTracks=esd->GetNumberOfTracks();
2130 nV0s=esd->GetNumberOfV0s();
2132 (Form("Number of ESD tracks and V0s after cleaning %d %d",nTracks,nV0s));
2137 //_____________________________________________________________________________
2138 Bool_t AliReconstruction::FillESD(AliESDEvent*& esd, const TString& detectors)
2140 // fill the event summary data
2142 AliCodeTimerAuto("")
2143 static Int_t eventNr=0;
2144 TString detStr = detectors;
2146 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2147 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2148 AliReconstructor* reconstructor = GetReconstructor(iDet);
2149 if (!reconstructor) continue;
2150 AliDebug(1, Form("filling ESD for %s", fgkDetectorName[iDet]));
2151 TTree* clustersTree = NULL;
2152 if (fLoader[iDet]) {
2153 fLoader[iDet]->LoadRecPoints("read");
2154 clustersTree = fLoader[iDet]->TreeR();
2155 if (!clustersTree) {
2156 AliError(Form("Can't get the %s clusters tree",
2157 fgkDetectorName[iDet]));
2158 if (fStopOnError) return kFALSE;
2161 if (fRawReader && !reconstructor->HasDigitConversion()) {
2162 reconstructor->FillESD(fRawReader, clustersTree, esd);
2164 TTree* digitsTree = NULL;
2165 if (fLoader[iDet]) {
2166 fLoader[iDet]->LoadDigits("read");
2167 digitsTree = fLoader[iDet]->TreeD();
2169 AliError(Form("Can't get the %s digits tree",
2170 fgkDetectorName[iDet]));
2171 if (fStopOnError) return kFALSE;
2174 reconstructor->FillESD(digitsTree, clustersTree, esd);
2175 if (fLoader[iDet]) fLoader[iDet]->UnloadDigits();
2177 if (fLoader[iDet]) {
2178 fLoader[iDet]->UnloadRecPoints();
2182 if ((detStr.CompareTo("ALL") != 0) && !detStr.IsNull()) {
2183 AliError(Form("the following detectors were not found: %s",
2185 if (fStopOnError) return kFALSE;
2187 AliSysInfo::AddStamp(Form("FillESD%d",eventNr), 0,1, eventNr);
2192 //_____________________________________________________________________________
2193 Bool_t AliReconstruction::FillTriggerESD(AliESDEvent*& esd)
2195 // Reads the trigger decision which is
2196 // stored in Trigger.root file and fills
2197 // the corresponding esd entries
2199 AliCodeTimerAuto("")
2201 AliInfo("Filling trigger information into the ESD");
2204 AliCTPRawStream input(fRawReader);
2205 if (!input.Next()) {
2206 AliWarning("No valid CTP (trigger) DDL raw data is found ! The trigger info is taken from the event header!");
2209 if (esd->GetTriggerMask() != input.GetClassMask())
2210 AliError(Form("Invalid trigger pattern found in CTP raw-data: %llx %llx",
2211 input.GetClassMask(),esd->GetTriggerMask()));
2212 if (esd->GetOrbitNumber() != input.GetOrbitID())
2213 AliError(Form("Invalid orbit id found in CTP raw-data: %x %x",
2214 input.GetOrbitID(),esd->GetOrbitNumber()));
2215 if (esd->GetBunchCrossNumber() != input.GetBCID())
2216 AliError(Form("Invalid bunch-crossing id found in CTP raw-data: %x %x",
2217 input.GetBCID(),esd->GetBunchCrossNumber()));
2220 // Here one has to add the filling of trigger inputs and
2221 // interaction records
2231 //_____________________________________________________________________________
2232 Bool_t AliReconstruction::FillRawEventHeaderESD(AliESDEvent*& esd)
2235 // Filling information from RawReader Header
2238 if (!fRawReader) return kFALSE;
2240 AliInfo("Filling information from RawReader Header");
2242 esd->SetBunchCrossNumber(fRawReader->GetBCID());
2243 esd->SetOrbitNumber(fRawReader->GetOrbitID());
2244 esd->SetPeriodNumber(fRawReader->GetPeriod());
2246 esd->SetTimeStamp(fRawReader->GetTimestamp());
2247 esd->SetEventType(fRawReader->GetType());
2253 //_____________________________________________________________________________
2254 Bool_t AliReconstruction::IsSelected(TString detName, TString& detectors) const
2256 // check whether detName is contained in detectors
2257 // if yes, it is removed from detectors
2259 // check if all detectors are selected
2260 if ((detectors.CompareTo("ALL") == 0) ||
2261 detectors.BeginsWith("ALL ") ||
2262 detectors.EndsWith(" ALL") ||
2263 detectors.Contains(" ALL ")) {
2268 // search for the given detector
2269 Bool_t result = kFALSE;
2270 if ((detectors.CompareTo(detName) == 0) ||
2271 detectors.BeginsWith(detName+" ") ||
2272 detectors.EndsWith(" "+detName) ||
2273 detectors.Contains(" "+detName+" ")) {
2274 detectors.ReplaceAll(detName, "");
2278 // clean up the detectors string
2279 while (detectors.Contains(" ")) detectors.ReplaceAll(" ", " ");
2280 while (detectors.BeginsWith(" ")) detectors.Remove(0, 1);
2281 while (detectors.EndsWith(" ")) detectors.Remove(detectors.Length()-1, 1);
2286 //_____________________________________________________________________________
2287 Bool_t AliReconstruction::InitRunLoader()
2289 // get or create the run loader
2291 if (gAlice) delete gAlice;
2294 if (!gSystem->AccessPathName(fGAliceFileName.Data())) { // galice.root exists
2295 // load all base libraries to get the loader classes
2296 TString libs = gSystem->GetLibraries();
2297 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2298 TString detName = fgkDetectorName[iDet];
2299 if (detName == "HLT") continue;
2300 if (libs.Contains("lib" + detName + "base.so")) continue;
2301 gSystem->Load("lib" + detName + "base.so");
2303 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data());
2305 AliError(Form("no run loader found in file %s", fGAliceFileName.Data()));
2310 fRunLoader->CdGAFile();
2311 fRunLoader->LoadgAlice();
2313 //PH This is a temporary fix to give access to the kinematics
2314 //PH that is needed for the labels of ITS clusters
2315 fRunLoader->LoadHeader();
2316 fRunLoader->LoadKinematics();
2318 } else { // galice.root does not exist
2320 AliError(Form("the file %s does not exist", fGAliceFileName.Data()));
2322 fRunLoader = AliRunLoader::Open(fGAliceFileName.Data(),
2323 AliConfig::GetDefaultEventFolderName(),
2326 AliError(Form("could not create run loader in file %s",
2327 fGAliceFileName.Data()));
2331 fIsNewRunLoader = kTRUE;
2332 fRunLoader->MakeTree("E");
2334 if (fNumberOfEventsPerFile > 0)
2335 fRunLoader->SetNumberOfEventsPerFile(fNumberOfEventsPerFile);
2337 fRunLoader->SetNumberOfEventsPerFile((UInt_t)-1);
2343 //_____________________________________________________________________________
2344 AliReconstructor* AliReconstruction::GetReconstructor(Int_t iDet)
2346 // get the reconstructor object and the loader for a detector
2348 if (fReconstructor[iDet]) {
2349 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2350 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2351 fReconstructor[iDet]->SetRecoParam(par);
2353 return fReconstructor[iDet];
2356 // load the reconstructor object
2357 TPluginManager* pluginManager = gROOT->GetPluginManager();
2358 TString detName = fgkDetectorName[iDet];
2359 TString recName = "Ali" + detName + "Reconstructor";
2361 if (!fIsNewRunLoader && !fRunLoader->GetLoader(detName+"Loader") && (detName != "HLT")) return NULL;
2363 AliReconstructor* reconstructor = NULL;
2364 // first check if a plugin is defined for the reconstructor
2365 TPluginHandler* pluginHandler =
2366 pluginManager->FindHandler("AliReconstructor", detName);
2367 // if not, add a plugin for it
2368 if (!pluginHandler) {
2369 AliDebug(1, Form("defining plugin for %s", recName.Data()));
2370 TString libs = gSystem->GetLibraries();
2371 if (libs.Contains("lib" + detName + "base.so") ||
2372 (gSystem->Load("lib" + detName + "base.so") >= 0)) {
2373 pluginManager->AddHandler("AliReconstructor", detName,
2374 recName, detName + "rec", recName + "()");
2376 pluginManager->AddHandler("AliReconstructor", detName,
2377 recName, detName, recName + "()");
2379 pluginHandler = pluginManager->FindHandler("AliReconstructor", detName);
2381 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2382 reconstructor = (AliReconstructor*) pluginHandler->ExecPlugin(0);
2384 if (reconstructor) {
2385 TObject* obj = fOptions.FindObject(detName.Data());
2386 if (obj) reconstructor->SetOption(obj->GetTitle());
2387 reconstructor->Init();
2388 fReconstructor[iDet] = reconstructor;
2391 // get or create the loader
2392 if (detName != "HLT") {
2393 fLoader[iDet] = fRunLoader->GetLoader(detName + "Loader");
2394 if (!fLoader[iDet]) {
2395 AliConfig::Instance()
2396 ->CreateDetectorFolders(fRunLoader->GetEventFolder(),
2398 // first check if a plugin is defined for the loader
2400 pluginManager->FindHandler("AliLoader", detName);
2401 // if not, add a plugin for it
2402 if (!pluginHandler) {
2403 TString loaderName = "Ali" + detName + "Loader";
2404 AliDebug(1, Form("defining plugin for %s", loaderName.Data()));
2405 pluginManager->AddHandler("AliLoader", detName,
2406 loaderName, detName + "base",
2407 loaderName + "(const char*, TFolder*)");
2408 pluginHandler = pluginManager->FindHandler("AliLoader", detName);
2410 if (pluginHandler && (pluginHandler->LoadPlugin() == 0)) {
2412 (AliLoader*) pluginHandler->ExecPlugin(2, detName.Data(),
2413 fRunLoader->GetEventFolder());
2415 if (!fLoader[iDet]) { // use default loader
2416 fLoader[iDet] = new AliLoader(detName, fRunLoader->GetEventFolder());
2418 if (!fLoader[iDet]) {
2419 AliWarning(Form("couldn't get loader for %s", detName.Data()));
2420 if (fStopOnError) return NULL;
2422 fRunLoader->AddLoader(fLoader[iDet]);
2423 fRunLoader->CdGAFile();
2424 if (gFile && !gFile->IsWritable()) gFile->ReOpen("UPDATE");
2425 fRunLoader->Write(0, TObject::kOverwrite);
2430 if (fRecoParam.GetDetRecoParamArray(iDet) && !AliReconstructor::GetRecoParam(iDet)) {
2431 const AliDetectorRecoParam *par = fRecoParam.GetDetRecoParam(iDet);
2432 reconstructor->SetRecoParam(par);
2434 return reconstructor;
2437 //_____________________________________________________________________________
2438 Bool_t AliReconstruction::CreateVertexer()
2440 // create the vertexer
2443 AliReconstructor* itsReconstructor = GetReconstructor(0);
2444 if (itsReconstructor) {
2445 fVertexer = itsReconstructor->CreateVertexer();
2448 AliWarning("couldn't create a vertexer for ITS");
2449 if (fStopOnError) return kFALSE;
2455 //_____________________________________________________________________________
2456 Bool_t AliReconstruction::CreateTrackers(const TString& detectors)
2458 // create the trackers
2460 TString detStr = detectors;
2461 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2462 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2463 AliReconstructor* reconstructor = GetReconstructor(iDet);
2464 if (!reconstructor) continue;
2465 TString detName = fgkDetectorName[iDet];
2466 if (detName == "HLT") {
2467 fRunHLTTracking = kTRUE;
2470 if (detName == "MUON") {
2471 fRunMuonTracking = kTRUE;
2476 fTracker[iDet] = reconstructor->CreateTracker();
2477 if (!fTracker[iDet] && (iDet < 7)) {
2478 AliWarning(Form("couldn't create a tracker for %s", detName.Data()));
2479 if (fStopOnError) return kFALSE;
2481 AliSysInfo::AddStamp(Form("LTracker%s",fgkDetectorName[iDet]), iDet,0);
2487 //_____________________________________________________________________________
2488 void AliReconstruction::CleanUp()
2490 // delete trackers and the run loader and close and delete the file
2492 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2493 delete fReconstructor[iDet];
2494 fReconstructor[iDet] = NULL;
2495 fLoader[iDet] = NULL;
2496 delete fTracker[iDet];
2497 fTracker[iDet] = NULL;
2508 if(!(AliCDBManager::Instance()->GetCacheFlag())) {
2509 delete fDiamondProfile;
2510 fDiamondProfile = NULL;
2511 delete fDiamondProfileTPC;
2512 fDiamondProfileTPC = NULL;
2521 delete fParentRawReader;
2522 fParentRawReader=NULL;
2531 void AliReconstruction::WriteAlignmentData(AliESDEvent* esd)
2533 // Write space-points which are then used in the alignment procedures
2534 // For the moment only ITS, TPC, TRD and TOF
2536 Int_t ntracks = esd->GetNumberOfTracks();
2537 for (Int_t itrack = 0; itrack < ntracks; itrack++)
2539 AliESDtrack *track = esd->GetTrack(itrack);
2542 for (Int_t iDet = 3; iDet >= 0; iDet--) {// TOF, TRD, TPC, ITS clusters
2543 nsp += track->GetNcls(iDet);
2545 if (iDet==0) { // ITS "extra" clusters
2546 track->GetClusters(iDet,idx);
2547 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nsp++;
2552 AliTrackPointArray *sp = new AliTrackPointArray(nsp);
2553 track->SetTrackPointArray(sp);
2555 for (Int_t iDet = 3; iDet >= 0; iDet--) {
2556 AliTracker *tracker = fTracker[iDet];
2557 if (!tracker) continue;
2558 Int_t nspdet = track->GetClusters(iDet,idx);
2560 if (iDet==0) // ITS "extra" clusters
2561 for (Int_t i=6; i<12; i++) if(idx[i] >= 0) nspdet++;
2563 if (nspdet <= 0) continue;
2567 while (isp2 < nspdet) {
2568 Bool_t isvalid=kTRUE;
2570 Int_t index=idx[isp++];
2571 if (index < 0) continue;
2573 TString dets = fgkDetectorName[iDet];
2574 if ((fUseTrackingErrorsForAlignment.CompareTo(dets) == 0) ||
2575 fUseTrackingErrorsForAlignment.BeginsWith(dets+" ") ||
2576 fUseTrackingErrorsForAlignment.EndsWith(" "+dets) ||
2577 fUseTrackingErrorsForAlignment.Contains(" "+dets+" ")) {
2578 isvalid = tracker->GetTrackPointTrackingError(index,p,track);
2580 isvalid = tracker->GetTrackPoint(index,p);
2583 if (!isvalid) continue;
2584 sp->AddPoint(isptrack,&p); isptrack++;
2591 //_____________________________________________________________________________
2592 void AliReconstruction::FillRawDataErrorLog(Int_t iEvent, AliESDEvent* esd)
2594 // The method reads the raw-data error log
2595 // accumulated within the rawReader.
2596 // It extracts the raw-data errors related to
2597 // the current event and stores them into
2598 // a TClonesArray inside the esd object.
2600 if (!fRawReader) return;
2602 for(Int_t i = 0; i < fRawReader->GetNumberOfErrorLogs(); i++) {
2604 AliRawDataErrorLog *log = fRawReader->GetErrorLog(i);
2606 if (iEvent != log->GetEventNumber()) continue;
2608 esd->AddRawDataErrorLog(log);
2613 TNamed* AliReconstruction::CopyFileToTNamed(TString fPath,TString pName){
2614 // Dump a file content into a char in TNamed
2616 in.open(fPath.Data(),ios::in | ios::binary|ios::ate);
2617 Int_t kBytes = (Int_t)in.tellg();
2618 printf("Size: %d \n",kBytes);
2621 char* memblock = new char [kBytes];
2622 in.seekg (0, ios::beg);
2623 in.read (memblock, kBytes);
2625 TString fData(memblock,kBytes);
2626 fn = new TNamed(pName,fData);
2627 printf("fData Size: %d \n",fData.Sizeof());
2628 printf("pName Size: %d \n",pName.Sizeof());
2629 printf("fn Size: %d \n",fn->Sizeof());
2633 AliInfo(Form("Could not Open %s\n",fPath.Data()));
2639 void AliReconstruction::TNamedToFile(TTree* fTree, TString pName){
2640 // This is not really needed in AliReconstruction at the moment
2641 // but can serve as a template
2643 TList *fList = fTree->GetUserInfo();
2644 TNamed *fn = (TNamed*)fList->FindObject(pName.Data());
2645 printf("fn Size: %d \n",fn->Sizeof());
2647 TString fTmp(fn->GetName()); // to be 100% sure in principle pName also works
2648 const char* cdata = fn->GetTitle();
2649 printf("fTmp Size %d\n",fTmp.Sizeof());
2651 int size = fn->Sizeof()-fTmp.Sizeof()-sizeof(UChar_t)-sizeof(Int_t); // see dfinition of TString::SizeOf()...
2652 printf("calculated size %d\n",size);
2653 ofstream out(pName.Data(),ios::out | ios::binary);
2654 out.write(cdata,size);
2659 //_____________________________________________________________________________
2660 void AliReconstruction::CheckQA()
2662 // check the QA of SIM for this run and remove the detectors
2663 // with status Fatal
2665 TString newRunLocalReconstruction ;
2666 TString newRunTracking ;
2667 TString newFillESD ;
2669 for (Int_t iDet = 0; iDet < AliQA::kNDET; iDet++) {
2670 TString detName(AliQA::GetDetName(iDet)) ;
2671 AliQA * qa = AliQA::Instance(AliQA::DETECTORINDEX_t(iDet)) ;
2672 if ( qa->IsSet(AliQA::DETECTORINDEX_t(iDet), AliQA::kSIM, AliQA::kFATAL)) {
2673 AliInfo(Form("QA status for %s in Hits and/or SDIGITS and/or Digits was Fatal; No reconstruction performed", detName.Data())) ;
2675 if ( fRunLocalReconstruction.Contains(AliQA::GetDetName(iDet)) ||
2676 fRunLocalReconstruction.Contains("ALL") ) {
2677 newRunLocalReconstruction += detName ;
2678 newRunLocalReconstruction += " " ;
2680 if ( fRunTracking.Contains(AliQA::GetDetName(iDet)) ||
2681 fRunTracking.Contains("ALL") ) {
2682 newRunTracking += detName ;
2683 newRunTracking += " " ;
2685 if ( fFillESD.Contains(AliQA::GetDetName(iDet)) ||
2686 fFillESD.Contains("ALL") ) {
2687 newFillESD += detName ;
2692 fRunLocalReconstruction = newRunLocalReconstruction ;
2693 fRunTracking = newRunTracking ;
2694 fFillESD = newFillESD ;
2697 //_____________________________________________________________________________
2698 Int_t AliReconstruction::GetDetIndex(const char* detector)
2700 // return the detector index corresponding to detector
2702 for (index = 0; index < fgkNDetectors ; index++) {
2703 if ( strcmp(detector, fgkDetectorName[index]) == 0 )
2708 //_____________________________________________________________________________
2709 Bool_t AliReconstruction::FinishPlaneEff() {
2711 // Here execute all the necessary operationis, at the end of the tracking phase,
2712 // in case that evaluation of PlaneEfficiencies was required for some detector.
2713 // E.g., write into a DataBase file the PlaneEfficiency which have been evaluated.
2715 // This Preliminary version works only FOR ITS !!!!!
2716 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2719 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2722 //for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2723 for (Int_t iDet = 0; iDet < 1; iDet++) { // for the time being only ITS
2724 //if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2725 if(fTracker[iDet]) {
2726 AliPlaneEff *planeeff=fTracker[iDet]->GetPlaneEff();
2727 TString name=planeeff->GetName();
2729 TFile* pefile = TFile::Open(name, "RECREATE");
2730 ret=(Bool_t)planeeff->Write();
2732 if(planeeff->GetCreateHistos()) {
2733 TString hname=planeeff->GetName();
2734 hname+="Histo.root";
2735 ret*=planeeff->WriteHistosToFile(hname,"RECREATE");
2741 //_____________________________________________________________________________
2742 Bool_t AliReconstruction::InitPlaneEff() {
2744 // Here execute all the necessary operations, before of the tracking phase,
2745 // for the evaluation of PlaneEfficiencies, in case required for some detectors.
2746 // E.g., read from a DataBase file a first evaluation of the PlaneEfficiency
2747 // which should be updated/recalculated.
2749 // This Preliminary version will work only FOR ITS !!!!!
2750 // other detectors (TOF,TRD, etc. have to develop their specific codes)
2753 // Return: kTRUE if all operations have been done properly, kFALSE otherwise
2755 AliWarning(Form("Implementation of this method not yet done !! Method return kTRUE"));
2759 //_____________________________________________________________________________
2760 Bool_t AliReconstruction::InitAliEVE()
2762 // This method should be called only in case
2763 // AliReconstruction is run
2764 // within the alieve environment.
2765 // It will initialize AliEVE in a way
2766 // so that it can visualize event processed
2767 // by AliReconstruction.
2768 // The return flag shows whenever the
2769 // AliEVE initialization was successful or not.
2772 macroStr.Form("%s/EVE/macros/alieve_online.C",gSystem->ExpandPathName("$ALICE_ROOT"));
2773 AliInfo(Form("Loading AliEVE macro: %s",macroStr.Data()));
2774 if (gROOT->LoadMacro(macroStr.Data()) != 0) return kFALSE;
2776 gROOT->ProcessLine("if (!gAliEveEvent) {gAliEveEvent = new AliEveEventManager();gAliEveEvent->AddNewEventCommand(\"alieve_online_on_new_event()\");gEve->AddEvent(gAliEveEvent);};");
2777 gROOT->ProcessLine("alieve_online_init()");
2782 //_____________________________________________________________________________
2783 void AliReconstruction::RunAliEVE()
2785 // Runs AliEVE visualisation of
2786 // the current event.
2787 // Should be executed only after
2788 // successful initialization of AliEVE.
2790 AliInfo("Running AliEVE...");
2791 gROOT->ProcessLine(Form("gAliEveEvent->SetEvent((AliRunLoader*)%p,(AliRawReader*)%p,(AliESDEvent*)%p);",fRunLoader,fRawReader,fesd));
2795 //_____________________________________________________________________________
2796 Bool_t AliReconstruction::SetRunQA(TString detAndAction)
2798 // Allows to run QA for a selected set of detectors
2799 // and a selected set of tasks among RAWS, RECPOINTS and ESDS
2800 // all selected detectors run the same selected tasks
2802 if (!detAndAction.Contains(":")) {
2803 AliError( Form("%s is a wrong syntax, use \"DetectorList:ActionList\" \n", detAndAction.Data()) ) ;
2807 Int_t colon = detAndAction.Index(":") ;
2808 fQADetectors = detAndAction(0, colon) ;
2809 if (fQADetectors.Contains("ALL") )
2810 fQADetectors = fFillESD ;
2811 fQATasks = detAndAction(colon+1, detAndAction.Sizeof() ) ;
2812 if (fQATasks.Contains("ALL") ) {
2813 fQATasks = Form("%d %d %d", AliQA::kRAWS, AliQA::kRECPOINTS, AliQA::kESDS) ;
2815 fQATasks.ToUpper() ;
2817 if ( fQATasks.Contains("RAW") )
2818 tempo = Form("%d ", AliQA::kRAWS) ;
2819 if ( fQATasks.Contains("RECPOINT") )
2820 tempo += Form("%d ", AliQA::kRECPOINTS) ;
2821 if ( fQATasks.Contains("ESD") )
2822 tempo += Form("%d ", AliQA::kESDS) ;
2824 if (fQATasks.IsNull()) {
2825 AliInfo("No QA requested\n") ;
2830 TString tempo(fQATasks) ;
2831 tempo.ReplaceAll(Form("%d", AliQA::kRAWS), AliQA::GetTaskName(AliQA::kRAWS)) ;
2832 tempo.ReplaceAll(Form("%d", AliQA::kRECPOINTS), AliQA::GetTaskName(AliQA::kRECPOINTS)) ;
2833 tempo.ReplaceAll(Form("%d", AliQA::kESDS), AliQA::GetTaskName(AliQA::kESDS)) ;
2834 AliInfo( Form("QA will be done on \"%s\" for \"%s\"\n", fQADetectors.Data(), tempo.Data()) ) ;
2839 //_____________________________________________________________________________
2840 Bool_t AliReconstruction::InitRecoParams()
2842 // The method accesses OCDB and retrieves all
2843 // the available reco-param objects from there.
2845 Bool_t isOK = kTRUE;
2847 TString detStr = fLoadCDB;
2848 for (Int_t iDet = 0; iDet < fgkNDetectors; iDet++) {
2850 if (!IsSelected(fgkDetectorName[iDet], detStr)) continue;
2852 if (fRecoParam.GetDetRecoParamArray(iDet)) {
2853 AliInfo(Form("Using custom reconstruction parameters for detector %s",fgkDetectorName[iDet]));
2857 AliInfo(Form("Loading reconstruction parameter objects for detector %s",fgkDetectorName[iDet]));
2859 AliCDBPath path(fgkDetectorName[iDet],"Calib","RecoParam");
2860 AliCDBEntry *entry=AliCDBManager::Instance()->Get(path.GetPath());
2862 AliWarning(Form("Couldn't find RecoParam entry in OCDB for detector %s",fgkDetectorName[iDet]));
2866 TObject *recoParamObj = entry->GetObject();
2867 if (dynamic_cast<TObjArray*>(recoParamObj)) {
2868 // The detector has a normal TobjArray of AliDetectorRecoParam objects
2869 // Registering them in AliRecoParam
2870 fRecoParam.AddDetRecoParamArray(iDet,dynamic_cast<TObjArray*>(recoParamObj));
2872 else if (dynamic_cast<AliDetectorRecoParam*>(recoParamObj)) {
2873 // The detector has only onse set of reco parameters
2874 // Registering it in AliRecoParam
2875 AliInfo(Form("Single set of reconstruction parameters found for detector %s",fgkDetectorName[iDet]));
2876 dynamic_cast<AliDetectorRecoParam*>(recoParamObj)->SetAsDefault();
2877 fRecoParam.AddDetRecoParam(iDet,dynamic_cast<AliDetectorRecoParam*>(recoParamObj));
2880 AliError(Form("No valid RecoParam object found in the OCDB for detector %s",fgkDetectorName[iDet]));
2884 AliCDBManager::Instance()->UnloadFromCache(path.GetPath());
2888 if (AliDebugLevel() > 0) fRecoParam.Print();
2893 //_____________________________________________________________________________
2894 Bool_t AliReconstruction::GetEventInfo()
2896 // Fill the event info object
2898 AliCodeTimerAuto("")
2900 AliCentralTrigger *aCTP = NULL;
2902 fEventInfo.SetEventType(fRawReader->GetType());
2904 ULong64_t mask = fRawReader->GetClassMask();
2905 fEventInfo.SetTriggerMask(mask);
2906 UInt_t clmask = fRawReader->GetDetectorPattern()[0];
2907 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(clmask));
2909 aCTP = new AliCentralTrigger();
2910 TString configstr("");
2911 if (!aCTP->LoadConfiguration(configstr)) { // Load CTP config from OCDB
2912 AliError("No trigger configuration found in OCDB! The trigger configuration information will not be used!");
2916 aCTP->SetClassMask(mask);
2917 aCTP->SetClusterMask(clmask);
2920 fEventInfo.SetEventType(AliRawEventHeaderBase::kPhysicsEvent);
2922 if (fRunLoader && (!fRunLoader->LoadTrigger())) {
2923 aCTP = fRunLoader->GetTrigger();
2924 fEventInfo.SetTriggerMask(aCTP->GetClassMask());
2925 fEventInfo.SetTriggerCluster(AliDAQ::ListOfTriggeredDetectors(aCTP->GetClusterMask()));
2928 AliWarning("No trigger can be loaded! The trigger information will not be used!");
2933 AliTriggerConfiguration *config = aCTP->GetConfiguration();
2935 AliError("No trigger configuration has been found! The trigger configuration information will not be used!");
2936 if (fRawReader) delete aCTP;
2940 UChar_t clustmask = 0;
2942 ULong64_t trmask = fEventInfo.GetTriggerMask();
2943 const TObjArray& classesArray = config->GetClasses();
2944 Int_t nclasses = classesArray.GetEntriesFast();
2945 for( Int_t iclass=0; iclass < nclasses; iclass++ ) {
2946 AliTriggerClass* trclass = (AliTriggerClass*)classesArray.At(iclass);
2948 Int_t trindex = (Int_t)TMath::Log2(trclass->GetMask());
2949 fesd->SetTriggerClass(trclass->GetName(),trindex);
2950 if (trmask & (1 << trindex)) {
2952 trclasses += trclass->GetName();
2954 clustmask |= trclass->GetCluster()->GetClusterMask();
2958 fEventInfo.SetTriggerClasses(trclasses);
2960 // Set the information in ESD
2961 fesd->SetTriggerMask(trmask);
2962 fesd->SetTriggerCluster(clustmask);
2964 if (!aCTP->CheckTriggeredDetectors()) {
2965 if (fRawReader) delete aCTP;
2969 if (fRawReader) delete aCTP;
2971 // We have to fill also the HLT decision here!!
2977 const char *AliReconstruction::MatchDetectorList(const char *detectorList, UInt_t detectorMask)
2979 // Match the detector list found in the rec.C or the default 'ALL'
2980 // to the list found in the GRP (stored there by the shuttle PP which
2981 // gets the information from ECS)
2982 static TString resultList;
2983 TString detList = detectorList;
2987 for(Int_t iDet = 0; iDet < (AliDAQ::kNDetectors-1); iDet++) {
2988 if ((detectorMask >> iDet) & 0x1) {
2989 TString det = AliDAQ::OfflineModuleName(iDet);
2990 if ((detList.CompareTo("ALL") == 0) ||
2991 detList.BeginsWith("ALL ") ||
2992 detList.EndsWith(" ALL") ||
2993 detList.Contains(" ALL ") ||
2994 (detList.CompareTo(det) == 0) ||
2995 detList.BeginsWith(det) ||
2996 detList.EndsWith(det) ||
2997 detList.Contains( " "+det+" " )) {
2998 if (!resultList.EndsWith(det + " ")) {
3007 if ((detectorMask >> AliDAQ::kHLTId) & 0x1) {
3008 TString hltDet = AliDAQ::OfflineModuleName(AliDAQ::kNDetectors-1);
3009 if ((detList.CompareTo("ALL") == 0) ||
3010 detList.BeginsWith("ALL ") ||
3011 detList.EndsWith(" ALL") ||
3012 detList.Contains(" ALL ") ||
3013 (detList.CompareTo(hltDet) == 0) ||
3014 detList.BeginsWith(hltDet) ||
3015 detList.EndsWith(hltDet) ||
3016 detList.Contains( " "+hltDet+" " )) {
3017 resultList += hltDet;
3021 return resultList.Data();
3025 //______________________________________________________________________________
3026 void AliReconstruction::Abort(const char *method, EAbort what)
3028 // Abort processing. If what = kAbortProcess, the Process() loop will be
3029 // aborted. If what = kAbortFile, the current file in a chain will be
3030 // aborted and the processing will continue with the next file, if there
3031 // is no next file then Process() will be aborted. Abort() can also be
3032 // called from Begin(), SlaveBegin(), Init() and Notify(). After abort
3033 // the SlaveTerminate() and Terminate() are always called. The abort flag
3034 // can be checked in these methods using GetAbort().
3036 // The method is overwritten in AliReconstruction for better handling of
3037 // reco specific errors
3039 if (!fStopOnError) return;
3043 TString whyMess = method;
3044 whyMess += " failed! Aborting...";
3046 AliError(whyMess.Data());
3049 TString mess = "Abort";
3050 if (fAbort == kAbortProcess)
3051 mess = "AbortProcess";
3052 else if (fAbort == kAbortFile)
3055 Info(mess, whyMess.Data());